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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillcorio (co, nk, kon, ipkon, lakon, ne, nodeboun, ndirboun, xboun, nboun, ipompc, nodempc, coefmpc, nmpc, nodeforc, ndirforc, xforc, nforc, nelemload, sideload, xload, nload, xbody, ipobody, nbody, cgr, ad, au, nactdof, icol, jq, irow, neq, nzl, nmethod, ikmpc, ilmpc, ikboun, ilboun, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, ielmat, ielorien, norien, orab, ntmat_, t0, t1, ithermal, prestr, iprestr, vold, iperturb, sti, nzs, stx, adb, aub, iexpl, plicon, nplicon, plkcon, nplkcon, xstiff, npmat_, dtime, matname, mi, ncmat_, ttime, time, istep, iinc, ibody, ielprop, prop)
 

Function/Subroutine Documentation

◆ mafillcorio()

subroutine mafillcorio ( real*8, dimension(3,*)  co,
integer  nk,
integer, dimension(*)  kon,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer  ne,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
real*8, dimension(*)  xboun,
integer  nboun,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer, dimension(2,*)  nodeforc,
integer, dimension(*)  ndirforc,
real*8, dimension(*)  xforc,
integer  nforc,
integer, dimension(2,*)  nelemload,
character*20, dimension(*)  sideload,
real*8, dimension(2,*)  xload,
integer  nload,
real*8, dimension(7,*)  xbody,
integer, dimension(2,*)  ipobody,
integer  nbody,
real*8, dimension(4,*)  cgr,
real*8, dimension(*)  ad,
real*8, dimension(*)  au,
integer, dimension(0:mi(2),*)  nactdof,
integer, dimension(*)  icol,
integer, dimension(*)  jq,
integer, dimension(*)  irow,
integer, dimension(2)  neq,
integer  nzl,
integer  nmethod,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
real*8, dimension(0:ncmat_,ntmat_,*)  elcon,
integer, dimension(2,*)  nelcon,
real*8, dimension(0:1,ntmat_,*)  rhcon,
integer, dimension(*)  nrhcon,
real*8, dimension(0:6,ntmat_,*)  alcon,
integer, dimension(2,*)  nalcon,
real*8, dimension(*)  alzero,
integer, dimension(mi(3),*)  ielmat,
integer, dimension(mi(3),*)  ielorien,
integer  norien,
real*8, dimension(7,*)  orab,
integer  ntmat_,
real*8, dimension(*)  t0,
real*8, dimension(*)  t1,
integer  ithermal,
real*8, dimension(6,mi(1),*)  prestr,
integer  iprestr,
real*8, dimension(0:mi(2),*)  vold,
integer  iperturb,
real*8, dimension(6,mi(1),*)  sti,
integer, dimension(3)  nzs,
real*8, dimension(6,mi(1),*)  stx,
real*8, dimension(*)  adb,
real*8, dimension(*)  aub,
integer  iexpl,
real*8, dimension(0:2*npmat_,ntmat_,*)  plicon,
integer, dimension(0:ntmat_,*)  nplicon,
real*8, dimension(0:2*npmat_,ntmat_,*)  plkcon,
integer, dimension(0:ntmat_,*)  nplkcon,
real*8, dimension(27,mi(1),*)  xstiff,
integer  npmat_,
real*8  dtime,
character*80, dimension(*)  matname,
integer, dimension(*)  mi,
integer  ncmat_,
real*8  ttime,
real*8  time,
integer  istep,
integer  iinc,
integer, dimension(3,*)  ibody,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop 
)
30 !
31 ! filling the damping matrix in spare matrix format (sm)
32 !
33  implicit none
34 !
35  character*8 lakon(*)
36  character*20 sideload(*)
37  character*80 matname(*)
38 !
39  integer kon(*),nodeboun(*),ndirboun(*),ipompc(*),nodempc(3,*),
40  & nodeforc(2,*),ndirforc(*),nelemload(2,*),icol(*),jq(*),ikmpc(*),
41  & ilmpc(*),ikboun(*),ilboun(*),mi(*),nactdof(0:mi(2),*),konl(20),
42  & irow(*),nelcon(2,*),nrhcon(*),nalcon(2,*),ielmat(mi(3),*),
43  & ielorien(mi(3),*),ipkon(*),ipobody(2,*),nbody,ibody(3,*),
44  & nk,ne,nboun,nmpc,nforc,nload,neq(2),nzl,nmethod,icolumn,
45  & ithermal,iprestr,iperturb,nzs(3),i,j,k,l,m,idist,jj,
46  & ll,id,id1,id2,ist,ist1,ist2,index,jdof1,jdof2,idof1,idof2,
47  & mpc1,mpc2,index1,index2,node1,node2,ielprop(*),
48  & ntmat_,indexe,nope,norien,iexpl,i0,ncmat_,istep,iinc,
49  & nplicon(0:ntmat_,*),nplkcon(0:ntmat_,*),npmat_
50 !
51  real*8 co(3,*),xboun(*),coefmpc(*),xforc(*),xload(2,*),p1(3),
52  & p2(3),ad(*),au(*),bodyf(3),t0(*),t1(*),prestr(6,mi(1),*),
53  & vold(0:mi(2),*),s(60,60),ff(60),prop(*),
54  & sti(6,mi(1),*),sm(60,60),stx(6,mi(1),*),adb(*),aub(*),
55  & elcon(0:ncmat_,ntmat_,*),rhcon(0:1,ntmat_,*),
56  & alcon(0:6,ntmat_,*),alzero(*),orab(7,*),xbody(7,*),cgr(4,*)
57 !
58  real*8 plicon(0:2*npmat_,ntmat_,*),plkcon(0:2*npmat_,ntmat_,*),
59  & xstiff(27,mi(1),*)
60 !
61  real*8 om,value,dtime,ttime,time
62 !
63  i0=0
64 !
65 ! determining nzl
66 !
67  nzl=0
68  do i=neq(2),1,-1
69  if(icol(i).gt.0) then
70  nzl=i
71  exit
72  endif
73  enddo
74 !
75 ! initializing the matrices
76 !
77  do i=1,neq(2)
78  ad(i)=0.d0
79  enddo
80  do i=1,nzs(2)
81  au(i)=0.d0
82  enddo
83 !
84 ! mechanical analysis: loop over all elements
85 !
86  loop: do i=1,ne
87 !
88  if(ipkon(i).lt.0) cycle
89  indexe=ipkon(i)
90 c Bernhardi start
91  if(lakon(i)(1:5).eq.'C3D8I') then
92  nope=11
93  elseif(lakon(i)(4:4).eq.'2') then
94 c Bernhardi end
95  nope=20
96  elseif(lakon(i)(4:4).eq.'8') then
97  nope=8
98  elseif(lakon(i)(4:5).eq.'10') then
99  nope=10
100  elseif(lakon(i)(4:4).eq.'4') then
101  nope=4
102  elseif(lakon(i)(4:5).eq.'15') then
103  nope=15
104  elseif(lakon(i)(4:4).eq.'6') then
105  nope=6
106  else
107  cycle
108  endif
109 !
110  do j=1,nope
111  konl(j)=kon(indexe+j)
112  enddo
113 !
114  om=0.d0
115 !
116  index=i
117  do
118  j=ipobody(1,index)
119  if(j.eq.0) exit
120  if(ibody(1,j).eq.1) then
121  om=xbody(1,j)
122  p1(1)=xbody(2,j)
123  p1(2)=xbody(3,j)
124  p1(3)=xbody(4,j)
125  p2(1)=xbody(5,j)
126  p2(2)=xbody(6,j)
127  p2(3)=xbody(7,j)
128  exit
129  endif
130  index=ipobody(2,index)
131  if(index.eq.0) cycle loop
132  enddo
133 !
134  call e_corio(co,nk,konl,lakon(i),p1,p2,om,bodyf,nbody,s,sm,ff,i,
135  & elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
136  & alzero,ielmat,ielorien,norien,orab,ntmat_,
137  & t0,t1,ithermal,vold,iperturb,
138  & nelemload,sideload,xload,nload,idist,sti,stx,
139  & iexpl,plicon,nplicon,plkcon,nplkcon,xstiff,npmat_,
140  & dtime,matname,mi(1),ncmat_,ttime,time,istep,iinc,
141  & nmethod,ielprop,prop)
142 !
143  do jj=1,3*nope
144 !
145  j=(jj-1)/3+1
146  k=jj-3*(j-1)
147 !
148  node1=kon(indexe+j)
149  jdof1=nactdof(k,node1)
150 !
151  do ll=jj,3*nope
152 !
153  l=(ll-1)/3+1
154  m=ll-3*(l-1)
155 !
156  node2=kon(indexe+l)
157  jdof2=nactdof(m,node2)
158 !
159 ! check whether one of the DOF belongs to a SPC or MPC
160 !
161  if((jdof1.gt.0).and.(jdof2.gt.0)) then
162  call add_sm_st_corio(au,ad,jq,irow,jdof1,jdof2,
163  & s(jj,ll),jj,ll)
164  elseif((jdof1.gt.0).or.(jdof2.gt.0)) then
165 !
166 ! idof1: genuine DOF
167 ! idof2: nominal DOF of the SPC/MPC
168 !
169  if(jdof1.le.0) then
170  idof1=jdof2
171 c idof2=(node1-1)*8+k
172  idof2=jdof1
173  else
174  idof1=jdof1
175 c idof2=(node2-1)*8+m
176  idof2=jdof2
177  endif
178  if(nmpc.gt.0) then
179 c call nident(ikmpc,idof2,nmpc,id)
180 c if((id.gt.0).and.(ikmpc(id).eq.idof2)) then
181  if(idof2.ne.2*(idof2/2)) then
182 !
183 ! regular DOF / MPC
184 !
185 c id=ilmpc(id)
186  id=(-idof2+1)/2
187  ist=ipompc(id)
188  index=nodempc(3,ist)
189  if(index.eq.0) cycle
190  do
191  idof2=nactdof(nodempc(2,index),nodempc(1,index))
192  value=-coefmpc(index)*s(jj,ll)/coefmpc(ist)
193  if(idof1.eq.idof2) value=2.d0*value
194  if(idof2.gt.0) then
195  call add_sm_st_corio(au,ad,jq,irow,idof1,
196  & idof2,value,i0,i0)
197  endif
198  index=nodempc(3,index)
199  if(index.eq.0) exit
200  enddo
201  cycle
202  endif
203  endif
204  else
205 c idof1=(node1-1)*8+k
206 c idof2=(node2-1)*8+m
207  idof1=jdof1
208  idof2=jdof2
209  mpc1=0
210  mpc2=0
211  if(nmpc.gt.0) then
212 c call nident(ikmpc,idof1,nmpc,id1)
213 c if((id1.gt.0).and.(ikmpc(id1).eq.idof1)) mpc1=1
214 c call nident(ikmpc,idof2,nmpc,id2)
215 c if((id2.gt.0).and.(ikmpc(id2).eq.idof2)) mpc2=1
216  if(idof1.ne.2*(idof1/2)) mpc1=1
217  if(idof2.ne.2*(idof2/2)) mpc2=1
218  endif
219  if((mpc1.eq.1).and.(mpc2.eq.1)) then
220 c id1=ilmpc(id1)
221 c id2=ilmpc(id2)
222  id1=(-idof1+1)/2
223  id2=(-idof2+1)/2
224  if(id1.eq.id2) then
225 !
226 ! MPC id1 / MPC id1
227 !
228  ist=ipompc(id1)
229  index1=nodempc(3,ist)
230  if(index1.eq.0) cycle
231  do
232  idof1=nactdof(nodempc(2,index1),
233  & nodempc(1,index1))
234  index2=index1
235  do
236  idof2=nactdof(nodempc(2,index2),
237  & nodempc(1,index2))
238  value=coefmpc(index1)*coefmpc(index2)*
239  & s(jj,ll)/coefmpc(ist)/coefmpc(ist)
240  if((idof1.gt.0).and.(idof2.gt.0)) then
241  call add_sm_st_corio(au,ad,jq,irow,
242  & idof1,idof2,value,i0,i0)
243  endif
244 !
245  index2=nodempc(3,index2)
246  if(index2.eq.0) exit
247  enddo
248  index1=nodempc(3,index1)
249  if(index1.eq.0) exit
250  enddo
251  else
252 !
253 ! MPC id1 / MPC id2
254 !
255  ist1=ipompc(id1)
256  index1=nodempc(3,ist1)
257  if(index1.eq.0) cycle
258  do
259  idof1=nactdof(nodempc(2,index1),
260  & nodempc(1,index1))
261  ist2=ipompc(id2)
262  index2=nodempc(3,ist2)
263  if(index2.eq.0) then
264  index1=nodempc(3,index1)
265  if(index1.eq.0) then
266  exit
267  else
268  cycle
269  endif
270  endif
271  do
272  idof2=nactdof(nodempc(2,index2),
273  & nodempc(1,index2))
274  value=coefmpc(index1)*coefmpc(index2)*
275  & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2)
276  if(idof1.eq.idof2) value=2.d0*value
277  if((idof1.gt.0).and.(idof2.gt.0)) then
278  call add_sm_st_corio(au,ad,jq,irow,
279  & idof1,idof2,value,i0,i0)
280  endif
281 !
282  index2=nodempc(3,index2)
283  if(index2.eq.0) exit
284  enddo
285  index1=nodempc(3,index1)
286  if(index1.eq.0) exit
287  enddo
288  endif
289  endif
290  endif
291  enddo
292 !
293  enddo
294  enddo loop
295 !
296  return
subroutine e_corio(co, nk, konl, lakonl, p1, p2, omx, bodyfx, nbody, s, sm, ff, nelem, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, ielmat, ielorien, norien, orab, ntmat_, t0, t1, ithermal, vold, iperturb, nelemload, sideload, xload, nload, idist, sti, stx, iexpl, plicon, nplicon, plkcon, nplkcon, xstiff, npmat_, dtime, matname, mi, ncmat_, ttime, time, istep, iinc, nmethod, ielprop, prop)
Definition: e_corio.f:26
static ITG * idist
Definition: radflowload.c:39
subroutine add_sm_st_corio(au, ad, jq, irow, i, j, value, i0, i1)
Definition: add_sm_st_corio.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)