CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
resultsforc_em.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine resultsforc_em (nk, f, fn, nactdof, ipompc, nodempc, coefmpc, labmpc, nmpc, mi, fmpc, calcul_fn, calcul_f, inomat)
 

Function/Subroutine Documentation

◆ resultsforc_em()

subroutine resultsforc_em ( integer  nk,
real*8, dimension(*)  f,
real*8, dimension(0:mi(2),*)  fn,
integer, dimension(0:mi(2),*)  nactdof,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
character*20, dimension(*)  labmpc,
integer  nmpc,
integer, dimension(*)  mi,
real*8, dimension(*)  fmpc,
integer  calcul_fn,
integer  calcul_f,
integer, dimension(*)  inomat 
)
21 !
22 ! calculating the equation system internal force vector
23 ! (one entry for each active degree of freedom)
24 !
25 ! This routine is based on resultsforc. Whereas resultsforc
26 ! applies to mechanical and thermal degrees of freedom,
27 ! resultsforcem is limited to thermal degrees of freedom.
28 ! (application for electromagnetic calculations)
29 !
30  implicit none
31 !
32  character*20 labmpc(*)
33 !
34  integer mi(*),nactdof(0:mi(2),*),ipompc(*),nodempc(3,*),nk,i,j,
35  & nmpc,ist,ndir,node,index,calcul_fn,calcul_f,inomat(*)
36 !
37  real*8 f(*),fn(0:mi(2),*),coefmpc(*),fmpc(*),forcempc
38 !
39 ! subtracting the mpc force (for each linear mpc there is one
40 ! force; the actual force in a node belonging to the mpc is
41 ! obtained by multiplying this force with the nodal coefficient.
42 ! The force has to be subtracted from f, since it does not
43 ! appear on the rhs of the equation system
44 !
45  if(calcul_fn.eq.1)then
46  do i=1,nmpc
47  ist=ipompc(i)
48  node=nodempc(1,ist)
49  if(inomat(node).eq.0) cycle
50  ndir=nodempc(2,ist)
51 c if(ndir.ne.0) cycle
52  forcempc=fn(ndir,node)/coefmpc(ist)
53  fmpc(i)=forcempc
54  fn(ndir,node)=0.d0
55  index=nodempc(3,ist)
56  if(index.eq.0) cycle
57  do
58  node=nodempc(1,index)
59  ndir=nodempc(2,index)
60  fn(ndir,node)=fn(ndir,node)-coefmpc(index)*forcempc
61  index=nodempc(3,index)
62  if(index.eq.0) exit
63  enddo
64  enddo
65  endif
66 !
67 ! calculating the system force vector
68 !
69  if(calcul_f.eq.1) then
70  do i=1,nk
71  if(inomat(i).eq.0) cycle
72  if(nactdof(0,i).gt.0) then
73  f(nactdof(0,i))=fn(0,i)
74  endif
75 c do j=0,mi(2)
76 c if(nactdof(j,i).ne.0) then
77 c f(nactdof(j,i))=fn(j,i)
78 c endif
79 c enddo
80  enddo
81  endif
82 !
83 ! adding the mpc force again to fn
84 !
85  if(calcul_fn.eq.1)then
86  do i=1,nmpc
87  ist=ipompc(i)
88  node=nodempc(1,ist)
89  if(inomat(node).eq.0) cycle
90  ndir=nodempc(2,ist)
91 c if(ndir.ne.0) cycle
92  forcempc=fmpc(i)
93  fn(ndir,node)=forcempc*coefmpc(ist)
94  index=nodempc(3,ist)
95 !
96 ! nodes not belonging to the structure have to be
97 ! taken out
98 !
99  if(labmpc(i)(1:7).eq.'MEANROT') then
100  if(nodempc(3,nodempc(3,index)).eq.0) cycle
101  elseif(labmpc(i)(1:10).eq.'PRETENSION') then
102  if(nodempc(3,index).eq.0) cycle
103  elseif(labmpc(i)(1:5).eq.'RIGID') then
104  if(nodempc(3,nodempc(3,nodempc(3,nodempc(3,nodempc(3,inde
105  &x))))).eq.0) cycle
106  else
107  if(index.eq.0) cycle
108  endif
109  do
110  node=nodempc(1,index)
111  ndir=nodempc(2,index)
112  fn(ndir,node)=fn(ndir,node)+coefmpc(index)*forcempc
113  index=nodempc(3,index)
114  if(labmpc(i)(1:7).eq.'MEANROT') then
115  if(nodempc(3,nodempc(3,index)).eq.0) exit
116  elseif(labmpc(i)(1:10).eq.'PRETENSION') then
117  if(nodempc(3,index).eq.0) exit
118  elseif(labmpc(i)(1:5).eq.'RIGID') then
119  if(nodempc(3,nodempc(3,nodempc(3,nodempc(3,nodempc(3,i
120  &ndex))))).eq.0) exit
121  else
122  if(index.eq.0) exit
123  endif
124  enddo
125  enddo
126  endif
127 !
128  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)