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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillsmcsas (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, fext, 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, ics, cs, nm, ncmat_, labmpc, mass, stiffness, buckling, rhsi, intscheme, mcs, coriolis, ibody, xloadold, reltime, ielcs, veold, springarea, thicke, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nasym, nstate_, xstateini, xstate, pslavsurf, pmastsurf, mortar, clearini, ielprop, prop, ne0, kscale)
 

Function/Subroutine Documentation

◆ mafillsmcsas()

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