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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillsmas (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, bb, 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_, mass, stiffness, buckling, rhsi, intscheme, physcon, shcon, nshcon, cocon, ncocon, ttime, time, istep, iinc, coriolis, ibody, xloadold, reltime, veold, springarea, nstate_, xstateini, xstate, thicke, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nasym, pslavsurf, pmastsurf, mortar, clearini, ielprop, prop, ne0, kscale, iponoel, inoel, network)
 

Function/Subroutine Documentation

◆ mafillsmas()

subroutine mafillsmas ( 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,
real*8, dimension(*)  bb,
integer, dimension(0:mi(2),*)  nactdof,
integer, dimension(*)  icol,
integer, dimension(*)  jq,
integer, dimension(*)  irow,
integer  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, dimension(2)  ithermal,
real*8, dimension(6,mi(1),*)  prestr,
integer  iprestr,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  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_,
logical, dimension(2)  mass,
logical  stiffness,
logical  buckling,
logical  rhsi,
integer  intscheme,
real*8, dimension(*)  physcon,
real*8, dimension(0:3,ntmat_,*)  shcon,
integer, dimension(*)  nshcon,
real*8, dimension(0:6,ntmat_,*)  cocon,
integer, dimension(2,*)  ncocon,
real*8  ttime,
real*8  time,
integer  istep,
integer  iinc,
logical  coriolis,
integer, dimension(3,*)  ibody,
real*8, dimension(2,*)  xloadold,
real*8  reltime,
real*8, dimension(0:mi(2),*)  veold,
real*8, dimension(2,*)  springarea,
integer  nstate_,
real*8, dimension(nstate_,mi(1),*)  xstateini,
real*8, dimension(nstate_,mi(1),*)  xstate,
real*8, dimension(mi(3),*)  thicke,
integer, dimension(*)  integerglob,
real*8, dimension(*)  doubleglob,
character*81, dimension(3,*)  tieset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  ntie,
integer  nasym,
real*8, dimension(3,*)  pslavsurf,
real*8, dimension(6,*)  pmastsurf,
integer  mortar,
real*8, dimension(3,9,*)  clearini,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
integer  ne0,
integer  kscale,
integer, dimension(*)  iponoel,
integer, dimension(2,*)  inoel,
integer  network 
)
36 !
37 ! filling the stiffness matrix in spare matrix format (sm)
38 ! asymmetric contributions
39 !
40  implicit none
41 !
42  logical mass(2),stiffness,buckling,rhsi,coriolis
43 !
44  character*8 lakon(*)
45  character*20 sideload(*)
46  character*80 matname(*)
47  character*81 tieset(3,*)
48 !
49  integer kon(*),nodeboun(*),ndirboun(*),ipompc(*),nodempc(3,*),
50  & nodeforc(2,*),ndirforc(*),nelemload(2,*),icol(*),jq(*),ikmpc(*),
51  & ilmpc(*),ikboun(*),ilboun(*),mi(*),integerglob(*),ist,mpc1,
52  & nactdof(0:mi(2),*),irow(*),istartset(*),iendset(*),
53  & nelcon(2,*),nrhcon(*),nalcon(2,*),ielmat(mi(3),*),index,
54  & ielorien(mi(3),*),ialset(*),ntie,ne0,nstate_,index1,index2,
55  & ipkon(*),intscheme,ncocon(2,*),nshcon(*),ipobody(2,*),nbody,
56  & ibody(3,*),nk,ne,nboun,nmpc,nforc,nload,neq,nzl,nmethod,
57  & ithermal(2),iprestr,iperturb(*),nzs(3),i,j,k,l,m,idist,jj,
58  & ll,jdof1,jdof2,node1,node2,id,i0,id1,id2,idof1,idof2,nasym,
59  & ntmat_,indexe,nope,norien,iexpl,ncmat_,istep,iinc,mpc2,
60  & nplicon(0:ntmat_,*),nplkcon(0:ntmat_,*),npmat_,ist1,ist2,
61  & mortar,ielprop(*),kscale,iponoel(*),inoel(2,*),network
62 !
63  real*8 co(3,*),xboun(*),coefmpc(*),xforc(*),xload(2,*),p1(3),
64  & p2(3),ad(*),au(*),bodyf(3),bb(*),xloadold(2,*),value,
65  & t0(*),t1(*),prestr(6,mi(1),*),vold(0:mi(2),*),s(60,60),
66  & ff(60),
67  & sti(6,mi(1),*),sm(60,60),stx(6,mi(1),*),adb(*),aub(*),
68  & elcon(0:ncmat_,ntmat_,*),rhcon(0:1,ntmat_,*),reltime,
69  & alcon(0:6,ntmat_,*),physcon(*),cocon(0:6,ntmat_,*),
70  & shcon(0:3,ntmat_,*),alzero(*),orab(7,*),xbody(7,*),cgr(4,*),
71  & xstate(nstate_,mi(1),*),xstateini(nstate_,mi(1),*),
72  & springarea(2,*),thicke(mi(3),*),clearini(3,9,*),prop(*),
73  & plicon(0:2*npmat_,ntmat_,*),plkcon(0:2*npmat_,ntmat_,*),
74  & xstiff(27,mi(1),*),veold(0:mi(2),*),doubleglob(*),
75  & om,dtime,ttime,time,pslavsurf(3,*),pmastsurf(6,*)
76 !
77  i0=0
78 !
79 ! storing the symmetric matrix in asymmetric format
80 !
81  do i=1,nzs(2)
82  au(nzs(3)+i)=au(i)
83  enddo
84 !
85 ! dynamic applications
86 !
87  if((nmethod.eq.2).or.(nmethod.eq.4)) then
88  do i=1,nzs(2)
89  aub(nzs(3)+i)=aub(i)
90  enddo
91  endif
92 !
93  if(rhsi) then
94 !
95 ! distributed forces (body forces or thermal loads or
96 ! residual stresses or distributed face loads)
97 !
98  if((nbody.ne.0).or.(ithermal(1).ne.0).or.
99  & (iprestr.ne.0).or.(nload.ne.0)) then
100  idist=1
101  else
102  idist=0
103  endif
104 !
105  endif
106 !
107  if((ithermal(1).le.1).or.(ithermal(1).eq.3)) then
108 !
109 ! mechanical analysis: asymmetric contributions
110 ! only contact friction
111 !
112  do i=ne0+1,ne
113 !
114  if(ipkon(i).lt.0) cycle
115  if(lakon(i)(1:2).ne.'ES') cycle
116  indexe=ipkon(i)
117  read(lakon(i)(8:8),'(i1)') nope
118  nope=nope+1
119 !
120 ! local contact spring number
121 !
122  if(lakon(i)(7:7).eq.'C') then
123  if(mortar.eq.1) nope=kon(indexe)
124  else
125  cycle
126  endif
127 !
128  call e_c3d(co,kon,lakon(i),p1,p2,om,bodyf,nbody,s,sm,ff,i,
129  & nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
130  & alzero,ielmat,ielorien,norien,orab,ntmat_,
131  & t0,t1,ithermal,vold,iperturb,nelemload,sideload,xload,
132  & nload,idist,sti,stx,iexpl,plicon,
133  & nplicon,plkcon,nplkcon,xstiff,npmat_,
134  & dtime,matname,mi(1),ncmat_,mass,stiffness,buckling,rhsi,
135  & intscheme,ttime,time,istep,iinc,coriolis,xloadold,
136  & reltime,ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
137  & springarea,nstate_,xstateini,xstate,ne0,ipkon,thicke,
138  & integerglob,
139  & doubleglob,tieset,istartset,iendset,ialset,ntie,nasym,
140  & pslavsurf,pmastsurf,mortar,clearini,ielprop,prop,kscale)
141 !
142  do jj=1,3*nope
143 !
144  j=(jj-1)/3+1
145  k=jj-3*(j-1)
146 !
147  node1=kon(indexe+j)
148  jdof1=nactdof(k,node1)
149 !
150  do ll=1,3*nope
151 !
152  l=(ll-1)/3+1
153  m=ll-3*(l-1)
154 !
155  node2=kon(indexe+l)
156  jdof2=nactdof(m,node2)
157 !
158 ! check whether one of the DOF belongs to a SPC or MPC
159 !
160  if((jdof1.gt.0).and.(jdof2.gt.0)) then
161  call add_sm_st_as(au,ad,jq,irow,jdof1,jdof2,
162  & s(jj,ll),jj,ll,nzs)
163  elseif((jdof1.gt.0).or.(jdof2.gt.0)) then
164 !
165 ! idof1: genuine DOF
166 ! idof2: nominal DOF of the SPC/MPC
167 !
168  if(jdof1.le.0) then
169 c idof1=(node1-1)*8+k
170  idof1=jdof1
171  idof2=jdof2
172  if(nmpc.gt.0) then
173 c call nident(ikmpc,idof1,nmpc,id)
174 c if((id.gt.0).and.(ikmpc(id).eq.idof1)) then
175  if(idof1.ne.2*(idof1/2)) then
176 !
177 ! regular DOF / MPC
178 !
179 c id=ilmpc(id)
180  id=(-idof1+1)/2
181  ist=ipompc(id)
182  index=nodempc(3,ist)
183  if(index.eq.0) cycle
184  do
185  idof1=nactdof(nodempc(2,index),
186  & nodempc(1,index))
187  value=-coefmpc(index)*s(jj,ll)/coefmpc(ist)
188  if(idof1.gt.0) then
189  call add_sm_st_as(au,ad,jq,irow,idof1,
190  & idof2,value,i0,i0,nzs)
191  endif
192  index=nodempc(3,index)
193  if(index.eq.0) exit
194  enddo
195  cycle
196  endif
197  endif
198  else
199  idof1=jdof1
200 c idof2=(node2-1)*8+m
201  idof2=jdof2
202  if(nmpc.gt.0) then
203 c call nident(ikmpc,idof2,nmpc,id)
204 c if((id.gt.0).and.(ikmpc(id).eq.idof2)) then
205  if(idof2.ne.2*(idof2/2)) then
206 !
207 ! regular DOF / MPC
208 !
209 c id=ilmpc(id)
210  id=(-idof2+1)/2
211  ist=ipompc(id)
212  index=nodempc(3,ist)
213  if(index.eq.0) cycle
214  do
215  idof2=nactdof(nodempc(2,index),
216  & nodempc(1,index))
217  value=-coefmpc(index)*s(jj,ll)/coefmpc(ist)
218  if(idof2.gt.0) then
219  call add_sm_st_as(au,ad,jq,irow,idof1,
220  & idof2,value,i0,i0,nzs)
221  endif
222  index=nodempc(3,index)
223  if(index.eq.0) exit
224  enddo
225  cycle
226  endif
227  endif
228  endif
229  else
230 c idof1=(node1-1)*8+k
231 c idof2=(node2-1)*8+m
232  idof1=jdof1
233  idof2=jdof2
234  mpc1=0
235  mpc2=0
236  if(nmpc.gt.0) then
237 c call nident(ikmpc,idof1,nmpc,id1)
238 c if((id1.gt.0).and.(ikmpc(id1).eq.idof1)) mpc1=1
239 c call nident(ikmpc,idof2,nmpc,id2)
240 c if((id2.gt.0).and.(ikmpc(id2).eq.idof2)) mpc2=1
241  if(idof1.ne.2*(idof1/2)) mpc1=1
242  if(idof2.ne.2*(idof2/2)) mpc2=1
243  endif
244  if((mpc1.eq.1).and.(mpc2.eq.1)) then
245 c id1=ilmpc(id1)
246 c id2=ilmpc(id2)
247  id1=(-idof1+1)/2
248  id2=(-idof2+1)/2
249  if(id1.eq.id2) then
250 !
251 ! MPC id1 / MPC id1
252 !
253  ist1=ipompc(id1)
254  index1=nodempc(3,ist1)
255  if(index1.eq.0) cycle
256  do
257  idof1=nactdof(nodempc(2,index1),
258  & nodempc(1,index1))
259 c index2=index1
260  ist2=ipompc(id2)
261  index2=nodempc(3,ist2)
262  if(index2.eq.0) then
263  index1=nodempc(3,index1)
264  if(index1.eq.0) then
265  exit
266  else
267  cycle
268  endif
269  endif
270  do
271  idof2=nactdof(nodempc(2,index2),
272  & nodempc(1,index2))
273 c value=coefmpc(index1)*coefmpc(index2)*
274 c & s(jj,ll)/coefmpc(ist)/coefmpc(ist)
275  value=coefmpc(index1)*coefmpc(index2)*
276  & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2)
277  if((idof1.gt.0).and.(idof2.gt.0)) then
278  call add_sm_st_as(au,ad,jq,irow,
279  & idof1,idof2,value,i0,i0,nzs)
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  else
289 !
290 ! MPC id1 / MPC id2
291 !
292  ist1=ipompc(id1)
293  index1=nodempc(3,ist1)
294  if(index1.eq.0) cycle
295  do
296  idof1=nactdof(nodempc(2,index1),
297  & nodempc(1,index1))
298  ist2=ipompc(id2)
299  index2=nodempc(3,ist2)
300  if(index2.eq.0) then
301  index1=nodempc(3,index1)
302  if(index1.eq.0) then
303  exit
304  else
305  cycle
306  endif
307  endif
308  do
309  idof2=nactdof(nodempc(2,index2),
310  & nodempc(1,index2))
311  value=coefmpc(index1)*coefmpc(index2)*
312  & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2)
313  if((idof1.gt.0).and.(idof2.gt.0)) then
314  call add_sm_st_as(au,ad,jq,irow,
315  & idof1,idof2,value,i0,i0,nzs)
316  endif
317 !
318  index2=nodempc(3,index2)
319  if(index2.eq.0) exit
320  enddo
321  index1=nodempc(3,index1)
322  if(index1.eq.0) exit
323  enddo
324  endif
325  endif
326  endif
327  enddo
328  enddo
329  enddo
330 !
331  endif
332 !
333  if(ithermal(1).gt.1) then
334 !
335 ! thermal analysis: asymmetric contributions
336 !
337  do i=1,ne
338 !
339  if(ipkon(i).lt.0) cycle
340  if((lakon(i)(1:2).ne.'D ').and.
341  & ((lakon(i)(1:1).ne.'D').or.(network.ne.1))) cycle
342  indexe=ipkon(i)
343 !
344 ! no entry or exit elements
345 !
346  if((kon(indexe+1).eq.0).or.(kon(indexe+3).eq.0)) cycle
347  nope=3
348 !
349  call e_c3d_th(co,nk,kon,lakon(i),s,sm,
350  & ff,i,nmethod,rhcon,nrhcon,ielmat,ielorien,norien,orab,
351  & ntmat_,t0,t1,ithermal,vold,iperturb,nelemload,
352  & sideload,xload,nload,idist,iexpl,dtime,
353  & matname,mi(1),mass(2),stiffness,buckling,rhsi,intscheme,
354  & physcon,shcon,nshcon,cocon,ncocon,ttime,time,istep,iinc,
355  & xstiff,xloadold,reltime,ipompc,nodempc,coefmpc,nmpc,ikmpc,
356  & ilmpc,springarea,plkcon,nplkcon,npmat_,ncmat_,elcon,nelcon,
357  & lakon,pslavsurf,pmastsurf,mortar,clearini,plicon,nplicon,ipkon,
358  & ielprop,prop,iponoel,inoel,sti,xstateini,xstate,nstate_,
359  & network,ipobody,xbody,ibody)
360 !
361  do jj=1,nope
362 !
363  j=jj
364 !
365  node1=kon(indexe+j)
366  jdof1=nactdof(0,node1)
367 !
368  do ll=1,nope
369 !
370  l=ll
371 !
372  node2=kon(indexe+l)
373  jdof2=nactdof(0,node2)
374 !
375 ! check whether one of the DOF belongs to a SPC or MPC
376 !
377  if((jdof1.gt.0).and.(jdof2.gt.0)) then
378  call add_sm_st_as(au,ad,jq,irow,jdof1,jdof2,
379  & s(jj,ll),jj,ll,nzs)
380  elseif((jdof1.gt.0).or.(jdof2.gt.0)) then
381 !
382 ! idof1: genuine DOF
383 ! idof2: nominal DOF of the SPC/MPC
384 !
385  if(jdof1.le.0) then
386 c idof1=(node1-1)*8
387  idof1=jdof1
388  idof2=jdof2
389  if(nmpc.gt.0) then
390 c call nident(ikmpc,idof1,nmpc,id)
391 c if((id.gt.0).and.(ikmpc(id).eq.idof1)) then
392  if(idof1.ne.2*(idof1/2)) then
393 !
394 ! regular DOF / MPC
395 !
396 c id=ilmpc(id)
397  id=(-idof1+1)/2
398  ist=ipompc(id)
399  index=nodempc(3,ist)
400  if(index.eq.0) cycle
401  do
402  idof1=nactdof(nodempc(2,index),
403  & nodempc(1,index))
404  value=-coefmpc(index)*s(jj,ll)/coefmpc(ist)
405  if(idof1.gt.0) then
406  call add_sm_st_as(au,ad,jq,irow,idof1,
407  & idof2,value,i0,i0,nzs)
408  endif
409  index=nodempc(3,index)
410  if(index.eq.0) exit
411  enddo
412  cycle
413  endif
414  endif
415  else
416  idof1=jdof1
417 c idof2=(node2-1)*8
418  idof2=jdof2
419  if(nmpc.gt.0) then
420 c call nident(ikmpc,idof2,nmpc,id)
421 c if((id.gt.0).and.(ikmpc(id).eq.idof2)) then
422  if(idof2.ne.2*(idof2/2)) then
423 !
424 ! regular DOF / MPC
425 !
426 c id=ilmpc(id)
427  id=(-idof2+1)/2
428  ist=ipompc(id)
429  index=nodempc(3,ist)
430  if(index.eq.0) cycle
431  do
432  idof2=nactdof(nodempc(2,index),
433  & nodempc(1,index))
434  value=-coefmpc(index)*s(jj,ll)/coefmpc(ist)
435  if(idof2.gt.0) then
436  call add_sm_st_as(au,ad,jq,irow,idof1,
437  & idof2,value,i0,i0,nzs)
438  endif
439  index=nodempc(3,index)
440  if(index.eq.0) exit
441  enddo
442  cycle
443  endif
444  endif
445  endif
446  else
447 c idof1=(node1-1)*8
448 c idof2=(node2-1)*8
449  idof1=jdof1
450  idof2=jdof2
451  mpc1=0
452  mpc2=0
453  if(nmpc.gt.0) then
454 c call nident(ikmpc,idof1,nmpc,id1)
455 c if((id1.gt.0).and.(ikmpc(id1).eq.idof1)) mpc1=1
456 c call nident(ikmpc,idof2,nmpc,id2)
457 c if((id2.gt.0).and.(ikmpc(id2).eq.idof2)) mpc2=1
458  if(idof1.ne.2*(idof1/2)) mpc1=1
459  if(idof2.ne.2*(idof2/2)) mpc2=1
460 c write(*,*) 'mafillsmas ',mpc1,mpc2
461  endif
462  if((mpc1.eq.1).and.(mpc2.eq.1)) then
463 c id1=ilmpc(id1)
464 c id2=ilmpc(id2)
465  id1=(-idof1+1)/2
466  id2=(-idof2+1)/2
467 c write(*,*) 'mafillsmas ',id1,id2
468  if(id1.eq.id2) then
469 !
470 ! MPC id1 / MPC id1
471 !
472  ist1=ipompc(id1)
473  index1=nodempc(3,ist1)
474  if(index1.eq.0) cycle
475  do
476  idof1=nactdof(nodempc(2,index1),
477  & nodempc(1,index1))
478 c index2=index1
479  ist2=ipompc(id2)
480  index2=nodempc(3,ist2)
481  if(index2.eq.0) then
482  index1=nodempc(3,index1)
483  if(index1.eq.0) then
484  exit
485  else
486  cycle
487  endif
488  endif
489  do
490  idof2=nactdof(nodempc(2,index2),
491  & nodempc(1,index2))
492 c value=coefmpc(index1)*coefmpc(index2)*
493 c & s(jj,ll)/coefmpc(ist)/coefmpc(ist)
494  value=coefmpc(index1)*coefmpc(index2)*
495  & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2)
496  if((idof1.gt.0).and.(idof2.gt.0)) then
497  call add_sm_st_as(au,ad,jq,irow,
498  & idof1,idof2,value,i0,i0,nzs)
499  endif
500 !
501  index2=nodempc(3,index2)
502  if(index2.eq.0) exit
503  enddo
504  index1=nodempc(3,index1)
505  if(index1.eq.0) exit
506  enddo
507  else
508 !
509 ! MPC id1 / MPC id2
510 !
511  ist1=ipompc(id1)
512  index1=nodempc(3,ist1)
513  if(index1.eq.0) cycle
514  do
515  idof1=nactdof(nodempc(2,index1),
516  & nodempc(1,index1))
517  ist2=ipompc(id2)
518  index2=nodempc(3,ist2)
519  if(index2.eq.0) then
520  index1=nodempc(3,index1)
521  if(index1.eq.0) then
522  exit
523  else
524  cycle
525  endif
526  endif
527  do
528  idof2=nactdof(nodempc(2,index2),
529  & nodempc(1,index2))
530  value=coefmpc(index1)*coefmpc(index2)*
531  & s(jj,ll)/coefmpc(ist1)/coefmpc(ist2)
532  if((idof1.gt.0).and.(idof2.gt.0)) then
533  call add_sm_st_as(au,ad,jq,irow,
534  & idof1,idof2,value,i0,i0,nzs)
535  endif
536 !
537  index2=nodempc(3,index2)
538  if(index2.eq.0) exit
539  enddo
540  index1=nodempc(3,index1)
541  if(index1.eq.0) exit
542  enddo
543  endif
544  endif
545  endif
546  enddo
547  enddo
548  enddo
549 !
550  endif
551 !
552  return
subroutine e_c3d(co, kon, lakonl, p1, p2, omx, bodyfx, nbody, s, sm, ff, nelem, nmethod, 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_, mass, stiffness, buckling, rhsi, intscheme, ttime, time, istep, iinc, coriolis, xloadold, reltime, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, veold, springarea, nstate_, xstateini, xstate, ne0, ipkon, thicke, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nasym, pslavsurf, pmastsurf, mortar, clearini, ielprop, prop, kscale)
Definition: e_c3d.f:31
static ITG * idist
Definition: radflowload.c:39
subroutine e_c3d_th(co, nk, kon, lakonl, s, sm, ff, nelem, nmethod, rhcon, nrhcon, ielmat, ielorien, norien, orab, ntmat_, t0, t1, ithermal, vold, iperturb, nelemload, sideload, xload, nload, idist, iexpl, dtime, matname, mi, mass, stiffness, buckling, rhsi, intscheme, physcon, shcon, nshcon, cocon, ncocon, ttime, time, istep, iinc, xstiff, xloadold, reltime, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, springarea, plkcon, nplkcon, npmat_, ncmat_, elcon, nelcon, lakon, pslavsurf, pmastsurf, mortar, clearini, plicon, nplicon, ipkon, ielprop, prop, iponoel, inoel, sti, xstateini, xstate, nstate_, network, ipobody, xbody, ibody)
Definition: e_c3d_th.f:31
subroutine add_sm_st_as(au, ad, jq, irow, i, j, value, i0, i1, nzs)
Definition: add_sm_st_as.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)