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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillsmcsse (co, kon, ipkon, lakon, ne, ipompc, nodempc, coefmpc, nmpc, nelemload, sideload, xload, nload, xbody, ipobody, nbody, cgr, nactdof, neq, nmethod, ikmpc, ilmpc, elcon, nelcon, rhcon, nrhcon, alcon, nalcon, alzero, ielmat, ielorien, norien, orab, ntmat_, t0, t1, ithermal, iprestr, vold, iperturb, sti, stx, iexpl, plicon, nplicon, plkcon, nplkcon, xstiff, npmat_, dtime, matname, mi, ncmat_, mass, stiffness, buckling, rhsi, intscheme, physcon, 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, nea, neb, distmin, ndesi, nodedesi, df, jqs, irows, dfl, icoordinate, dxstiff, xdesi, istartelem, ialelem, v, sigma, labmpc, ics, cs, mcs, nk, nzss)
 

Function/Subroutine Documentation

◆ mafillsmcsse()

subroutine mafillsmcsse ( real*8, dimension(3,*), intent(in)  co,
integer, dimension(*), intent(in)  kon,
integer, dimension(*), intent(in)  ipkon,
character*8, dimension(*), intent(in)  lakon,
integer, intent(in)  ne,
integer, dimension(*), intent(in)  ipompc,
integer, dimension(3,*), intent(in)  nodempc,
real*8, dimension(*), intent(in)  coefmpc,
integer, intent(in)  nmpc,
integer, dimension(2,*), intent(in)  nelemload,
character*20, dimension(*), intent(in)  sideload,
real*8, dimension(2,*), intent(inout)  xload,
integer, intent(in)  nload,
real*8, dimension(7,*), intent(in)  xbody,
integer, dimension(2,*), intent(in)  ipobody,
integer, intent(in)  nbody,
real*8, dimension(4,*), intent(inout)  cgr,
integer, dimension(0:mi(2),*), intent(in)  nactdof,
integer, dimension(2), intent(in)  neq,
integer, intent(inout)  nmethod,
integer, dimension(*), intent(in)  ikmpc,
integer, dimension(*), intent(in)  ilmpc,
real*8, dimension(0:ncmat_,ntmat_,*), intent(in)  elcon,
integer, dimension(2,*), intent(in)  nelcon,
real*8, dimension(0:1,ntmat_,*), intent(in)  rhcon,
integer, dimension(*), intent(in)  nrhcon,
real*8, dimension(0:6,ntmat_,*), intent(in)  alcon,
integer, dimension(2,*), intent(in)  nalcon,
real*8, dimension(*), intent(in)  alzero,
integer, dimension(mi(3),*), intent(in)  ielmat,
integer, dimension(mi(3),*), intent(in)  ielorien,
integer, intent(in)  norien,
real*8, dimension(7,*), intent(in)  orab,
integer, intent(in)  ntmat_,
real*8, dimension(*), intent(in)  t0,
real*8, dimension(*), intent(in)  t1,
integer, dimension(2), intent(in)  ithermal,
integer, intent(in)  iprestr,
real*8, dimension(0:mi(2),*), intent(in)  vold,
integer, dimension(*), intent(in)  iperturb,
real*8, dimension(6,mi(1),*), intent(in)  sti,
real*8, dimension(6,mi(1),*), intent(in)  stx,
integer, intent(in)  iexpl,
real*8, dimension(0:2*npmat_,ntmat_,*), intent(in)  plicon,
integer, dimension(0:ntmat_,*), intent(in)  nplicon,
real*8, dimension(0:2*npmat_,ntmat_,*), intent(in)  plkcon,
integer, dimension(0:ntmat_,*), intent(in)  nplkcon,
real*8, dimension(27,mi(1),*), intent(in)  xstiff,
integer, intent(in)  npmat_,
real*8, intent(in)  dtime,
character*80, dimension(*), intent(in)  matname,
integer, dimension(*), intent(in)  mi,
integer, intent(in)  ncmat_,
integer, dimension(2), intent(in)  mass,
integer, intent(in)  stiffness,
integer, intent(in)  buckling,
integer, intent(in)  rhsi,
integer, intent(in)  intscheme,
real*8, dimension(*), intent(in)  physcon,
real*8, intent(in)  ttime,
real*8, intent(in)  time,
integer, intent(in)  istep,
integer, intent(in)  iinc,
integer, intent(in)  coriolis,
integer, dimension(3,*), intent(in)  ibody,
real*8, dimension(2,*), intent(in)  xloadold,
real*8, intent(in)  reltime,
real*8, dimension(0:mi(2),*), intent(in)  veold,
real*8, dimension(2,*), intent(inout)  springarea,
integer, intent(in)  nstate_,
real*8, dimension(nstate_,mi(1),*), intent(in)  xstateini,
real*8, dimension(nstate_,mi(1),*), intent(inout)  xstate,
real*8, dimension(mi(3),*), intent(in)  thicke,
integer, dimension(*), intent(in)  integerglob,
real*8, dimension(*), intent(in)  doubleglob,
character*81, dimension(3,*), intent(in)  tieset,
integer, dimension(*), intent(in)  istartset,
integer, dimension(*), intent(in)  iendset,
integer, dimension(*), intent(in)  ialset,
integer, intent(in)  ntie,
integer, intent(in)  nasym,
real*8, dimension(3,*), intent(in)  pslavsurf,
real*8, dimension(6,*), intent(in)  pmastsurf,
integer, intent(in)  mortar,
real*8, dimension(3,9,*), intent(in)  clearini,
integer, dimension(*), intent(in)  ielprop,
real*8, dimension(*), intent(in)  prop,
integer, intent(in)  ne0,
integer, intent(in)  nea,
integer, intent(in)  neb,
real*8, intent(in)  distmin,
integer, intent(in)  ndesi,
integer, dimension(*), intent(in)  nodedesi,
real*8, dimension(*), intent(inout)  df,
integer, dimension(*), intent(in)  jqs,
integer, dimension(*), intent(in)  irows,
real*8, dimension(ndesi,60), intent(inout)  dfl,
integer, intent(in)  icoordinate,
real*8, dimension(27,mi(1),ne,*), intent(in)  dxstiff,
real*8, dimension(3,*), intent(in)  xdesi,
integer, dimension(*), intent(in)  istartelem,
integer, dimension(*), intent(in)  ialelem,
real*8, dimension(0:mi(2),*), intent(in)  v,
real*8  sigma,
character*20, dimension(*)  labmpc,
integer, dimension(*)  ics,
real*8, dimension(17,*)  cs,
integer  mcs,
integer  nk,
integer  nzss 
)
33 !
34  implicit none
35 !
36  integer mass(2),stiffness,buckling,rhsi,stiffonly(2),coriolis
37 !
38  character*8 lakon(*)
39  character*20 sideload(*),labmpc(*)
40  character*80 matname(*)
41  character*81 tieset(3,*)
42 !
43  integer kon(*),ipompc(*),nodempc(3,*),nelemload(2,*),ikmpc(*),
44  & ilmpc(*),mi(*),nstate_,ne0,nasym,nactdof(0:mi(2),*),ialset(*),
45  & ielprop(*),nelcon(2,*),nrhcon(*),nalcon(2,*),ielmat(mi(3),*),
46  & ntie,ielorien(mi(3),*),integerglob(*),istartset(*),iendset(*),
47  & ipkon(*),intscheme,ipobody(2,*),nbody,ibody(3,*),ne,nmpc,
48  & nload,neq(2),nmethod,ithermal(2),iprestr,iperturb(*),i,j,k,
49  & idist,jj,id,ist,index,jdof1,node1,kflag,icalccg,ntmat_,nk,
50  & indexe,nope,norien,iexpl,i0,ncmat_,istep,iinc,jqs(*),irows(*),
51  & nplicon(0:ntmat_,*),nplkcon(0:ntmat_,*),npmat_,mortar,nea,
52  & neb,ndesi,nodedesi(*),idesvar,istartelem(*),ialelem(*),
53  & icoordinate,ii,inode,lprev,ij,ilength,idof2,id1,
54  & icomplex,ics(*),mcs,nzss
55 !
56  real*8 co(3,*),coefmpc(*),xload(2,*),p1(3),p2(3),bodyf(3),
57  & xloadold(2,*),reltime,t0(*),t1(*),vold(0:mi(2),*),
58  & s(60,60),ff(60),sti(6,mi(1),*),sm(60,60),xdesi(3,*),
59  & stx(6,mi(1),*),elcon(0:ncmat_,ntmat_,*),sigma,
60  & rhcon(0:1,ntmat_,*),springarea(2,*),alcon(0:6,ntmat_,*),
61  & physcon(*),prop(*),xstate(nstate_,mi(1),*),xcos,xsin,
62  & xstateini(nstate_,mi(1),*),alzero(*),orab(7,*),valuer,
63  & xbody(7,*),cgr(4,*),plicon(0:2*npmat_,ntmat_,*),
64  & plkcon(0:2*npmat_,ntmat_,*),xstiff(27,mi(1),*),cs(17,*),
65  & veold(0:mi(2),*),om,dtime,ttime,time,thicke(mi(3),*),
66  & doubleglob(*),clearini(3,9,*),pslavsurf(3,*),valuei,
67  & pmastsurf(6,*),distmin,dfl(ndesi,60),waluer,waluei,
68  & df(*),dxstiff(27,mi(1),ne,*),v(0:mi(2),*)
69 !
70  intent(in) co,kon,ipkon,lakon,ne,ipompc,nodempc,coefmpc,nmpc,
71  & nelemload,sideload,nload,xbody,ipobody,nbody,nactdof,neq,
72  & ikmpc,ilmpc,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,alzero,
73  & ielmat,ielorien,norien,orab,ntmat_,t0,t1,ithermal,iprestr,
74  & vold,iperturb,sti,stx,iexpl,plicon,nplicon,plkcon,nplkcon,
75  & xstiff,npmat_,dtime,matname,mi,ncmat_,mass,stiffness,
76  & buckling,rhsi,intscheme,physcon,ttime,time,istep,iinc,
77  & coriolis,ibody,xloadold,reltime,veold,nstate_,xstateini,
78  & thicke,integerglob,doubleglob,tieset,istartset,iendset,
79  & ialset,ntie,nasym,pslavsurf,pmastsurf,mortar,clearini,
80  & ielprop,prop,ne0,nea,neb,ndesi,nodedesi,distmin,
81  & icoordinate,jqs,irows,dxstiff,xdesi,istartelem,ialelem,v
82 !
83  intent(inout) xload,nmethod,cgr,springarea,xstate,df,dfl
84 !
85  kflag=2
86  i0=0
87  icalccg=0
88 !
89  if((stiffness.eq.1).and.(mass(1).eq.0).and.(buckling.eq.0)) then
90  stiffonly(1)=1
91  else
92  stiffonly(1)=0
93  endif
94  if((stiffness.eq.1).and.(mass(2).eq.0).and.(buckling.eq.0)) then
95  stiffonly(2)=1
96  else
97  stiffonly(2)=0
98  endif
99 !
100  if(rhsi.eq.1) then
101 !
102 ! distributed forces (body forces or thermal loads or
103 ! residual stresses or distributed face loads)
104 !
105  if((nbody.ne.0).or.(ithermal(1).ne.0).or.
106  & (iprestr.ne.0).or.(nload.ne.0)) then
107  idist=1
108  else
109  idist=0
110  endif
111 !
112  endif
113 !
114  if((ithermal(1).le.1).or.(ithermal(1).eq.3)) then
115 !
116 ! mechanical analysis: loop over all elements
117 !
118  do i=nea,neb
119 !
120  if(ipkon(i).lt.0) cycle
121  indexe=ipkon(i)
122 c Bernhardi start
123  if(lakon(i)(1:5).eq.'C3D8I') then
124  nope=11
125  elseif(lakon(i)(4:5).eq.'20') then
126 c Bernhardi end
127  nope=20
128  elseif(lakon(i)(4:4).eq.'8') then
129  nope=8
130  elseif(lakon(i)(4:5).eq.'10') then
131  nope=10
132  elseif(lakon(i)(4:4).eq.'4') then
133  nope=4
134  elseif(lakon(i)(4:5).eq.'15') then
135  nope=15
136  elseif(lakon(i)(4:4).eq.'6') then
137  nope=6
138  elseif((lakon(i)(1:2).eq.'ES').and.(lakon(i)(7:7).ne.'F'))
139  & then
140 !
141 ! spring and contact spring elements (NO dashpot elements
142 ! = ED... elements)
143 !
144  nope=ichar(lakon(i)(8:8))-47
145 !
146 ! local contact spring number
147 ! if friction is involved, the contact spring element
148 ! matrices are determined in mafillsmas.f
149 !
150  if(lakon(i)(7:7).eq.'C') then
151  if(nasym.eq.1) cycle
152  if(mortar.eq.1) nope=kon(indexe)
153  endif
154  else
155  cycle
156  endif
157 !
158  om=0.d0
159 !
160  if((nbody.gt.0).and.(lakon(i)(1:1).ne.'E')) then
161 !
162 ! assigning centrifugal forces
163 !
164  bodyf(1)=0.
165  bodyf(2)=0.
166  bodyf(3)=0.
167 !
168  index=i
169  do
170  j=ipobody(1,index)
171  if(j.eq.0) exit
172  if(ibody(1,j).eq.1) then
173  om=xbody(1,j)
174  p1(1)=xbody(2,j)
175  p1(2)=xbody(3,j)
176  p1(3)=xbody(4,j)
177  p2(1)=xbody(5,j)
178  p2(2)=xbody(6,j)
179  p2(3)=xbody(7,j)
180 !
181 ! assigning gravity forces
182 !
183  elseif(ibody(1,j).eq.2) then
184  bodyf(1)=bodyf(1)+xbody(1,j)*xbody(2,j)
185  bodyf(2)=bodyf(2)+xbody(1,j)*xbody(3,j)
186  bodyf(3)=bodyf(3)+xbody(1,j)*xbody(4,j)
187 !
188 ! assigning newton gravity forces
189 !
190  elseif(ibody(1,j).eq.3) then
191  call newton(icalccg,ne,ipkon,lakon,kon,t0,co,rhcon,
192  & nrhcon,ntmat_,physcon,i,cgr,bodyf,ielmat,
193  & ithermal,vold,mi)
194  endif
195  index=ipobody(2,index)
196  if(index.eq.0) exit
197  enddo
198  endif
199 !
200  call e_c3d_cs_se(co,kon,lakon(i),p1,p2,om,bodyf,nbody,s,sm,
201  & ff,i,nmethod,elcon,nelcon,rhcon,nrhcon,alcon,nalcon,
202  & alzero,ielmat,ielorien,norien,orab,ntmat_,
203  & t0,t1,ithermal,vold,iperturb,nelemload,sideload,xload,
204  & nload,idist,sti,stx,iexpl,plicon,
205  & nplicon,plkcon,nplkcon,xstiff,npmat_,
206  & dtime,matname,mi(1),ncmat_,mass(1),stiffness,buckling,
207  & rhsi,intscheme,ttime,time,istep,iinc,coriolis,xloadold,
208  & reltime,ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
209  & springarea,nstate_,xstateini,xstate,ne0,ipkon,thicke,
210  & integerglob,doubleglob,tieset,istartset,
211  & iendset,ialset,ntie,nasym,pslavsurf,pmastsurf,mortar,
212  & clearini,ielprop,prop,distmin,ndesi,nodedesi,
213  & dfl,icoordinate,dxstiff,ne,xdesi,istartelem,
214  & ialelem,v,sigma,nk)
215 !
216  do ii=istartelem(i),istartelem(i+1)-1
217  idesvar=ialelem(ii)
218  if(idesvar.eq.0) cycle
219 !
220  do jj=1,3*nope
221 !
222  j=(jj-1)/3+1
223  k=jj-3*(j-1)
224 !
225  node1=kon(indexe+j)
226  jdof1=nactdof(k,node1)
227 !
228  if(jdof1.le.0) then
229  idof2=jdof1
230  if(nmpc.gt.0) then
231  if(idof2.ne.2*(idof2/2)) then
232  id1=(-idof2+1)/2
233  ist=ipompc(id1)
234  index=nodempc(3,ist)
235  if(index.eq.0) cycle
236  do
237  inode=nodempc(1,index)
238  icomplex=0
239  if(labmpc(id1)(1:6).eq.'CYCLIC') then
240  read(labmpc(id1)(7:20),'(i14)')icomplex
241  elseif(labmpc(id1)(1:9).eq.'SUBCYCLIC')
242  & then
243  do ij=1,mcs
244  ilength=int(cs(4,ij))
245  lprev=int(cs(14,ij))
246  call nident(ics(lprev+1),inode,
247  & ilength,id)
248  if(id.gt.0) then
249  if(ics(lprev+id).eq.inode) then
250  icomplex=ij
251  exit
252  endif
253  endif
254  enddo
255  endif
256  idof2=nactdof(nodempc(2,index),inode)
257  if(idof2.gt.0) then
258  valuer=-coefmpc(index)*dfl(idesvar,jj)/
259  & coefmpc(ist)
260  valuei=-coefmpc(index)*
261  & dfl(idesvar,60+jj)/coefmpc(ist)
262  if(icomplex.eq.0) then
263  call add_bo_st(df,jqs,irows,idof2,
264  & idesvar,valuer)
265  call add_bo_st(df(nzss+1),jqs,irows,
266  & idof2,idesvar,valuei)
267  else
268  xcos=cs(15,icomplex)
269  xsin=cs(16,icomplex)
270  waluer=valuer*xcos+valuei*xsin
271  waluei=-valuer*xsin+valuei*xcos
272  call add_bo_st(df,jqs,irows,idof2,
273  & idesvar,waluer)
274  call add_bo_st(df(nzss+1),jqs,irows,
275  & idof2,idesvar,waluei)
276  endif
277  endif
278  index=nodempc(3,index)
279  if(index.eq.0) exit
280  enddo
281  endif
282  endif
283  cycle
284  endif
285  call add_bo_st(df,jqs,irows,jdof1,
286  & idesvar,dfl(idesvar,jj))
287  call add_bo_st(df(nzss+1),jqs,irows,jdof1,
288  & idesvar,dfl(idesvar,60+jj))
289  enddo
290  enddo
291  enddo
292  endif
293 !
294  return
subroutine newton(icalccg, ne, ipkon, lakon, kon, t0, co, rhcon, nrhcon, ntmat_, physcon, nelem, cgr, bodyf, ielmat, ithermal, vold, mi)
Definition: newton.f:22
subroutine e_c3d_cs_se(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, distmin, ndesi, nodedesi, dfl, icoordinate, dxstiff, ne, xdesi, istartelem, ialelem, v, sigma, nk)
Definition: e_c3d_cs_se.f:33
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
subroutine add_bo_st(au, jq, irow, i, j, value)
Definition: add_bo_st.f:20
static ITG * idist
Definition: radflowload.c:39
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)