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

Go to the source code of this file.

Functions/Subroutines

subroutine advecstiff (nope, voldl, ithermal, xl, nelemload, nelemadvec, nload, lakon, xload, istep, time, ttime, dtime, sideload, vold, mi, xloadold, reltime, nmethod, s, iinc, iponoel, inoel, ielprop, prop, ielmat, shcon, nshcon, rhcon, nrhcon, ntmat_, ipkon, kon, cocon, ncocon, ipobody, xbody, ibody)
 

Function/Subroutine Documentation

◆ advecstiff()

subroutine advecstiff ( integer, intent(in)  nope,
real*8, dimension(0:mi(2),9), 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(60,60), intent(inout)  s,
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,*), intent(in)  ipobody,
real*8, dimension(7,*), intent(in)  xbody,
integer, dimension(3,*), intent(in)  ibody 
)
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),9),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,*),s(60,60),sref,sref2,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,
51  & ncocon,ipobody,ibody,xbody
52 !
53  intent(inout) s,xload
54 !
55  include "gauss.f"
56 !
57  iflag=2
58 !
59  timeend(1)=time
60  timeend(2)=ttime+time
61 !
62 ! number of nodes in the advective face
63 !
64  nopes=nope-1
65 !
66 ! temperature and displacements in the element's nodes
67 !
68  do i=1,nope
69  tl2(i)=voldl(0,i)
70  enddo
71  do i=1,nopes
72  if(ithermal(2).eq.2) then
73  do j=1,3
74  xl2(j,i)=xl(j,i)
75  enddo
76  else
77  do j=1,3
78  xl2(j,i)=xl(j,i)+voldl(j,i)
79  enddo
80  endif
81  enddo
82 !
83  call nident2(nelemload,nelemadvec,nload,id)
84 !
85 ! the second entry in nelemload points to the original
86 ! film loading
87 !
88  id=nelemload(2,id)
89 !
90 ! number of the original element
91 !
92  nelem=nelemload(1,id)
93  lakonl=lakon(nelem)
94  read(sideload(id)(2:2),'(i1)') ig
95 !
96 ! number of integration points
97 !
98  if(lakonl(4:5).eq.'8R') then
99  mint2d=1
100  elseif((lakonl(4:4).eq.'8').or.(lakonl(4:6).eq.'20R')) then
101  if((lakonl(7:7).eq.'A').or.(lakonl(7:7).eq.'E')) then
102  mint2d=2
103  else
104  mint2d=4
105  endif
106  elseif(lakonl(4:4).eq.'2') then
107  mint2d=9
108  elseif(lakonl(4:5).eq.'10') then
109  mint2d=3
110  elseif(lakonl(4:4).eq.'4') then
111  mint2d=1
112  elseif(lakonl(4:5).eq.'15') then
113  if(ig.le.2) then
114  mint2d=3
115  else
116  mint2d=4
117  endif
118  elseif(lakonl(4:4).eq.'6') then
119  mint2d=1
120  endif
121 !
122  do i=1,mint2d
123 !
124 ! copying the sink temperature to ensure the same
125 ! value in each integration point (sinktemp can be
126 ! changed in subroutine film: requirement from the
127 ! thermal people)
128 !
129  sinktemp=tl2(nope)
130 !
131  if((lakonl(4:5).eq.'8R').or.
132  & ((lakonl(4:4).eq.'6').and.(nopes.eq.4))) then
133  xi=gauss2d1(1,i)
134  et=gauss2d1(2,i)
135  weight=weight2d1(i)
136  elseif((lakonl(4:4).eq.'8').or.
137  & (lakonl(4:6).eq.'20R').or.
138  & ((lakonl(4:5).eq.'15').and.(nopes.eq.8))) then
139  xi=gauss2d2(1,i)
140  et=gauss2d2(2,i)
141  weight=weight2d2(i)
142  elseif(lakonl(4:4).eq.'2') then
143  xi=gauss2d3(1,i)
144  et=gauss2d3(2,i)
145  weight=weight2d3(i)
146  elseif((lakonl(4:5).eq.'10').or.
147  & ((lakonl(4:5).eq.'15').and.(nopes.eq.6))) then
148  xi=gauss2d5(1,i)
149  et=gauss2d5(2,i)
150  weight=weight2d5(i)
151  elseif((lakonl(4:4).eq.'4').or.
152  & ((lakonl(4:4).eq.'6').and.(nopes.eq.3))) then
153  xi=gauss2d4(1,i)
154  et=gauss2d4(2,i)
155  weight=weight2d4(i)
156  endif
157 !
158  if(nopes.eq.8) then
159  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
160  elseif(nopes.eq.4) then
161  call shape4q(xi,et,xl2,xsj2,xs2,shp2,iflag)
162  elseif(nopes.eq.6) then
163  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
164  else
165  call shape3tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
166  endif
167 !
168  dxsj2=dsqrt(xsj2(1)*xsj2(1)+xsj2(2)*xsj2(2)+
169  & xsj2(3)*xsj2(3))
170  areaj=dxsj2*weight
171 !
172  temp=0.d0
173  do j=1,nopes
174  temp=temp+tl2(j)*shp2(4,j)
175  enddo
176 !
177 ! for nonuniform load: determine the coordinates of the
178 ! point (transferred into the user subroutine)
179 !
180  if((sideload(id)(3:4).eq.'NU').or.
181  & (sideload(id)(5:6).eq.'NU')) then
182  do k=1,3
183  coords(k)=0.d0
184  do j=1,nopes
185  coords(k)=coords(k)+xl2(k,j)*shp2(4,j)
186  enddo
187  enddo
188  read(sideload(id)(2:2),'(i1)') jltyp
189  jltyp=jltyp+10
190  node=nelemload(2,id)
191  sideloadl=sideload(id)
192  sideloadl(1:1)='F'
193  call film(xload(1,id),sinktemp,temp,istep,
194  & iinc,timeend,nelem,i,coords,jltyp,field,nfield,
195  & sideloadl,node,areaj,vold,mi,
196  & ipkon,kon,lakon,iponoel,inoel,ielprop,prop,ielmat,
197  & shcon,nshcon,rhcon,nrhcon,ntmat_,cocon,ncocon,
198  & ipobody,xbody,ibody,heatnod,heatfac)
199 c if(nmethod.eq.1) xload(1,id)=xloadold(1,id)+
200 c & (xload(1,id)-xloadold(1,id))*reltime
201  endif
202 !
203  sref=xload(1,id)*areaj
204  shp2(4,nope)=-1.d0
205 !
206  do j=1,nope
207  sref2=sref*shp2(4,j)
208  s(nope,j)=s(nope,j)-sref2
209  do k=1,nopes
210  s(k,j)=s(k,j)+sref2*shp2(4,k)
211  enddo
212  enddo
213  enddo
214 !
215 ! for some axisymmetric and plane strain elements only half the
216 ! number of integration points was used
217 !
218  if(((lakonl(4:4).eq.'8').or.(lakonl(4:6).eq.'20R')).and.
219  & ((lakonl(7:7).eq.'A').or.(lakonl(7:7).eq.'E'))) then
220  do i=1,nope
221  do j=1,nope
222  s(i,j)=2.d0*s(i,j)
223  enddo
224  enddo
225  endif
226 !
227  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)