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

Go to the source code of this file.

Functions/Subroutines

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)
 

Function/Subroutine Documentation

◆ thermmodel()

subroutine thermmodel ( character*80  amat,
integer  iel,
integer  iint,
integer  kode,
real*8, dimension(*)  coconloc,
real*8, dimension(0:3,3)  vkl,
real*8  dtime,
real*8  time,
real*8  ttime,
integer, dimension(*)  mi,
integer  nstate_,
real*8, dimension(nstate_,mi(1),*)  xstateini,
real*8, dimension(nstate_,mi(1),*)  xstate,
real*8, dimension(3)  qflux,
real*8, dimension(27,mi(1),*)  xstiff,
integer  iorien,
real*8, dimension(3)  pgauss,
real*8, dimension(7,*)  orab,
real*8  t1l,
real*8  t1lold,
real*8, dimension(0:mi(2),*)  vold,
real*8, dimension(3,*)  co,
character*8  lakonl,
integer, dimension(20)  konl,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer  nmethod,
integer, dimension(*)  iperturb 
)
23 !
24  implicit none
25 !
26  character*8 lakonl
27  character*80 amat
28 !
29  integer iel,iint,kode,mi(*),nstate_,iorien,ntgrd,ncoconst,
30  & layer,kspt,kstep,kinc,kal(2,6),konl(20),ipompc(*),nmethod,
31  & nodempc(3,*),nmpc,ikmpc(*),ilmpc(*),j3,j2,j4,jj,j,i,j1,
32  & iperturb(*)
33 !
34  real*8 coconloc(*),vkl(0:3,3),dtime,time,ttime,cond(6),
35  & xstateini(nstate_,mi(1),*),xstate(nstate_,mi(1),*),qflux(3),
36  & pgauss(3),orab(7,*),abqtime(2),u,dudt,dudg(3),dfdt(3),
37  & dfdg(3,3),dtemp,dtemdx(3),predef(1),dpred(1),pnewdt,
38  & skl(3,3),t1lold,xstiff(27,mi(1),*),xa(3,3),vold(0:mi(2),*),
39  & co(3,*),coefmpc(*),t1l
40 !
41  kal=reshape((/1,1,2,2,3,3,1,2,1,3,2,3/),(/2,6/))
42 !
43  if(kode.eq.1) then
44 !
45 ! linear isotropic
46 !
47  do i=1,3
48  cond(i)=coconloc(1)
49  enddo
50  do i=4,6
51  cond(i)=0.d0
52  enddo
53 !
54  do i=1,3
55  qflux(i)=-coconloc(1)*vkl(0,i)
56  enddo
57 !
58  elseif((kode.eq.3).or.(kode.eq.6)) then
59  if((kode.eq.3).and.(iorien.eq.0)) then
60 !
61 ! orthotropic
62 !
63  do i=1,3
64  cond(i)=coconloc(i)
65  enddo
66  do i=4,6
67  cond(i)=0.d0
68  enddo
69 !
70  do i=1,3
71  qflux(i)=-coconloc(i)*vkl(0,i)
72  enddo
73 !
74  else
75  if(iorien.ne.0) then
76 !
77 ! transformation due to special orientation
78 !
79 ! calculating the transformation matrix
80 !
81  call transformatrix(orab(1,iorien),pgauss,skl)
82 !
83 ! modifying the conductivity constants
84 !
85  if(kode.eq.3) then
86  do j=4,6
87  coconloc(j)=0.d0
88  enddo
89  endif
90 !
91  xa(1,1)=coconloc(1)
92  xa(1,2)=coconloc(4)
93  xa(1,3)=coconloc(5)
94  xa(2,1)=coconloc(4)
95  xa(2,2)=coconloc(2)
96  xa(2,3)=coconloc(6)
97  xa(3,1)=coconloc(5)
98  xa(3,2)=coconloc(6)
99  xa(3,3)=coconloc(3)
100 !
101  do jj=1,6
102  coconloc(jj)=0.d0
103  j1=kal(1,jj)
104  j2=kal(2,jj)
105  do j3=1,3
106  do j4=1,3
107  coconloc(jj)=coconloc(jj)+
108  & xa(j3,j4)*skl(j1,j3)*skl(j2,j4)
109  enddo
110  enddo
111  enddo
112  endif
113 !
114 ! anisotropy
115 !
116  do i=1,6
117  cond(i)=coconloc(i)
118  enddo
119 !
120  qflux(1)=-coconloc(1)*vkl(0,1)-coconloc(4)*vkl(0,2)-
121  & coconloc(5)*vkl(0,3)
122  qflux(2)=-coconloc(4)*vkl(0,1)-coconloc(2)*vkl(0,2)-
123  & coconloc(6)*vkl(0,3)
124  qflux(3)=-coconloc(5)*vkl(0,1)-coconloc(6)*vkl(0,2)-
125  & coconloc(3)*vkl(0,3)
126 !
127  endif
128  else
129 !
130 ! user material
131 !
132  ncoconst=-kode-100
133 !
134  do i=1,nstate_
135  xstate(i,iint,iel)=xstateini(i,iint,iel)
136  enddo
137 !
138  abqtime(1)=time-dtime
139  abqtime(2)=ttime+time-dtime
140 !
141  ntgrd=3
142  dtemp=t1l-t1lold
143  do i=1,3
144  dtemdx(i)=vkl(0,i)
145  enddo
146 !
147  call umatht(u,dudt,dudg,qflux,dfdt,dfdg,xstate,t1lold,dtemp,
148  & dtemdx,abqtime,dtime,predef,dpred,amat,ntgrd,nstate_,
149  & coconloc,ncoconst,pgauss,pnewdt,iel,iint,layer,kspt,
150  & kstep,kinc,vold,co,lakonl,konl,
151  & ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,mi)
152 !
153  cond(1)=dfdg(1,1)
154  cond(2)=dfdg(2,2)
155  cond(3)=dfdg(3,3)
156  cond(4)=dfdg(1,2)
157  cond(5)=dfdg(1,3)
158  cond(6)=dfdg(2,3)
159 !
160  endif
161 !
162  if(((nmethod.ne.4).or.(iperturb(1).ne.0)).and.
163  & (nmethod.ne.5)) then
164  do i=1,6
165  xstiff(21+i,iint,iel)=cond(i)
166  enddo
167  endif
168 !
169  return
subroutine umatht(u, dudt, dudg, flux, dfdt, dfdg, statev, temp, dtemp, dtemdx, time, dtime, predef, dpred, cmname, ntgrd, nstatv, props, nprops, coords, pnewdt, noel, npt, layer, kspt, kstep, kinc, vold, co, lakonl, konl, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, mi)
Definition: umatht.f:24
subroutine transformatrix(xab, p, a)
Definition: transformatrix.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)