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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillsmse (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, ieigenfrequency)
 

Function/Subroutine Documentation

◆ mafillsmse()

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