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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ resultsforc()

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