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

Go to the source code of this file.

Functions/Subroutines

subroutine mafillt (nef, ipnei, neifa, neiel, vfa, xxn, area, au, ad, jq, irow, nzs, b, vel, umel, xlet, xle, gradtfa, xxi, body, volume, ielfa, lakonf, ifabou, nbody, neq, dtimef, velo, veloo, cvfa, hcfa, cvel, gradvel, xload, gamma, xrlfa, xxj, nactdohinv, a1, a2, a3, flux, nefa, nefb, iau6, xxni, xxnj, iturbulent)
 

Function/Subroutine Documentation

◆ mafillt()

subroutine mafillt ( integer, intent(in)  nef,
integer, dimension(*), intent(in)  ipnei,
integer, dimension(*), intent(in)  neifa,
integer, dimension(*), intent(in)  neiel,
real*8, dimension(0:7,*), intent(in)  vfa,
real*8, dimension(3,*), intent(in)  xxn,
real*8, dimension(*), intent(in)  area,
real*8, dimension(*), intent(inout)  au,
real*8, dimension(*), intent(inout)  ad,
integer, dimension(*), intent(in)  jq,
integer, dimension(*), intent(in)  irow,
integer, intent(in)  nzs,
real*8, dimension(neq), intent(inout)  b,
real*8, dimension(nef,0:7), intent(in)  vel,
real*8, dimension(*), intent(in)  umel,
real*8, dimension(*), intent(in)  xlet,
real*8, dimension(*), intent(in)  xle,
real*8, dimension(3,*), intent(in)  gradtfa,
real*8, dimension(3,*), intent(in)  xxi,
real*8, dimension(0:3,*), intent(in)  body,
real*8, dimension(*), intent(in)  volume,
integer, dimension(4,*), intent(in)  ielfa,
character*8, dimension(*), intent(in)  lakonf,
integer, dimension(*), intent(in)  ifabou,
integer, intent(in)  nbody,
integer, intent(in)  neq,
real*8, intent(in)  dtimef,
real*8, dimension(nef,0:7), intent(in)  velo,
real*8, dimension(nef,0:7), intent(in)  veloo,
real*8, dimension(*), intent(in)  cvfa,
real*8, dimension(*), intent(in)  hcfa,
real*8, dimension(*), intent(in)  cvel,
real*8, dimension(3,3,*), intent(in)  gradvel,
real*8, dimension(2,*), intent(in)  xload,
real*8, dimension(*), intent(in)  gamma,
real*8, dimension(3,*), intent(in)  xrlfa,
real*8, dimension(3,*), intent(in)  xxj,
integer, dimension(*), intent(in)  nactdohinv,
real*8, intent(in)  a1,
real*8, intent(in)  a2,
real*8, intent(in)  a3,
real*8, dimension(*), intent(in)  flux,
integer, intent(in)  nefa,
integer, intent(in)  nefb,
integer, dimension(6,*)  iau6,
real*8, dimension(3,*)  xxni,
real*8, dimension(3,*)  xxnj,
integer  iturbulent 
)
25 !
26 ! filling the matrix for the conservation of energy
27 !
28  implicit none
29 !
30  logical knownflux
31 !
32  character*8 lakonf(*)
33 !
34  integer i,nef,indexf,ipnei(*),j,ifa,iel,neifa(*),
35  & neiel(*),jq(*),irow(*),nzs,ielfa(4,*),nefa,nefb,
36  & ipointer,ifabou(*),nbody,neq,indexb,nactdohinv(*),
37  & iau6(6,*),iturbulent
38 !
39  real*8 xflux,vfa(0:7,*),xxn(3,*),area(*),au(*),ad(*),b(neq),
40  & vel(nef,0:7),umel(*),xlet(*),xle(*),coef,gradtfa(3,*),
41  & xxi(3,*),body(0:3,*),volume(*),dtimef,velo(nef,0:7),
42  & veloo(nef,0:7),rhovol,constant,cvel(*),gradvel(3,3,*),
43  & cvfa(*),hcfa(*),div,xload(2,*),gamma(*),xrlfa(3,*),
44  & xxj(3,*),a1,a2,a3,flux(*),xxnj(3,*),xxni(3,*),visdis
45 !
46  intent(in) nef,ipnei,neifa,neiel,vfa,xxn,area,
47  & jq,irow,nzs,vel,umel,xlet,xle,gradtfa,xxi,
48  & body,volume,ielfa,lakonf,ifabou,nbody,neq,
49  & dtimef,velo,veloo,cvfa,hcfa,cvel,gradvel,xload,gamma,xrlfa,
50  & xxj,nactdohinv,a1,a2,a3,flux,nefa,nefb
51 !
52  intent(inout) au,ad,b
53 !
54  do i=nefa,nefb
55  do indexf=ipnei(i)+1,ipnei(i+1)
56 !
57 ! convection
58 !
59  ifa=neifa(indexf)
60  iel=neiel(indexf)
61  xflux=flux(indexf)*cvfa(ifa)
62 !
63  if(xflux.ge.0.d0) then
64 !
65 ! outflowing flux
66 !
67  ad(i)=ad(i)+xflux
68 !
69  b(i)=b(i)-gamma(ifa)*(vfa(0,ifa)-vel(i,0))*xflux
70 !
71  else
72  if(iel.gt.0) then
73 !
74 ! incoming flux from neighboring element
75 !
76  au(indexf)=au(indexf)+xflux
77 !
78  b(i)=b(i)-gamma(ifa)*(vfa(0,ifa)-vel(iel,0))*xflux
79 !
80  else
81 !
82 ! incoming flux through boundary
83 !
84  if(ielfa(2,ifa).lt.0) then
85  indexb=-ielfa(2,ifa)
86  if((ifabou(indexb).ne.0).or.
87  & (dabs(xflux).lt.1.d-10)) then
88  b(i)=b(i)-vfa(0,ifa)*xflux
89  else
90 c write(*,*) '*ERROR in mafillt: the tempera-'
91 c write(*,*) ' ture of an incoming flux'
92 c write(*,*) ' through face ',
93 c & indexf-ipnei(i),'of'
94 c write(*,*)' element ',nactdohinv(i),
95 c & ' is not given'
96  endif
97  else
98 c write(*,*) '*ERROR in mafillt: the tempera-'
99 c write(*,*) ' ture of an incoming flux'
100 c write(*,*) ' through face ',
101 c & indexf-ipnei(i),'of'
102 c write(*,*)' element ',nactdohinv(i),
103 c & ' is not given'
104  endif
105  endif
106  endif
107 !
108 ! diffusion
109 !
110  if(iel.ne.0) then
111 !
112 ! neighboring element
113 !
114  coef=hcfa(ifa)*area(ifa)/xlet(indexf)
115  ad(i)=ad(i)+coef
116  au(indexf)=au(indexf)-coef
117 !
118 ! correction for non-orthogonal grid
119 !
120  b(i)=b(i)+hcfa(ifa)*area(ifa)*
121  & (gradtfa(1,ifa)*xxnj(1,indexf)+
122  & gradtfa(2,ifa)*xxnj(2,indexf)+
123  & gradtfa(3,ifa)*xxnj(3,indexf))
124  else
125 !
126 ! boundary; either temperature given or adiabatic
127 ! or outlet
128 !
129  knownflux=.false.
130  ipointer=abs(ielfa(2,ifa))
131  if(ipointer.gt.0) then
132  if(ifabou(ipointer+6).gt.0) then
133 !
134 ! heat flux is known
135 !
136  b(i)=b(i)-xload(1,ifabou(ipointer+6))
137  knownflux=.true.
138 !
139  elseif((ifabou(ipointer).ne.0).or.
140  & (ifabou(ipointer+1).ne.0).or.
141  & (ifabou(ipointer+2).ne.0).or.
142  & (ifabou(ipointer+3).ne.0)) then
143 !
144 ! temperature given or no outlet:
145 ! temperature is assumed fixed
146 !
147  coef=hcfa(ifa)*area(ifa)/xle(indexf)
148  ad(i)=ad(i)+coef
149  b(i)=b(i)+coef*vfa(0,ifa)
150  else
151 !
152 ! outlet: no diffusion
153 !
154  endif
155  endif
156 !
157 ! correction for non-orthogonal grid
158 !
159  if(.not.knownflux) then
160  b(i)=b(i)+hcfa(ifa)*area(ifa)*
161  & (gradtfa(1,ifa)*xxni(1,indexf)+
162  & gradtfa(2,ifa)*xxni(2,indexf)+
163  & gradtfa(3,ifa)*xxni(3,indexf))
164  endif
165  endif
166  enddo
167 !
168 ! viscous dissipation
169 !
170  if(iturbulent.eq.0) then
171  visdis=umel(i)
172  else
173  visdis=umel(i)+vel(i,5)*vel(i,6)/vel(i,7)
174  endif
175 !
176  b(i)=b(i)+visdis*volume(i)*
177  & (2.d0*(gradvel(1,1,i)**2+gradvel(2,2,i)**2+
178  & gradvel(3,3,i)**2)+
179  & (gradvel(1,2,i)+gradvel(2,1,i))**2+
180  & (gradvel(1,3,i)+gradvel(3,1,i))**2+
181  & (gradvel(2,3,i)+gradvel(3,2,i))**2)
182 !
183 ! body heat source and body sources
184 !
185  rhovol=vel(i,5)*volume(i)
186 !
187  if(nbody.gt.0) then
188  b(i)=b(i)+rhovol*body(0,i)
189  endif
190 !
191 ! transient term
192 !
193  constant=rhovol*cvel(i)
194  b(i)=b(i)-(a2*velo(i,0)+a3*veloo(i,0))*constant
195  constant=a1*constant
196  ad(i)=ad(i)+constant
197 !
198  enddo
199 !
200  return
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
Hosted by OpenAircraft.com, (Michigan UAV, LLC)