36 character*8 lakon(*),lakonl
37 character*20 sideload(*)
38 character*80 amat,matname(*)
40 integer kon(*),konl(26),iperm(20),ikmpc(*),ilmpc(*),mi(*),
41 & nelcon(2,*),nrhcon(*),ielmat(mi(3),*),ielorien(mi(3),*),
42 & ntmat_,ipkon(*),ipompc(*),nodempc(3,*),mortar,igauss,
43 & ncocon(2,*),iflag,nshcon(*),istep,iinc,mt,mattyp,
44 & i,j,k,m1,kk,i1,m3,indexe,nope,norien,iperturb(*),iout,
45 &
nal,nmpc,kode,imat,mint3d,iorien,istiff,ncmat_,nstate_,
46 & nplkcon(0:ntmat_,*),npmat_,calcul_fn,calcul_qa,nea,neb,
47 & nelemload(2,*),nload,ithermal(2),nmethod,nopered,iloc,
48 & jfaces,node,nplicon(0:ntmat_,*),null,ielprop(*),
49 & iponoel(*),inoel(2,*),network,ipobody(2,*),ibody(3,*)
51 real*8 co(3,*),v(0:mi(2),*),shp(4,26),reltime,
52 & xl(3,26),vl(0:mi(2),26),elcon(0:ncmat_,ntmat_,*),
53 & rhcon(0:1,ntmat_,*),qfx(3,mi(1),*),orab(7,*),
54 & rho,fn(0:mi(2),*),tnl(19),timeend(2),q(0:mi(2),26),
55 & vkl(0:3,3),t0(*),vold(0:mi(2),*),coefmpc(*),
56 & springarea(2,*),elconloc(21),cocon(0:6,ntmat_,*),
57 & shcon(0:3,ntmat_,*),sph,
c1,xi,et,ze,xsj,qa(*),t0l,t1l,dtime,
58 & weight,pgauss(3),coconloc(6),qflux(3),time,ttime,
59 & t1lold,plkcon(0:2*npmat_,ntmat_,*),xstiff(27,mi(1),*),
60 & xstate(nstate_,mi(1),*),xstateini(nstate_,mi(1),*),
61 & xload(2,*),xloadold(2,*),clearini(3,9,*),pslavsurf(3,*),
62 & pmastsurf(6,*),plicon(0:2*npmat_,ntmat_,*),prop(*),
69 iperm=(/5,6,7,8,1,2,3,4,13,14,15,16,9,10,11,12,17,18,19,20/)
79 if(ipkon(i).lt.0) cycle
83 if(lakon(i)(1:1).eq.
'F') cycle
84 if(lakon(i)(1:7).eq.
'DCOUP3D') cycle
95 if(lakon(i)(4:5).eq.
'20')
then 97 elseif(lakon(i)(4:4).eq.
'8')
then 99 elseif(lakon(i)(4:5).eq.
'10')
then 101 elseif(lakon(i)(4:4).eq.
'4')
then 103 elseif(lakon(i)(4:5).eq.
'15')
then 105 elseif(lakon(i)(4:4).eq.
'6')
then 107 elseif((lakon(i)(1:2).eq.
'ES').and.(lakon(i)(7:7).ne.
'A'))
then 112 if(lakon(i)(7:7).eq.
'C')
then 121 elseif(mortar.eq.0)
then 125 nope=ichar(lakon(i)(8:8))-47
126 konl(nope+1)=kon(indexe+nope+1)
132 nope=ichar(lakon(i)(8:8))-47
139 elseif((lakon(i)(1:2).eq.
'D ').or.
140 & ((lakon(i)(1:1).eq.
'D').and.(network.eq.1)))
then 144 if((kon(indexe+1).eq.0).or.(kon(indexe+3).eq.0)) cycle
150 if(lakon(i)(4:5).eq.
'8R')
then 152 elseif(lakon(i)(4:7).eq.
'20RB')
then 153 if((lakon(i)(8:8).eq.
'R').or.(lakon(i)(8:8).eq.
'C'))
then 157 & null,xi,et,ze,weight)
159 elseif((lakon(i)(4:4).eq.
'8').or.
160 & (lakon(i)(4:6).eq.
'20R'))
then 161 if(lakon(i)(6:7).eq.
'RA')
then 166 elseif(lakon(i)(4:4).eq.
'2')
then 168 elseif(lakon(i)(4:5).eq.
'10')
then 170 elseif(lakon(i)(4:4).eq.
'4')
then 172 elseif(lakon(i)(4:5).eq.
'15')
then 174 elseif(lakon(i)(4:4).eq.
'6')
then 181 konl(j)=kon(indexe+j)
183 xl(k,j)=co(k,konl(j))
191 if((iperturb(1).ge.2).or.((iperturb(1).le.0).and.(iout.lt.1)))
194 q(0,m1)=fn(0,konl(m1))
212 if(lakonl(2:2).eq.
'S')
then 213 if(lakonl(7:7).eq.
'C')
then 220 timeend(2)=ttime+time
223 & tnl,ncmat_,ntmat_,nope,kode,elconloc,
224 & plkcon,nplkcon,npmat_,mi,
225 & springarea(1,konl(nope+1)),timeend,matname,
226 & konl(nope),i,istep,iinc,iperturb)
227 elseif(mortar.eq.1)
then 228 iloc=kon(indexe+nope+1)
229 jfaces=kon(indexe+nope+2)
230 igauss=kon(indexe+nope+1)
233 & tnl,ncmat_,ntmat_,nope,lakonl,kode,elconloc,
234 & plicon,nplicon,npmat_,mi,springarea(1,iloc),
235 & nmethod,reltime,jfaces,igauss,
236 & pslavsurf,pmastsurf,clearini,timeend,istep,
237 & iinc,plkcon,nplkcon,node,i,matname)
240 elseif(lakonl(7:7).eq.
'F')
then 244 call advecforc(nope,vl,ithermal,xl,nelemload,
245 & i,nload,lakon,xload,istep,time,ttime,
246 & dtime,sideload,v,mi,xloadold,reltime,nmethod,
247 & tnl,iinc,iponoel,inoel,ielprop,prop,ielmat,shcon,
248 & nshcon,rhcon,nrhcon,ntmat_,ipkon,kon,cocon,
249 & ncocon,ipobody,xbody,ibody)
252 elseif((lakonl(1:2).eq.
'D ').or.
253 & ((lakonl(1:1).eq.
'D').and.(network.eq.1)))
then 258 & nshcon,rhcon,nrhcon)
263 fn(0,konl(j))=fn(0,konl(j))+tnl(j)
268 if(lakon(i)(4:5).eq.
'8R')
then 273 elseif(lakon(i)(4:7).eq.
'20RB')
then 274 if((lakon(i)(8:8).eq.
'R').or.(lakon(i)(8:8).eq.
'C'))
then 278 weight=weight3d13(kk)
281 & prop,kk,xi,et,ze,weight)
283 elseif((lakon(i)(4:4).eq.
'8').or.
284 & (lakon(i)(4:6).eq.
'20R'))
290 elseif(lakon(i)(4:4).eq.
'2')
then 295 elseif(lakon(i)(4:5).eq.
'10')
then 300 elseif(lakon(i)(4:4).eq.
'4')
then 305 elseif(lakon(i)(4:5).eq.
'15')
then 310 elseif(lakon(i)(4:4).eq.
'6')
then 318 if(lakon(i)(7:7).eq.
'A')
then 320 elseif((lakon(i)(7:7).eq.
'E').or.
321 & (lakon(i)(7:7).eq.
'S'))
then 324 call shape20h(xi,et,ze,xl,xsj,shp,iflag)
326 elseif(nope.eq.8)
then 327 call shape8h(xi,et,ze,xl,xsj,shp,iflag)
328 elseif(nope.eq.10)
then 330 elseif(nope.eq.4)
then 331 call shape4tet(xi,et,ze,xl,xsj,shp,iflag)
332 elseif(nope.eq.15)
then 333 call shape15w(xi,et,ze,xl,xsj,shp,iflag)
335 call shape6w(xi,et,ze,xl,xsj,shp,iflag)
349 vkl(0,m3)=vkl(0,m3)+shp(m3,m1)*vl(0,m1)
360 if((lakon(i)(4:5).eq.
'8 ').or.
361 & (lakon(i)(4:5).eq.
'8I'))
then 363 t1lold=t1lold+vold(0,konl(i1))/8.d0
364 t1l=t1l+v(0,konl(i1))/8.d0
366 elseif(lakon(i)(4:6).eq.
'20 ')
then 368 call lintemp_th(t0,vold,konl,nopered,kk,t0l,t1lold,mi)
369 call lintemp_th(t0,v,konl,nopered,kk,t0l,t1l,mi)
370 elseif(lakon(i)(4:6).eq.
'10T')
then 371 call linscal10(vold,konl,t1lold,mi(2),shp)
375 t1lold=t1lold+shp(4,i1)*vold(0,konl(i1))
376 t1l=t1l+shp(4,i1)*v(0,konl(i1))
384 if((iorien.gt.0).or.(kode.le.-100))
then 388 pgauss(j)=pgauss(j)+shp(4,i1)*co(j,konl(i1))
400 & ntmat_,coconloc,mattyp,t1l,rhcon,nrhcon,rho,shcon,
401 & nshcon,sph,xstiff,kk,i,istiff,mi(1))
403 call thermmodel(amat,i,kk,kode,coconloc,vkl,dtime,
404 & time,ttime,mi(1),nstate_,xstateini,xstate,qflux,xstiff,
405 & iorien,pgauss,orab,t1l,t1lold,vold,co,lakon(i),konl,
406 & ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,nmethod,
412 if(lakon(i)(6:7).eq.
'RA')
then 413 qfx(1,kk+4,i)=qflux(1)
414 qfx(2,kk+4,i)=qflux(2)
415 qfx(3,kk+4,i)=qflux(3)
420 if(calcul_fn.eq.1)
then 424 if(lakon(i)(6:7).eq.
'RA')
then 426 fn(0,konl(m1))=fn(0,konl(m1))
427 & -
c1*(qflux(1)*(shp(1,m1)+shp(1,iperm(m1)))
428 & +qflux(2)*(shp(2,m1)+shp(2,iperm(m1)))
429 & +qflux(3)*(shp(3,m1)+shp(3,iperm(m1))))
434 fn(0,konl(m1))=fn(0,konl(m1))-
435 &
c1*qflux(m3)*shp(m3,m1)
448 if(calcul_qa.eq.1)
then 450 qa(2)=qa(2)+dabs(fn(0,konl(m1))-q(0,m1))
subroutine shape20h_pl(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h_pl.f:20
subroutine shape6w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape6w.f:20
subroutine linscal10(scal, konl, scall, idim, shp)
Definition: linscal10.f:20
subroutine springforc_n2f_th(xl, vl, imat, elcon, nelcon, tnl, ncmat_, ntmat_, nope, kode, elconloc, plkcon, nplkcon, npmat_, mi, springarea, timeend, matname, node, noel, istep, iinc, iperturb)
Definition: springforc_n2f_th.f:23
subroutine networkforc(vl, tnl, imat, konl, mi, ntmat_, shcon, nshcon, rhcon, nrhcon)
Definition: networkforc.f:21
subroutine shape10tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape10tet.f:20
subroutine shape8h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape8h.f:20
subroutine springforc_f2f_th(xl, vl, imat, elcon, nelcon, tnl, ncmat_, ntmat_, nope, lakonl, kode, elconloc, plicon, nplicon, npmat_, mi, springarea, nmethod, reltime, jfaces, igauss, pslavsurf, pmastsurf, clearini, timeend, istep, iinc, plkcon, nplkcon, node, noel, matname)
Definition: springforc_f2f_th.f:25
subroutine advecforc(nope, voldl, ithermal, xl, nelemload, nelemadvec, nload, lakon, xload, istep, time, ttime, dtime, sideload, vold, mi, xloadold, reltime, nmethod, tnl, iinc, iponoel, inoel, ielprop, prop, ielmat, shcon, nshcon, rhcon, nrhcon, ntmat_, ipkon, kon, cocon, ncocon, ipobody, ibody, xbody)
Definition: advecforc.f:24
subroutine shape15w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape15w.f:20
static double * c1
Definition: mafillvcompmain.c:30
subroutine shape20h_ax(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h_ax.f:20
subroutine materialdata_th(cocon, ncocon, imat, iorien, pgauss, orab, ntmat_, coconloc, mattyp, t1l, rhcon, nrhcon, rho, shcon, nshcon, sph, xstiff, iint, iel, istiff, mi)
Definition: materialdata_th.f:21
subroutine lintemp_th(t0, vold, konl, nope, jj, t0l, t1l, mi)
Definition: lintemp_th.f:20
subroutine shape20h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h.f:20
subroutine thermmodel(amat, iel, iint, kode, coconloc, vkl, dtime, time, ttime, mi, nstate_, xstateini, xstate, qflux, xstiff, iorien, pgauss, orab, t1l, t1lold, vold, co, lakonl, konl, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, nmethod, iperturb)
Definition: thermmodel.f:23
subroutine beamintscheme(lakonl, mint3d, npropstart, prop, kk, xi, et, ze, weight)
Definition: beamintscheme.f:21
subroutine shape4tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape4tet.f:20
static ITG * nal
Definition: results.c:31