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

Go to the source code of this file.

Functions/Subroutines

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)
 

Function/Subroutine Documentation

◆ advecforc()

subroutine advecforc ( integer, intent(in)  nope,
real*8, dimension(0:mi(2),20), intent(in)  voldl,
integer, dimension(2), intent(in)  ithermal,
real*8, dimension(3,9), intent(in)  xl,
integer, dimension(2,*), intent(in)  nelemload,
integer, intent(in)  nelemadvec,
integer, intent(in)  nload,
character*8, dimension(*), intent(in)  lakon,
real*8, dimension(2,*), intent(inout)  xload,
integer, intent(in)  istep,
real*8, intent(in)  time,
real*8, intent(in)  ttime,
real*8, intent(in)  dtime,
character*20, dimension(*), intent(in)  sideload,
real*8, dimension(0:mi(2),*), intent(in)  vold,
integer, dimension(*), intent(in)  mi,
real*8, dimension(2,*), intent(in)  xloadold,
real*8, intent(in)  reltime,
integer, intent(in)  nmethod,
real*8, dimension(9), intent(inout)  tnl,
integer, intent(in)  iinc,
integer, dimension(*), intent(in)  iponoel,
integer, dimension(2,*), intent(in)  inoel,
integer, dimension(*), intent(in)  ielprop,
real*8, dimension(*), intent(in)  prop,
integer, dimension(*), intent(in)  ielmat,
real*8, dimension(0:3,ntmat_,*), intent(in)  shcon,
integer, dimension(*), intent(in)  nshcon,
real*8, dimension(0:1,ntmat_,*), intent(in)  rhcon,
integer, dimension(*), intent(in)  nrhcon,
integer, intent(in)  ntmat_,
integer, dimension(*), intent(in)  ipkon,
integer, dimension(*), intent(in)  kon,
real*8, dimension(0:6,ntmat_,*), intent(in)  cocon,
integer, dimension(2,*), intent(in)  ncocon,
integer, dimension(2,*)  ipobody,
integer, dimension(3,*)  ibody,
real*8, dimension(7,*)  xbody 
)
24 !
25 ! calculates the stiffness of an advective element.
26 ! An advective element consists of a face with a forced convection
27 ! film condition and a network node
28 !
29  implicit none
30 !
31  character*8 lakonl,lakon(*)
32  character*20 sideload(*),sideloadl
33 !
34  integer nope,i,ithermal(2),j,nelemload(2,*),nelemadvec,nload,id,
35  & nelem,ig,mint2d,iflag,istep,jltyp,nfield,mi(*),nmethod,k,iinc,
36  & node,nopes,iponoel(*),inoel(2,*),ielprop(*),ielmat(*),ipkon(*),
37  & nshcon(*),nrhcon(*),ntmat_,kon(*),ncocon(2,*),ipobody(2,*),
38  & ibody(3,*)
39 !
40  real*8 tl2(9),voldl(0:mi(2),20),xl(3,9),sinktemp,xi,et,weight,
41  & xl2(3,8),xsj2(3),shp2(7,9),coords(3),xs2(3,7),dxsj2,areaj,
42  & temp,xload(2,*),timeend(2),time,ttime,dtime,field,reltime,
43  & vold(0:mi(2),*),xloadold(2,*),tnl(9),tnlref,prop(*),
44  & shcon(0:3,ntmat_,*),rhcon(0:1,ntmat_,*),cocon(0:6,ntmat_,*),
45  & xbody(7,*),heatnod,heatfac
46 !
47  intent(in) nope,voldl,ithermal,xl,nelemload,nelemadvec,
48  & nload,lakon,istep,time,ttime,dtime,sideload,vold,mi,
49  & xloadold,reltime,nmethod,iinc,iponoel,inoel,ielprop,prop,
50  & ielmat,shcon,nshcon,rhcon,nrhcon,ntmat_,ipkon,kon,cocon,ncocon
51 !
52  intent(inout) xload,tnl
53 !
54  include "gauss.f"
55 !
56  data iflag /2/
57 !
58  timeend(1)=time
59  timeend(2)=ttime+time
60 !
61  heatfac=0.d0
62 !
63 ! number of nodes in the advective face
64 !
65  nopes=nope-1
66 !
67 ! temperature and displacements in the element's nodes
68 !
69  do i=1,nope
70  tl2(i)=voldl(0,i)
71  enddo
72  do i=1,nopes
73  if(ithermal(2).eq.2) then
74  do j=1,3
75  xl2(j,i)=xl(j,i)
76  enddo
77  else
78  do j=1,3
79  xl2(j,i)=xl(j,i)+voldl(j,i)
80  enddo
81  endif
82  enddo
83 !
84  call nident2(nelemload,nelemadvec,nload,id)
85 !
86 ! the second entry in nelemload points to the original
87 ! film loading
88 !
89  id=nelemload(2,id)
90 !
91 ! number of the original element
92 !
93  nelem=nelemload(1,id)
94  lakonl=lakon(nelem)
95 !
96 ! the next line originally ran:
97 ! read(sideload(id)(2:2),'(i1)') ig
98 ! it was replaced since the read statement might
99 ! cause problems in the multithreading version
100 !
101  ig=ichar(sideload(id)(2:2))-48
102 !
103 ! number of integration points
104 !
105  if(lakonl(4:5).eq.'8R') then
106  mint2d=1
107  elseif((lakonl(4:4).eq.'8').or.(lakonl(4:6).eq.'20R')) then
108  if((lakonl(7:7).eq.'A').or.(lakonl(7:7).eq.'E')) then
109  mint2d=2
110  else
111  mint2d=4
112  endif
113  elseif(lakonl(4:4).eq.'2') then
114  mint2d=9
115  elseif(lakonl(4:5).eq.'10') then
116  mint2d=3
117  elseif(lakonl(4:4).eq.'4') then
118  mint2d=1
119  elseif(lakonl(4:5).eq.'15') then
120  if(ig.le.2) then
121  mint2d=3
122  else
123  mint2d=4
124  endif
125  elseif(lakonl(4:4).eq.'6') then
126  mint2d=1
127  endif
128 !
129  do i=1,mint2d
130 !
131 ! copying the sink temperature to ensure the same
132 ! value in each integration point (sinktemp can be
133 ! changed in subroutine film: requirement from the
134 ! thermal people)
135 !
136  sinktemp=tl2(nope)
137 !
138  if((lakonl(4:5).eq.'8R').or.
139  & ((lakonl(4:4).eq.'6').and.(nopes.eq.4))) then
140  xi=gauss2d1(1,i)
141  et=gauss2d1(2,i)
142  weight=weight2d1(i)
143  elseif((lakonl(4:4).eq.'8').or.
144  & (lakonl(4:6).eq.'20R').or.
145  & ((lakonl(4:5).eq.'15').and.(nopes.eq.8))) then
146  xi=gauss2d2(1,i)
147  et=gauss2d2(2,i)
148  weight=weight2d2(i)
149  elseif(lakonl(4:4).eq.'2') then
150  xi=gauss2d3(1,i)
151  et=gauss2d3(2,i)
152  weight=weight2d3(i)
153  elseif((lakonl(4:5).eq.'10').or.
154  & ((lakonl(4:5).eq.'15').and.(nopes.eq.6))) then
155  xi=gauss2d5(1,i)
156  et=gauss2d5(2,i)
157  weight=weight2d5(i)
158  elseif((lakonl(4:4).eq.'4').or.
159  & ((lakonl(4:4).eq.'6').and.(nopes.eq.3))) then
160  xi=gauss2d4(1,i)
161  et=gauss2d4(2,i)
162  weight=weight2d4(i)
163  endif
164 !
165  if(nopes.eq.8) then
166  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
167  elseif(nopes.eq.4) then
168  call shape4q(xi,et,xl2,xsj2,xs2,shp2,iflag)
169  elseif(nopes.eq.6) then
170  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
171  else
172  call shape3tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
173  endif
174 !
175  dxsj2=dsqrt(xsj2(1)*xsj2(1)+xsj2(2)*xsj2(2)+
176  & xsj2(3)*xsj2(3))
177  areaj=dxsj2*weight
178 !
179  temp=0.d0
180  do j=1,nopes
181  temp=temp+tl2(j)*shp2(4,j)
182  enddo
183 !
184 ! for nonuniform load: determine the coordinates of the
185 ! point (transferred into the user subroutine)
186 !
187  if((sideload(id)(3:4).eq.'NU').or.
188  & (sideload(id)(5:6).eq.'NU')) then
189  do k=1,3
190  coords(k)=0.d0
191  do j=1,nopes
192  coords(k)=coords(k)+xl2(k,j)*shp2(4,j)
193  enddo
194  enddo
195 c read(sideload(id)(2:2),'(i1)') jltyp
196 c jltyp=jltyp+10
197  jltyp=ichar(sideload(id)(2:2))-38
198  node=nelemload(2,id)
199  sideloadl=sideload(id)
200  sideloadl(1:1)='F'
201  call film(xload(1,id),sinktemp,temp,istep,
202  & iinc,timeend,nelem,i,coords,jltyp,field,nfield,
203  & sideloadl,node,areaj,vold,mi,
204  & ipkon,kon,lakon,iponoel,inoel,ielprop,prop,ielmat,
205  & shcon,nshcon,rhcon,nrhcon,ntmat_,cocon,ncocon,
206  & ipobody,xbody,ibody,heatnod,heatfac)
207 c if(nmethod.eq.1) xload(1,id)=xloadold(1,id)+
208 c & (xload(1,id)-xloadold(1,id))*reltime
209  endif
210 !
211  tnlref=(xload(1,id)*(sinktemp-temp)+heatfac)*areaj
212  do j=1,nopes
213  tnl(j)=tnl(j)-tnlref*shp2(4,j)
214  enddo
215  tnl(nope)=tnl(nope)+tnlref
216  enddo
217 !
218 ! for some axisymmetric and plane strain elements only half the
219 ! number of integration points was used
220 !
221  if(((lakonl(4:4).eq.'8').or.(lakonl(4:6).eq.'20R')).and.
222  & ((lakonl(7:7).eq.'A').or.(lakonl(7:7).eq.'E'))) then
223  do i=1,nope
224  tnl(i)=2.d0*tnl(i)
225  enddo
226  endif
227 !
228  return
subroutine shape8q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape8q.f:20
subroutine shape3tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape3tri.f:20
subroutine nident2(x, px, n, id)
Definition: nident2.f:27
subroutine shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
subroutine film(h, sink, temp, kstep, kinc, time, noel, npt, coords, jltyp, field, nfield, loadtype, node, area, vold, mi, ipkon, kon, lakon, iponoel, inoel, ielprop, prop, ielmat, shcon, nshcon, rhcon, nrhcon, ntmat_, cocon, ncocon, ipobody, xbody, ibody, heatnod, heatfac)
Definition: film.f:24
Hosted by OpenAircraft.com, (Michigan UAV, LLC)