xref: /phasta/phSolver/incompressible/e3ivar.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine e3ivar (yl,          acl,       shpfun,
2*59599516SKenneth E. Jansen     &                   shgl,        xl,
3*59599516SKenneth E. Jansen     &                   aci,         g1yi,      g2yi,
4*59599516SKenneth E. Jansen     &                   g3yi,        shg,       dxidx,
5*59599516SKenneth E. Jansen     &                   WdetJ,       rho,       pres,
6*59599516SKenneth E. Jansen     &                   u1,          u2,        u3,
7*59599516SKenneth E. Jansen     &                   ql,          rLui,      src,
8*59599516SKenneth E. Jansen     &                   rerrl,       rlsl,      rlsli,
9*59599516SKenneth E. Jansen     &                   dwl)
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansenc----------------------------------------------------------------------
12*59599516SKenneth E. Jansenc
13*59599516SKenneth E. Jansenc  This routine computes the variables at integration point.
14*59599516SKenneth E. Jansenc
15*59599516SKenneth E. Jansenc input:
16*59599516SKenneth E. Jansenc  yl     (npro,nshl,ndof)      : primitive variables
17*59599516SKenneth E. Jansenc  acl    (npro,nshl,ndof)      : prim.var. accel.
18*59599516SKenneth E. Jansenc  shp    (nen)                 : element shape-functions
19*59599516SKenneth E. Jansenc  shgl   (nsd,nen)             : element local-grad-shape-functions
20*59599516SKenneth E. Jansenc  xl     (npro,nenl,nsd)       : nodal coordinates at current step
21*59599516SKenneth E. Jansenc  ql     (npro,nshl,nsd*nsd) : diffusive flux vector
22*59599516SKenneth E. Jansenc  rlsl   (npro,nshl,6)       : resolved Leonard stresses
23*59599516SKenneth E. Jansenc
24*59599516SKenneth E. Jansenc output:
25*59599516SKenneth E. Jansenc  aci    (npro,3)              : primvar accel. variables
26*59599516SKenneth E. Jansenc  g1yi   (npro,ndof)           : grad-y in direction 1
27*59599516SKenneth E. Jansenc  g2yi   (npro,ndof)           : grad-y in direction 2
28*59599516SKenneth E. Jansenc  g3yi   (npro,ndof)           : grad-y in direction 3
29*59599516SKenneth E. Jansenc  shg    (npro,nshl,nsd)       : element global grad-shape-functions
30*59599516SKenneth E. Jansenc  dxidx  (npro,nsd,nsd)        : inverse of deformation gradient
31*59599516SKenneth E. Jansenc  WdetJ  (npro)                : weighted Jacobian
32*59599516SKenneth E. Jansenc  rho    (npro)                : density
33*59599516SKenneth E. Jansenc  pres   (npro)                : pressure
34*59599516SKenneth E. Jansenc  u1     (npro)                : x1-velocity component
35*59599516SKenneth E. Jansenc  u2     (npro)                : x2-velocity component
36*59599516SKenneth E. Jansenc  u3     (npro)                : x3-velocity component
37*59599516SKenneth E. Jansenc  rLui   (npro,nsd)            : xi-momentum residual
38*59599516SKenneth E. Jansenc  src    (npro,nsd)            : body force term (not density weighted)
39*59599516SKenneth E. Jansenc  rlsli  (npro,6)              : resolved Leonard stresses at quad pt
40*59599516SKenneth E. Jansenc
41*59599516SKenneth E. Jansenc locally calculated and used
42*59599516SKenneth E. Jansenc  divqi  (npro,nsd+isurf)      : divergence of reconstructed quantity
43*59599516SKenneth E. Jansenc
44*59599516SKenneth E. Jansenc Zdenek Johan, Summer 1990. (Modified from e2ivar.f)
45*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991. (Fortran 90)
46*59599516SKenneth E. Jansenc Kenneth Jansen, Winter 1997. Primitive Variables
47*59599516SKenneth E. Jansenc Christian Whiting, Winter 1999. (uBar formulation)
48*59599516SKenneth E. Jansenc
49*59599516SKenneth E. Jansenc----------------------------------------------------------------------
50*59599516SKenneth E. Jansenc
51*59599516SKenneth E. Jansen      include "common.h"
52*59599516SKenneth E. Jansenc
53*59599516SKenneth E. Jansenc  passed arrays
54*59599516SKenneth E. Jansenc
55*59599516SKenneth E. Jansen      dimension yl(npro,nshl,ndof),        dwl(npro,nenl),
56*59599516SKenneth E. Jansen     &            acl(npro,nshl,ndof),       shpfun(npro,nshl),
57*59599516SKenneth E. Jansen     &            shgl(npro,nsd,nshl),       xl(npro,nenl,nsd),
58*59599516SKenneth E. Jansen     &            aci(npro,nsd),             g1yi(npro,ndof),
59*59599516SKenneth E. Jansen     &            g2yi(npro,ndof),           g3yi(npro,ndof),
60*59599516SKenneth E. Jansen     &            shg(npro,nshl,nsd),        dxidx(npro,nsd,nsd),
61*59599516SKenneth E. Jansen     &            WdetJ(npro),
62*59599516SKenneth E. Jansen     &            rho(npro),                 pres(npro),
63*59599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
64*59599516SKenneth E. Jansen     &            u3(npro),                  divqi(npro,nflow-1+isurf),
65*59599516SKenneth E. Jansen     &            ql(npro,nshl,idflx),       rLui(npro,nsd),
66*59599516SKenneth E. Jansen     &            src(npro,nsd), Temp(npro),xx(npro,nsd)
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansen        dimension tmp(npro), tmp1(npro), dkei(npro), dist2w(npro)
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansen        dimension rlsl(npro,nshl,6),         rlsli(npro,6)
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansen        real*8    rerrl(npro,nshl,6), omega(3), divu(npro)
73*59599516SKenneth E. Jansen        dimension gyti(npro,nsd),            gradh(npro,nsd),
74*59599516SKenneth E. Jansen     &            sforce(npro,3),            weber(npro),
75*59599516SKenneth E. Jansen     &            Sclr(npro)
76*59599516SKenneth E. Jansenc
77*59599516SKenneth E. Jansenc.... ------------->  Primitive variables at int. point  <--------------
78*59599516SKenneth E. Jansenc
79*59599516SKenneth E. Jansenc.... compute primitive variables
80*59599516SKenneth E. Jansenc
81*59599516SKenneth E. Jansen       pres = zero
82*59599516SKenneth E. Jansen       u1   = zero
83*59599516SKenneth E. Jansen       u2   = zero
84*59599516SKenneth E. Jansen       u3   = zero
85*59599516SKenneth E. Jansenc
86*59599516SKenneth E. Jansen       do n = 1, nshl
87*59599516SKenneth E. Jansen          pres = pres + shpfun(:,n) * yl(:,n,1)
88*59599516SKenneth E. Jansen          u1   = u1   + shpfun(:,n) * yl(:,n,2)
89*59599516SKenneth E. Jansen          u2   = u2   + shpfun(:,n) * yl(:,n,3)
90*59599516SKenneth E. Jansen          u3   = u3   + shpfun(:,n) * yl(:,n,4)
91*59599516SKenneth E. Jansen       enddo
92*59599516SKenneth E. Jansen       if(matflg(5,1).eq.2) then ! boussinesq body force
93*59599516SKenneth E. Jansen          Temp = zero
94*59599516SKenneth E. Jansen          do n = 1, nshl
95*59599516SKenneth E. Jansen             Temp = Temp + shpfun(:,n) * yl(:,n,5)
96*59599516SKenneth E. Jansen          enddo
97*59599516SKenneth E. Jansen       endif
98*59599516SKenneth E. Jansen       if(matflg(5,1).eq.3.or.matflg(6,1).eq.1) then
99*59599516SKenneth E. Jansenc         user-specified body force or coriolis force specified
100*59599516SKenneth E. Jansen                 xx = zero
101*59599516SKenneth E. Jansen          do n  = 1,nenl
102*59599516SKenneth E. Jansen             xx(:,1) = xx(:,1)  + shpfun(:,n) * xl(:,n,1)
103*59599516SKenneth E. Jansen             xx(:,2) = xx(:,2)  + shpfun(:,n) * xl(:,n,2)
104*59599516SKenneth E. Jansen             xx(:,3) = xx(:,3)  + shpfun(:,n) * xl(:,n,3)
105*59599516SKenneth E. Jansen          enddo
106*59599516SKenneth E. Jansen       endif
107*59599516SKenneth E. Jansenc
108*59599516SKenneth E. Jansen       if(iRANS.eq.-2) then ! kay-epsilon
109*59599516SKenneth E. Jansen          dist2w = zero
110*59599516SKenneth E. Jansen          do n = 1, nenl
111*59599516SKenneth E. Jansen             dist2w = dist2w + shpfun(:,n) * dwl(:,n)
112*59599516SKenneth E. Jansen          enddo
113*59599516SKenneth E. Jansen       endif
114*59599516SKenneth E. Jansenc
115*59599516SKenneth E. Jansen
116*59599516SKenneth E. Jansen       if( (iLES.gt.10).and.(iLES.lt.20))  then  ! bardina
117*59599516SKenneth E. Jansen       rlsli = zero
118*59599516SKenneth E. Jansen       do n = 1, nshl
119*59599516SKenneth E. Jansen
120*59599516SKenneth E. Jansen          rlsli(:,1) = rlsli(:,1) + shpfun(:,n) * rlsl(:,n,1)
121*59599516SKenneth E. Jansen          rlsli(:,2) = rlsli(:,2) + shpfun(:,n) * rlsl(:,n,2)
122*59599516SKenneth E. Jansen          rlsli(:,3) = rlsli(:,3) + shpfun(:,n) * rlsl(:,n,3)
123*59599516SKenneth E. Jansen          rlsli(:,4) = rlsli(:,4) + shpfun(:,n) * rlsl(:,n,4)
124*59599516SKenneth E. Jansen          rlsli(:,5) = rlsli(:,5) + shpfun(:,n) * rlsl(:,n,5)
125*59599516SKenneth E. Jansen          rlsli(:,6) = rlsli(:,6) + shpfun(:,n) * rlsl(:,n,6)
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen       enddo
128*59599516SKenneth E. Jansen       else
129*59599516SKenneth E. Jansen          rlsli = zero
130*59599516SKenneth E. Jansen       endif
131*59599516SKenneth E. Jansenc
132*59599516SKenneth E. Jansenc.... ----------------------->  accel. at int. point  <----------------------
133*59599516SKenneth E. Jansenc
134*59599516SKenneth E. Jansen       aci = zero
135*59599516SKenneth E. Jansen       do n = 1, nshl
136*59599516SKenneth E. Jansen          aci(:,1) = aci(:,1) + shpfun(:,n) * acl(:,n,2)
137*59599516SKenneth E. Jansen          aci(:,2) = aci(:,2) + shpfun(:,n) * acl(:,n,3)
138*59599516SKenneth E. Jansen          aci(:,3) = aci(:,3) + shpfun(:,n) * acl(:,n,4)
139*59599516SKenneth E. Jansen       enddo
140*59599516SKenneth E. Jansenc
141*59599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
142*59599516SKenneth E. Jansenc
143*59599516SKenneth E. Jansen       call e3metric( xl,         shgl,       dxidx,
144*59599516SKenneth E. Jansen     &                shg,        WdetJ)
145*59599516SKenneth E. Jansenc
146*59599516SKenneth E. Jansenc.... compute the global gradient of u and P
147*59599516SKenneth E. Jansenc
148*59599516SKenneth E. Jansenc
149*59599516SKenneth E. Jansen       g1yi = zero
150*59599516SKenneth E. Jansen       g2yi = zero
151*59599516SKenneth E. Jansen       g3yi = zero
152*59599516SKenneth E. Jansen       do n = 1, nshl
153*59599516SKenneth E. Jansen          g1yi(:,1) = g1yi(:,1) + shg(:,n,1) * yl(:,n,1)
154*59599516SKenneth E. Jansen          g1yi(:,2) = g1yi(:,2) + shg(:,n,1) * yl(:,n,2)
155*59599516SKenneth E. Jansen          g1yi(:,3) = g1yi(:,3) + shg(:,n,1) * yl(:,n,3)
156*59599516SKenneth E. Jansen          g1yi(:,4) = g1yi(:,4) + shg(:,n,1) * yl(:,n,4)
157*59599516SKenneth E. Jansenc
158*59599516SKenneth E. Jansen          g2yi(:,1) = g2yi(:,1) + shg(:,n,2) * yl(:,n,1)
159*59599516SKenneth E. Jansen          g2yi(:,2) = g2yi(:,2) + shg(:,n,2) * yl(:,n,2)
160*59599516SKenneth E. Jansen          g2yi(:,3) = g2yi(:,3) + shg(:,n,2) * yl(:,n,3)
161*59599516SKenneth E. Jansen          g2yi(:,4) = g2yi(:,4) + shg(:,n,2) * yl(:,n,4)
162*59599516SKenneth E. Jansenc
163*59599516SKenneth E. Jansen          g3yi(:,1) = g3yi(:,1) + shg(:,n,3) * yl(:,n,1)
164*59599516SKenneth E. Jansen          g3yi(:,2) = g3yi(:,2) + shg(:,n,3) * yl(:,n,2)
165*59599516SKenneth E. Jansen          g3yi(:,3) = g3yi(:,3) + shg(:,n,3) * yl(:,n,3)
166*59599516SKenneth E. Jansen          g3yi(:,4) = g3yi(:,4) + shg(:,n,3) * yl(:,n,4)
167*59599516SKenneth E. Jansen       enddo
168*59599516SKenneth E. Jansen
169*59599516SKenneth E. Jansen       divqi = zero
170*59599516SKenneth E. Jansen       idflow = 3
171*59599516SKenneth E. Jansen       if ( idiff >= 1 .or. isurf==1 ) then
172*59599516SKenneth E. Jansenc
173*59599516SKenneth E. Jansenc.... compute divergence of diffusive flux vector, qi,i
174*59599516SKenneth E. Jansenc
175*59599516SKenneth E. Jansen          if(idiff >= 1) then
176*59599516SKenneth E. Jansen             do n=1, nshl
177*59599516SKenneth E. Jansen                divqi(:,1) = divqi(:,1) + shg(:,n,1)*ql(:,n,1 )
178*59599516SKenneth E. Jansen     &                                  + shg(:,n,2)*ql(:,n,4 )
179*59599516SKenneth E. Jansen     &                                  + shg(:,n,3)*ql(:,n,7 )
180*59599516SKenneth E. Jansen
181*59599516SKenneth E. Jansen                divqi(:,2) = divqi(:,2) + shg(:,n,1)*ql(:,n,2 )
182*59599516SKenneth E. Jansen     &                                  + shg(:,n,2)*ql(:,n,5 )
183*59599516SKenneth E. Jansen     &                                  + shg(:,n,3)*ql(:,n,8)
184*59599516SKenneth E. Jansen
185*59599516SKenneth E. Jansen                divqi(:,3) = divqi(:,3) + shg(:,n,1)*ql(:,n,3 )
186*59599516SKenneth E. Jansen     &                                  + shg(:,n,2)*ql(:,n,6 )
187*59599516SKenneth E. Jansen     &                                  + shg(:,n,3)*ql(:,n,9 )
188*59599516SKenneth E. Jansen
189*59599516SKenneth E. Jansen          enddo
190*59599516SKenneth E. Jansen
191*59599516SKenneth E. Jansen          endif                 !end of idiff
192*59599516SKenneth E. Jansenc
193*59599516SKenneth E. Jansen          if (isurf .eq. 1) then
194*59599516SKenneth E. Jansenc     .... divergence of normal calculation (curvature)
195*59599516SKenneth E. Jansen             do n=1, nshl
196*59599516SKenneth E. Jansen                divqi(:,idflow+1) = divqi(:,idflow+1)
197*59599516SKenneth E. Jansen     &               + shg(:,n,1)*ql(:,n,idflx-2)
198*59599516SKenneth E. Jansen     &               + shg(:,n,2)*ql(:,n,idflx-1)
199*59599516SKenneth E. Jansen     &               + shg(:,n,3)*ql(:,n,idflx)
200*59599516SKenneth E. Jansen             enddo
201*59599516SKenneth E. Jansenc     .... initialization of some variables
202*59599516SKenneth E. Jansen             Sclr = zero
203*59599516SKenneth E. Jansen             gradh= zero
204*59599516SKenneth E. Jansen             gyti = zero
205*59599516SKenneth E. Jansen             sforce=zero
206*59599516SKenneth E. Jansen             do i = 1, npro
207*59599516SKenneth E. Jansen                do n = 1, nshl
208*59599516SKenneth E. Jansen                   Sclr(i) = Sclr(i) + shpfun(i,n) * yl(i,n,6) !scalar
209*59599516SKenneth E. Jansenc
210*59599516SKenneth E. Jansenc     .... compute the global gradient of Scalar variable
211*59599516SKenneth E. Jansenc
212*59599516SKenneth E. Jansen                   gyti(i,1) = gyti(i,1) + shg(i,n,1) * yl(i,n,6)
213*59599516SKenneth E. Jansen                   gyti(i,2) = gyti(i,2) + shg(i,n,2) * yl(i,n,6)
214*59599516SKenneth E. Jansen                   gyti(i,3) = gyti(i,3) + shg(i,n,3) * yl(i,n,6)
215*59599516SKenneth E. Jansenc
216*59599516SKenneth E. Jansen                enddo
217*59599516SKenneth E. Jansen
218*59599516SKenneth E. Jansen                if (abs (sclr(i)) .le. epsilon_ls) then
219*59599516SKenneth E. Jansen                   gradh(i,1) = 0.5/epsilon_ls * (1.0
220*59599516SKenneth E. Jansen     &                  + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,1)
221*59599516SKenneth E. Jansen                   gradh(i,2) = 0.5/epsilon_ls * (1.0
222*59599516SKenneth E. Jansen     &                  + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,2)
223*59599516SKenneth E. Jansen                   gradh(i,3) = 0.5/epsilon_ls * (1.0
224*59599516SKenneth E. Jansen     &                  + cos(pi*Sclr(i)/epsilon_ls)) * gyti(i,3)
225*59599516SKenneth E. Jansen                endif
226*59599516SKenneth E. Jansen             enddo              !end of the loop over npro
227*59599516SKenneth E. Jansenc
228*59599516SKenneth E. Jansenc .. surface tension force calculation
229*59599516SKenneth E. Jansenc .. divide by density now as it gets multiplied in e3res.f, as surface
230*59599516SKenneth E. Jansenc    tension force is already in the form of force per unti volume
231*59599516SKenneth E. Jansenc
232*59599516SKenneth E. Jansen             weber(:) = Bo
233*59599516SKenneth E. Jansen             sforce(:,1) = -(1.0/weber(:)) * divqi(:,idflow+1) !x-direction
234*59599516SKenneth E. Jansen     &            *gradh(:,1) /rho(:)
235*59599516SKenneth E. Jansen             sforce(:,2) = -(1.0/weber(:)) * divqi(:,idflow+1) !y-direction
236*59599516SKenneth E. Jansen     &            *gradh(:,2) /rho(:)
237*59599516SKenneth E. Jansen             sforce(:,3) = -(1.0/weber(:)) * divqi(:,idflow+1) !z-direction
238*59599516SKenneth E. Jansen     &            *gradh(:,3) /rho(:)
239*59599516SKenneth E. Jansenc
240*59599516SKenneth E. Jansen          endif        ! end of the surface tension force calculation
241*59599516SKenneth E. Jansen       endif           ! diffusive flux computation
242*59599516SKenneth E. Jansenc
243*59599516SKenneth E. Jansenc Calculate strong form of pde as well as the source term
244*59599516SKenneth E. Jansenc
245*59599516SKenneth E. Jansen       call e3resStrongPDE(
246*59599516SKenneth E. Jansen     &      aci,  u1,   u2,   u3,   Temp, rho,  xx,
247*59599516SKenneth E. Jansen     &            g1yi, g2yi, g3yi,
248*59599516SKenneth E. Jansen     &      rLui, src, divqi)
249*59599516SKenneth E. Jansenc
250*59599516SKenneth E. Jansenc.... take care of the surface tension force term here
251*59599516SKenneth E. Jansenc
252*59599516SKenneth E. Jansen       if (isurf .eq. 1) then  ! note multiplied by density in e3res.f
253*59599516SKenneth E. Jansen          src(:,1) = src(:,1) + sforce(:,1)
254*59599516SKenneth E. Jansen          src(:,2) = src(:,2) + sforce(:,2)
255*59599516SKenneth E. Jansen          src(:,3) = src(:,3) + sforce(:,3)
256*59599516SKenneth E. Jansen       endif
257*59599516SKenneth E. Jansenc
258*59599516SKenneth E. Jansenc.... -------------------> error calculation  <-----------------
259*59599516SKenneth E. Jansenc
260*59599516SKenneth E. Jansen       if((ierrcalc.eq.1).and.(nitr.eq.iter)) then
261*59599516SKenneth E. Jansen          do ia=1,nshl
262*59599516SKenneth E. Jansen             tmp=shpfun(:,ia)*WdetJ(:)
263*59599516SKenneth E. Jansen             tmp1=shpfun(:,ia)*Qwt(lcsyst,intp)
264*59599516SKenneth E. Jansen             rerrl(:,ia,1) = rerrl(:,ia,1) +
265*59599516SKenneth E. Jansen     &                       tmp1(:)*rLui(:,1)*rLui(:,1)
266*59599516SKenneth E. Jansen             rerrl(:,ia,2) = rerrl(:,ia,2) +
267*59599516SKenneth E. Jansen     &                       tmp1(:)*rLui(:,2)*rLui(:,2)
268*59599516SKenneth E. Jansen             rerrl(:,ia,3) = rerrl(:,ia,3) +
269*59599516SKenneth E. Jansen     &                       tmp1(:)*rLui(:,3)*rLui(:,3)
270*59599516SKenneth E. Jansen
271*59599516SKenneth E. Jansen             rerrl(:,ia,4) = rerrl(:,ia,4) +
272*59599516SKenneth E. Jansen     &                       tmp(:)*divqi(:,1)*divqi(:,1)
273*59599516SKenneth E. Jansen             rerrl(:,ia,5) = rerrl(:,ia,5) +
274*59599516SKenneth E. Jansen     &                       tmp(:)*divqi(:,2)*divqi(:,2)
275*59599516SKenneth E. Jansen             rerrl(:,ia,6) = rerrl(:,ia,6) +
276*59599516SKenneth E. Jansen     &                       tmp(:)*divqi(:,3)*divqi(:,3)
277*59599516SKenneth E. Jansen          enddo
278*59599516SKenneth E. Jansen       endif
279*59599516SKenneth E. Jansen       distcalc=0  ! return to 1 if you want to compute T-S instability
280*59599516SKenneth E. Jansen       if(distcalc.eq.1) then
281*59599516SKenneth E. Jansenc
282*59599516SKenneth E. Jansenc.... ----------------------->  dist. kin energy at int. point  <--------------
283*59599516SKenneth E. Jansenc
284*59599516SKenneth E. Jansen
285*59599516SKenneth E. Jansen       if (ires .ne. 2 .and. iter.eq.1)  then  !only do at beginning of step
286*59599516SKenneth E. Jansenc
287*59599516SKenneth E. Jansenc calc exact velocity for a channel at quadrature points.
288*59599516SKenneth E. Jansenc
289*59599516SKenneth E. Jansen       dkei=0.0
290*59599516SKenneth E. Jansenc
291*59599516SKenneth E. Jansen       do n = 1, nenl
292*59599516SKenneth E. Jansen          dkei = dkei + shpfun(:,n) * (1.0-xl(:,n,2)**2) !u_ex^~ (in FEM space)
293*59599516SKenneth E. Jansen       enddo
294*59599516SKenneth E. Jansen          dkei = (u1(:)-dkei)**2 +u2(:)**2  ! u'^2+v'^2
295*59599516SKenneth E. Jansen          dkei = dkei*WdetJ  ! mult function*W*det of jacobian to
296*59599516SKenneth E. Jansenc                              get this quadrature point contribution
297*59599516SKenneth E. Jansen          dke  = dke+sum(dkei) ! we move the sum over elements inside of the
298*59599516SKenneth E. Jansenc                              sum over quadrature to save memory (we want
299*59599516SKenneth E. Jansenc                              a scalar only)
300*59599516SKenneth E. Jansen       endif
301*59599516SKenneth E. Jansen       endif
302*59599516SKenneth E. Jansenc
303*59599516SKenneth E. Jansenc.... return
304*59599516SKenneth E. Jansenc
305*59599516SKenneth E. Jansen       return
306*59599516SKenneth E. Jansen       end
307*59599516SKenneth E. Jansen
308*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
309*59599516SKenneth E. Jansenc
310*59599516SKenneth E. Jansenc     Calculate the variables for the scalar advection-diffusion
311*59599516SKenneth E. Jansenc     equation.
312*59599516SKenneth E. Jansenc
313*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
314*59599516SKenneth E. Jansen      subroutine e3ivarSclr (yl,          acl,       shpfun,
315*59599516SKenneth E. Jansen     &                      shgl,        xl,        xmudmi,
316*59599516SKenneth E. Jansen     &                      Sclr,        Sdot,      gradS,
317*59599516SKenneth E. Jansen     &                      shg,         dxidx,     WdetJ,
318*59599516SKenneth E. Jansen     &                      u1,          u2,        u3,
319*59599516SKenneth E. Jansen     &                      ql,          rLS ,       SrcR,
320*59599516SKenneth E. Jansen     &                      SrcL,        uMod,      dwl,
321*59599516SKenneth E. Jansen     &                      diffus,      srcRat)
322*59599516SKenneth E. Jansenc
323*59599516SKenneth E. Jansen      include "common.h"
324*59599516SKenneth E. Jansenc
325*59599516SKenneth E. Jansenc  passed arrays
326*59599516SKenneth E. Jansenc
327*59599516SKenneth E. Jansen      dimension yl(npro,nshl,ndof),        acl(npro,nshl,ndof),
328*59599516SKenneth E. Jansen     &          Sclr(npro),                Sdot(npro),
329*59599516SKenneth E. Jansen     &          gradS(npro,nsd),           shpfun(npro,nshl),
330*59599516SKenneth E. Jansen     &          shgl(npro,nsd,nshl),       xl(npro,nenl,nsd),
331*59599516SKenneth E. Jansen     &          shg(npro,nshl,nsd),        dxidx(npro,nsd,nsd),
332*59599516SKenneth E. Jansen     &          WdetJ(npro),
333*59599516SKenneth E. Jansen     &          u1(npro),                  u2(npro),
334*59599516SKenneth E. Jansen     &          u3(npro),                  divS(npro),
335*59599516SKenneth E. Jansen     &          ql(npro,nshl,nsd),         rLS(npro),
336*59599516SKenneth E. Jansen     &          SrcR(npro),                 SrcL(npro),
337*59599516SKenneth E. Jansen     &          dwl(npro,nshl),            diffus(npro),
338*59599516SKenneth E. Jansen     &          umod(npro,nsd), Temp(npro),xx(npro,nsd),
339*59599516SKenneth E. Jansen     &          divqi(npro)
340*59599516SKenneth E. Jansenc
341*59599516SKenneth E. Jansen      dimension tmp(npro), srcRat(npro)
342*59599516SKenneth E. Jansen      real*8 rLui(npro,nsd),     aci(npro,nsd),
343*59599516SKenneth E. Jansen     &       g1yi(npro,nflow),   g2yi(npro,nflow),
344*59599516SKenneth E. Jansen     &       g3yi(npro,nflow),
345*59599516SKenneth E. Jansen     &       src(npro,nsd),      rho(npro),
346*59599516SKenneth E. Jansen     &       rmu(npro)
347*59599516SKenneth E. Jansen      real*8 uBar(npro,nsd), xmudmi(npro,ngauss)
348*59599516SKenneth E. Jansen
349*59599516SKenneth E. Jansenc
350*59599516SKenneth E. Jansenc.... ------------->  Primitive variables at int. point  <--------------
351*59599516SKenneth E. Jansenc
352*59599516SKenneth E. Jansenc.... compute primitive variables
353*59599516SKenneth E. Jansenc
354*59599516SKenneth E. Jansen      u1   = zero
355*59599516SKenneth E. Jansen      u2   = zero
356*59599516SKenneth E. Jansen      u3   = zero
357*59599516SKenneth E. Jansen      Sclr = zero
358*59599516SKenneth E. Jansenc
359*59599516SKenneth E. Jansen      id=isclr+5
360*59599516SKenneth E. Jansen      do n = 1, nshl
361*59599516SKenneth E. Jansen         u1   = u1   + shpfun(:,n) * yl(:,n,2)
362*59599516SKenneth E. Jansen         u2   = u2   + shpfun(:,n) * yl(:,n,3)
363*59599516SKenneth E. Jansen         u3   = u3   + shpfun(:,n) * yl(:,n,4)
364*59599516SKenneth E. Jansen         Sclr = Sclr + shpfun(:,n) * yl(:,n,id)
365*59599516SKenneth E. Jansen      enddo
366*59599516SKenneth E. Jansenc
367*59599516SKenneth E. Jansenc
368*59599516SKenneth E. Jansenc.... ----------------------->  dS/dt at int. point  <----------------------
369*59599516SKenneth E. Jansenc
370*59599516SKenneth E. Jansen      Sdot = zero
371*59599516SKenneth E. Jansen      do n = 1, nshl
372*59599516SKenneth E. Jansen         Sdot = Sdot + shpfun(:,n) * acl(:,n,id)
373*59599516SKenneth E. Jansen      enddo
374*59599516SKenneth E. Jansenc
375*59599516SKenneth E. Jansenc.... --------------------->  Element Metrics  <-----------------------
376*59599516SKenneth E. Jansenc
377*59599516SKenneth E. Jansen
378*59599516SKenneth E. Jansen      call e3metric( xl,         shgl,        dxidx,
379*59599516SKenneth E. Jansen     &               shg,        WdetJ)
380*59599516SKenneth E. Jansen
381*59599516SKenneth E. Jansenc
382*59599516SKenneth E. Jansenc.... compute the global gradient of u and P
383*59599516SKenneth E. Jansenc
384*59599516SKenneth E. Jansenc
385*59599516SKenneth E. Jansen       gradS = zero
386*59599516SKenneth E. Jansen       do n = 1, nshl
387*59599516SKenneth E. Jansen          gradS(:,1) = gradS(:,1) + shg(:,n,1) * yl(:,n,id)
388*59599516SKenneth E. Jansen          gradS(:,2) = gradS(:,2) + shg(:,n,2) * yl(:,n,id)
389*59599516SKenneth E. Jansen          gradS(:,3) = gradS(:,3) + shg(:,n,3) * yl(:,n,id)
390*59599516SKenneth E. Jansen       enddo
391*59599516SKenneth E. Jansen
392*59599516SKenneth E. Jansen       divS = zero
393*59599516SKenneth E. Jansen       if ( idiff >= 1 ) then
394*59599516SKenneth E. Jansenc
395*59599516SKenneth E. Jansenc.... compute divergence of diffusive flux vector, qi,i
396*59599516SKenneth E. Jansenc
397*59599516SKenneth E. Jansen          do n=1, nshl
398*59599516SKenneth E. Jansen             divS(:) = divS(:) + shg(:,n,1)*ql(:,n,1 )
399*59599516SKenneth E. Jansen     &                         + shg(:,n,2)*ql(:,n,2 )
400*59599516SKenneth E. Jansen     &                         + shg(:,n,3)*ql(:,n,3 )
401*59599516SKenneth E. Jansen          enddo
402*59599516SKenneth E. Jansen       endif                    ! diffusive flux computation
403*59599516SKenneth E. Jansen
404*59599516SKenneth E. Jansen       if(consrv_sclr_conv_vel.eq.1) then
405*59599516SKenneth E. Jansenc         Calculate uBar = u - TauM*L, where TauM is the momentum
406*59599516SKenneth E. Jansenc         stabilization factor and L is the momentum residual
407*59599516SKenneth E. Jansen
408*59599516SKenneth E. Jansen          if(matflg(5,1).eq.2) then ! boussinesq body force
409*59599516SKenneth E. Jansen             Temp = zero
410*59599516SKenneth E. Jansen             do n = 1, nshl
411*59599516SKenneth E. Jansen                Temp = Temp + shpfun(:,n) * yl(:,n,5)
412*59599516SKenneth E. Jansen             enddo
413*59599516SKenneth E. Jansen          endif
414*59599516SKenneth E. Jansen          if(matflg(5,1).eq.3.or.matflg(6,1).eq.1) then
415*59599516SKenneth E. Jansenc     user-specified body force or coriolis force specified
416*59599516SKenneth E. Jansen             xx = zero
417*59599516SKenneth E. Jansen             do n  = 1,nenl
418*59599516SKenneth E. Jansen                xx(:,1) = xx(:,1)  + shpfun(:,n) * xl(:,n,1)
419*59599516SKenneth E. Jansen                xx(:,2) = xx(:,2)  + shpfun(:,n) * xl(:,n,2)
420*59599516SKenneth E. Jansen                xx(:,3) = xx(:,3)  + shpfun(:,n) * xl(:,n,3)
421*59599516SKenneth E. Jansen             enddo
422*59599516SKenneth E. Jansen          endif
423*59599516SKenneth E. Jansen          aci = zero
424*59599516SKenneth E. Jansen          do n = 1, nshl
425*59599516SKenneth E. Jansen             aci(:,1) = aci(:,1) + shpfun(:,n) * acl(:,n,2)
426*59599516SKenneth E. Jansen             aci(:,2) = aci(:,2) + shpfun(:,n) * acl(:,n,3)
427*59599516SKenneth E. Jansen             aci(:,3) = aci(:,3) + shpfun(:,n) * acl(:,n,4)
428*59599516SKenneth E. Jansen          enddo
429*59599516SKenneth E. Jansen          g1yi = zero
430*59599516SKenneth E. Jansen          g2yi = zero
431*59599516SKenneth E. Jansen          g3yi = zero
432*59599516SKenneth E. Jansen          do n = 1, nshl
433*59599516SKenneth E. Jansen             g1yi(:,1) = g1yi(:,1) + shg(:,n,1) * yl(:,n,1)
434*59599516SKenneth E. Jansen             g1yi(:,2) = g1yi(:,2) + shg(:,n,1) * yl(:,n,2)
435*59599516SKenneth E. Jansen             g1yi(:,3) = g1yi(:,3) + shg(:,n,1) * yl(:,n,3)
436*59599516SKenneth E. Jansen             g1yi(:,4) = g1yi(:,4) + shg(:,n,1) * yl(:,n,4)
437*59599516SKenneth E. Jansenc
438*59599516SKenneth E. Jansen             g2yi(:,1) = g2yi(:,1) + shg(:,n,2) * yl(:,n,1)
439*59599516SKenneth E. Jansen             g2yi(:,2) = g2yi(:,2) + shg(:,n,2) * yl(:,n,2)
440*59599516SKenneth E. Jansen             g2yi(:,3) = g2yi(:,3) + shg(:,n,2) * yl(:,n,3)
441*59599516SKenneth E. Jansen             g2yi(:,4) = g2yi(:,4) + shg(:,n,2) * yl(:,n,4)
442*59599516SKenneth E. Jansenc
443*59599516SKenneth E. Jansen             g3yi(:,1) = g3yi(:,1) + shg(:,n,3) * yl(:,n,1)
444*59599516SKenneth E. Jansen             g3yi(:,2) = g3yi(:,2) + shg(:,n,3) * yl(:,n,2)
445*59599516SKenneth E. Jansen             g3yi(:,3) = g3yi(:,3) + shg(:,n,3) * yl(:,n,3)
446*59599516SKenneth E. Jansen             g3yi(:,4) = g3yi(:,4) + shg(:,n,3) * yl(:,n,4)
447*59599516SKenneth E. Jansen          enddo
448*59599516SKenneth E. Jansenc
449*59599516SKenneth E. Jansen          if (iLSet .eq. 0)then
450*59599516SKenneth E. Jansen             rho  = datmat(1,1,1)
451*59599516SKenneth E. Jansen             rmu = datmat(1,2,1)
452*59599516SKenneth E. Jansen          else
453*59599516SKenneth E. Jansen             write(*,*) 'Not sure if we can handle level set with K-E'
454*59599516SKenneth E. Jansen             write(*,*) '(different uMods? correct value of rho?)'
455*59599516SKenneth E. Jansen          endif
456*59599516SKenneth E. Jansen          divqi=zero  ! until we reconstruct q_flow for scalar solve
457*59599516SKenneth E. Jansen          call e3resStrongPDE(
458*59599516SKenneth E. Jansen     &         aci,  u1,   u2,   u3,   Temp, rho,  x,
459*59599516SKenneth E. Jansen     &               g1yi, g2yi, g3yi,
460*59599516SKenneth E. Jansen     &         rLui, src, divqi)
461*59599516SKenneth E. Jansen          src(:,1)=u1           !
462*59599516SKenneth E. Jansen          src(:,2)=u2           ! store u in src memory
463*59599516SKenneth E. Jansen          src(:,3)=u3           !
464*59599516SKenneth E. Jansenc         e3uBar calculates Tau_M and assembles uBar
465*59599516SKenneth E. Jansen          call getdiff(dwl, yl, shpfun, xmudmi, xl, rmu, rho)
466*59599516SKenneth E. Jansen          call e3uBar(rho, src, dxidx, rLui, rmu, uBar)
467*59599516SKenneth E. Jansen          u1=ubar(:,1)          ! the entire scalar residual
468*59599516SKenneth E. Jansen          u2=ubar(:,2)          ! is based on the modified
469*59599516SKenneth E. Jansen          u3=ubar(:,3)          ! velocity for conservation
470*59599516SKenneth E. Jansen       endif
471*59599516SKenneth E. Jansenc
472*59599516SKenneth E. Jansenc.... Initialize uMod, the modified velocity uMod
473*59599516SKenneth E. Jansenc      We initialize it to u_i and then calculate
474*59599516SKenneth E. Jansenc      the correction in e3sourcesclr
475*59599516SKenneth E. Jansenc
476*59599516SKenneth E. Jansen
477*59599516SKenneth E. Jansen       umod(:,1) = u1
478*59599516SKenneth E. Jansen       umod(:,2) = u2
479*59599516SKenneth E. Jansen       umod(:,3) = u3
480*59599516SKenneth E. Jansenc
481*59599516SKenneth E. Jansenc.... compute  source terms
482*59599516SKenneth E. Jansenc
483*59599516SKenneth E. Jansencad
484*59599516SKenneth E. Jansencad    if we are solving the redistancing equation, the umod(:,:) are
485*59599516SKenneth E. JansenCAD    modified in e3sourceSclr.
486*59599516SKenneth E. JansenCAD
487*59599516SKenneth E. JansenCAD  if we are redistancing levelset variable we want to use a use the
488*59599516SKenneth E. JansenCAD  convective term from the equation.
489*59599516SKenneth E. Jansen
490*59599516SKenneth E. Jansen
491*59599516SKenneth E. Jansen       if(nosource.ne.1) then
492*59599516SKenneth E. Jansen        call e3sourceSclr ( Sclr,         Sdot,      gradS,  dwl,
493*59599516SKenneth E. Jansen     &                      shpfun,       shg,       yl,     dxidx,
494*59599516SKenneth E. Jansen     &                      diffus,       u1,        u2,     u3,
495*59599516SKenneth E. Jansen     &                      xl,           srcR,      srcL,   uMod,
496*59599516SKenneth E. Jansen     &                      srcRat)
497*59599516SKenneth E. Jansen       else
498*59599516SKenneth E. Jansen        srcRat = zero
499*59599516SKenneth E. Jansen        srcR   = zero
500*59599516SKenneth E. Jansen        srcL   = zero
501*59599516SKenneth E. Jansen       endif
502*59599516SKenneth E. Jansenc
503*59599516SKenneth E. Jansenc.... -------------------> Scalar residual  <-----------------
504*59599516SKenneth E. Jansenc
505*59599516SKenneth E. Jansen
506*59599516SKenneth E. Jansen         rLS(:) = ( Sdot(:) +  (u1*gradS(:,1) +
507*59599516SKenneth E. Jansen     &                              u2*gradS(:,2) +
508*59599516SKenneth E. Jansen     &                              u3*gradS(:,3)) )
509*59599516SKenneth E. Jansen     &        - divS(:)
510*59599516SKenneth E. Jansen
511*59599516SKenneth E. Jansenc
512*59599516SKenneth E. Jansenc.... return
513*59599516SKenneth E. Jansenc
514*59599516SKenneth E. Jansen       return
515*59599516SKenneth E. Jansen       end
516*59599516SKenneth E. Jansen
517