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

Go to the source code of this file.

Functions/Subroutines

subroutine jouleheating (ipkon, lakon, kon, co, elcon, nelcon, mi, ne, sti, ielmat, nelemload, sideload, xload, nload, nload_, iamload, nam, idefload, ncmat_, ntmat_, alcon, nalcon, ithermal, vold, t1)
 

Function/Subroutine Documentation

◆ jouleheating()

subroutine jouleheating ( integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
real*8, dimension(3,*)  co,
real*8, dimension(0:ncmat_,ntmat_,*)  elcon,
integer, dimension(2,*)  nelcon,
integer, dimension(*)  mi,
integer  ne,
real*8, dimension(6,mi(1),*)  sti,
integer, dimension(mi(3),*)  ielmat,
integer, dimension(2,*)  nelemload,
character*20, dimension(*)  sideload,
real*8, dimension(2,*)  xload,
integer  nload,
integer  nload_,
integer, dimension(2,*)  iamload,
integer  nam,
integer, dimension(*)  idefload,
integer  ncmat_,
integer  ntmat_,
real*8, dimension(0:6,ntmat_,*)  alcon,
integer, dimension(2,*)  nalcon,
integer  ithermal,
real*8, dimension(0:mi(2),*)  vold,
real*8, dimension(*)  t1 
)
23 !
24 ! determines the effect of Joule heating
25 !
26  implicit none
27 !
28  character*8 lakon(*)
29  character*20 label,sideload(*)
30 !
31  integer ipkon(*),nelem,kon(*),mi(*),nope,indexe,j,k,null,
32  & mint3d,jj,iflag,ne,nelemload(2,*),iamload(2,*),nload,nload_,
33  & ielmat(mi(3),*),konl(20),idefload(*),iamplitude,isector,nam,
34  & one,nelcon(2,*),nalcon(2,*),ithermal,i1,ncmat_,ntmat_,imat
35 !
36  real*8 co(3,*),xl(3,20),xi,et,ze,xsj,shp(4,20),weight,xload(2,*),
37  & sti(6,mi(1),*),alpha(6),heat,elcon(0:ncmat_,ntmat_,*),volume,
38  & elconloc(21),t1l,alcon(0:6,ntmat_,*),vold(0:mi(2),*),t1(*)
39 !
40  include "gauss.f"
41 !
42  data iflag /2/
43 !
44  null=0
45  one=1
46 !
47  do nelem=1,ne
48  if(ipkon(nelem).lt.0) cycle
49 !
50 ! only elements belonging to the A,V-domain experience
51 ! Joule heating
52 !
53  if(int(elcon(2,1,ielmat(1,nelem))).ne.2) cycle
54 !
55  imat=ielmat(1,nelem)
56  indexe=ipkon(nelem)
57 !
58  if(lakon(nelem)(1:5).eq.'C3D8I') then
59  nope=11
60  elseif(lakon(nelem)(4:4).eq.'2') then
61  nope=20
62  elseif(lakon(nelem)(4:4).eq.'8') then
63  nope=8
64  elseif(lakon(nelem)(4:5).eq.'10') then
65  nope=10
66  elseif(lakon(nelem)(4:4).eq.'4') then
67  nope=4
68  elseif(lakon(nelem)(4:5).eq.'15') then
69  nope=15
70  elseif(lakon(nelem)(4:5).eq.'6') then
71  nope=6
72  endif
73 !
74  do j=1,nope
75  konl(j)=kon(indexe+j)
76  do k=1,3
77  xl(k,j)=co(k,konl(j))
78  enddo
79  enddo
80 !
81  heat=0.d0
82  volume=0.d0
83 !
84  if(lakon(nelem)(4:5).eq.'8R') then
85  mint3d=1
86  elseif((lakon(nelem)(4:4).eq.'8').or.
87  & (lakon(nelem)(4:6).eq.'20R')) then
88  mint3d=8
89  elseif(lakon(nelem)(4:4).eq.'2') then
90  mint3d=27
91  elseif(lakon(nelem)(4:5).eq.'10') then
92  mint3d=4
93  elseif(lakon(nelem)(4:4).eq.'4') then
94  mint3d=1
95  elseif(lakon(nelem)(4:5).eq.'15') then
96  mint3d=9
97  elseif(lakon(nelem)(4:5).eq.'6') then
98  mint3d=2
99  endif
100 !
101 ! loop over the integration points
102 !
103  do jj=1,mint3d
104  if(lakon(nelem)(4:5).eq.'8R') then
105  xi=gauss3d1(1,jj)
106  et=gauss3d1(2,jj)
107  ze=gauss3d1(3,jj)
108  weight=weight3d1(jj)
109  elseif((lakon(nelem)(4:4).eq.'8').or.
110  & (lakon(nelem)(4:6).eq.'20R'))
111  & then
112  xi=gauss3d2(1,jj)
113  et=gauss3d2(2,jj)
114  ze=gauss3d2(3,jj)
115  weight=weight3d2(jj)
116  elseif(lakon(nelem)(4:4).eq.'2') then
117  xi=gauss3d3(1,jj)
118  et=gauss3d3(2,jj)
119  ze=gauss3d3(3,jj)
120  weight=weight3d3(jj)
121  elseif(lakon(nelem)(4:5).eq.'10') then
122  xi=gauss3d5(1,jj)
123  et=gauss3d5(2,jj)
124  ze=gauss3d5(3,jj)
125  weight=weight3d5(jj)
126  elseif(lakon(nelem)(4:4).eq.'4') then
127  xi=gauss3d4(1,jj)
128  et=gauss3d4(2,jj)
129  ze=gauss3d4(3,jj)
130  weight=weight3d4(jj)
131  elseif(lakon(nelem)(4:5).eq.'15') then
132  xi=gauss3d8(1,jj)
133  et=gauss3d8(2,jj)
134  ze=gauss3d8(3,jj)
135  weight=weight3d8(jj)
136  else
137  xi=gauss3d7(1,jj)
138  et=gauss3d7(2,jj)
139  ze=gauss3d7(3,jj)
140  weight=weight3d7(jj)
141  endif
142 !
143  if(nope.eq.20) then
144  call shape20h(xi,et,ze,xl,xsj,shp,iflag)
145  elseif(nope.eq.8) then
146  call shape8h(xi,et,ze,xl,xsj,shp,iflag)
147  elseif(nope.eq.10) then
148  call shape10tet(xi,et,ze,xl,xsj,shp,iflag)
149  elseif(nope.eq.4) then
150  call shape4tet(xi,et,ze,xl,xsj,shp,iflag)
151  elseif(nope.eq.15) then
152  call shape15w(xi,et,ze,xl,xsj,shp,iflag)
153  else
154  call shape6w(xi,et,ze,xl,xsj,shp,iflag)
155  endif
156 !
157 ! calculating the temperature
158 !
159  t1l=0.d0
160  if(ithermal.eq.1) then
161  if(lakon(nelem)(4:5).eq.'8 ') then
162  do i1=1,nope
163  t1l=t1l+t1(konl(i1))/8.d0
164  enddo
165  elseif(lakon(nelem)(4:6).eq.'20 ')then
166  call linscal(t1,konl,nope,jj,t1l,one)
167  elseif(lakon(nelem)(4:6).eq.'10T') then
168  call linscal10(t1,konl,t1l,null,shp)
169  else
170  do i1=1,nope
171  t1l=t1l+shp(4,i1)*t1(konl(i1))
172  enddo
173  endif
174  elseif(ithermal.ge.2) then
175  if(lakon(nelem)(4:5).eq.'8 ') then
176  do i1=1,nope
177  t1l=t1l+vold(0,konl(i1))/8.d0
178  enddo
179  elseif(lakon(nelem)(4:6).eq.'20 ')then
180  call linscal(vold,konl,nope,jj,t1l,mi(2))
181  elseif(lakon(nelem)(4:6).eq.'10T') then
182  call linscal10(vold,konl,t1l,mi(2),shp)
183  else
184  do i1=1,nope
185  t1l=t1l+shp(4,i1)*vold(0,konl(i1))
186  enddo
187  endif
188  endif
189 !
190 ! material data (electric conductivity and
191 ! magnetic permeability)
192 !
193  call materialdata_em(elcon,nelcon,alcon,nalcon,
194  & imat,ntmat_,t1l,elconloc,ncmat_,alpha)
195 !
196  heat=heat+weight*xsj*alpha(1)*
197  & (sti(1,jj,nelem)*sti(1,jj,nelem)+
198  & sti(2,jj,nelem)*sti(2,jj,nelem)+
199  & sti(3,jj,nelem)*sti(3,jj,nelem))
200  volume=volume+weight*xsj
201  enddo
202 !
203  heat=heat/volume
204 !
205 ! adding the Joule heating to the distributed loading
206 !
207  label='BF '
208  iamplitude=0
209  isector=0
210  call loadadd(nelem,label,heat,nelemload,sideload,
211  & xload,nload,nload_,iamload,iamplitude,nam,isector,
212  & idefload)
213  enddo
214 !
215  return
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 loadadd(nelement, label, value, nelemload, sideload, xload, nload, nload_, iamload, iamplitude, nam, isector, idefload)
Definition: loadadd.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 shape15w(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape15w.f:20
subroutine linscal(scal, konl, nope, jj, scall, idim)
Definition: linscal.f:20
subroutine shape20h(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape20h.f:20
subroutine shape4tet(xi, et, ze, xl, xsj, shp, iflag)
Definition: shape4tet.f:20
subroutine materialdata_em(elcon, nelcon, alcon, nalcon, imat, ntmat_, t1l, elconloc, ncmat_, alpha)
Definition: materialdata_em.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)