xref: /phasta/phSolver/incompressible/e3source/e3source.fRHSonly (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen      subroutine e3source(xx, src)
2*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
3*59599516SKenneth E. Jansenc
4*59599516SKenneth E. Jansenc  this routine computes the body force term.
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc  currently this computes a swirl body with the axis alligned with
7*59599516SKenneth E. Jansenc  the z-coordinate
8*59599516SKenneth E. Jansenc
9*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
10*59599516SKenneth E. Jansen
11*59599516SKenneth E. Jansen      include "common.h"
12*59599516SKenneth E. Jansen
13*59599516SKenneth E. Jansen      real*8   xx(npro,nsd), src(npro,nsd)
14*59599516SKenneth E. Jansen
15*59599516SKenneth E. Jansen      real*8   nu
16*59599516SKenneth E. Jansen
17*59599516SKenneth E. Jansen      real*8   r, Stheta, dpdz, rP5
18*59599516SKenneth E. Jansen      if(datmat(2,5,1).eq. 0) then  ! calc swirl
19*59599516SKenneth E. Jansenc
20*59599516SKenneth E. Jansenc  This is the body force which will drive a swirl in a pipe flow
21*59599516SKenneth E. Jansenc
22*59599516SKenneth E. Jansen         bigR    = 0.5d0
23*59599516SKenneth E. Jansen         dpdz    = datmat(1,5,1)
24*59599516SKenneth E. Jansen         do iel = 1, npro
25*59599516SKenneth E. Jansen
26*59599516SKenneth E. Jansen            r   = sqrt( xx(iel,1)**2 + xx(iel,2)**2)
27*59599516SKenneth E. Jansen            rP5 = (r/bigR)**5
28*59599516SKenneth E. Jansen
29*59599516SKenneth E. Jansen            Stheta = dpdz * sin(0.5*pi*rP5)
30*59599516SKenneth E. Jansen
31*59599516SKenneth E. Jansen            src(iel,1) = -xx(iel,2)/r * Stheta
32*59599516SKenneth E. Jansen            src(iel,2) =  xx(iel,1)/r * Stheta
33*59599516SKenneth E. Jansen            src(iel,3) =  dpdz
34*59599516SKenneth E. Jansen         enddo
35*59599516SKenneth E. Jansen      else  ! contrived test problem
36*59599516SKenneth E. Jansen
37*59599516SKenneth E. Jansenc$$$         nu=datmat(1,2,1)
38*59599516SKenneth E. Jansenc$$$         omag=datmat(3,5,1) ! frame rotation rate
39*59599516SKenneth E. Jansenc$$$         e1 = one/sqrt(3.0d0) ! axis of rotation
40*59599516SKenneth E. Jansenc$$$         e2 = e1
41*59599516SKenneth E. Jansenc$$$         e3 = e1
42*59599516SKenneth E. Jansenc$$$         om1=omag*e1
43*59599516SKenneth E. Jansenc$$$         om2=omag*e2
44*59599516SKenneth E. Jansenc$$$         om3=omag*e3
45*59599516SKenneth E. Jansenc$$$         w=0
46*59599516SKenneth E. Jansenc$$$
47*59599516SKenneth E. Jansenc$$$         do iel = 1, npro
48*59599516SKenneth E. Jansenc$$$
49*59599516SKenneth E. Jansenc$$$            x = xx(iel,1)
50*59599516SKenneth E. Jansenc$$$            y = xx(iel,2)
51*59599516SKenneth E. Jansenc$$$            z = xx(iel,3)
52*59599516SKenneth E. Jansenc$$$
53*59599516SKenneth E. Jansenc$$$c
54*59599516SKenneth E. Jansenc$$$c  The following are MAPLE generated forcing functions for
55*59599516SKenneth E. Jansenc$$$c  a contrived flow in a rotating reference frame.
56*59599516SKenneth E. Jansenc$$$c
57*59599516SKenneth E. Jansenc$$$
58*59599516SKenneth E. Jansenc$$$      t1 = x**2
59*59599516SKenneth E. Jansenc$$$      t2 = t1**2
60*59599516SKenneth E. Jansenc$$$      t5 = dexp(14.D0*x)
61*59599516SKenneth E. Jansenc$$$      t6 = t2*x*t5
62*59599516SKenneth E. Jansenc$$$      t7 = y**2
63*59599516SKenneth E. Jansenc$$$      t8 = t7**2
64*59599516SKenneth E. Jansenc$$$      t9 = t8*t7
65*59599516SKenneth E. Jansenc$$$      t12 = t2*t5
66*59599516SKenneth E. Jansenc$$$      t15 = nu*t1
67*59599516SKenneth E. Jansenc$$$      t17 = dexp(7.D0*x)
68*59599516SKenneth E. Jansenc$$$      t18 = t7*y
69*59599516SKenneth E. Jansenc$$$      t19 = t17*t18
70*59599516SKenneth E. Jansenc$$$      t22 = t1*x
71*59599516SKenneth E. Jansenc$$$      t23 = nu*t22
72*59599516SKenneth E. Jansenc$$$      t24 = t17*y
73*59599516SKenneth E. Jansenc$$$      t29 = t17*t7
74*59599516SKenneth E. Jansenc$$$      t34 = nu*x
75*59599516SKenneth E. Jansenc$$$      t43 = nu*t2
76*59599516SKenneth E. Jansenc$$$      t46 = -40.D0*t6*t9-6.D0*t12*t7+92.D0*t15*t19+132.D0*t23*t24+80.D0*
77*59599516SKenneth E. Jansenc$$$     #t6*t18-252.D0*t23*t29+168.D0*t23*t19+96.D0*t34*t29-64.D0*t34*t19+2
78*59599516SKenneth E. Jansenc$$$     #2.D0*t15*t24-138.D0*t15*t29+294.D0*t43*t29
79*59599516SKenneth E. Jansenc$$$      t50 = t2*t22*t5
80*59599516SKenneth E. Jansenc$$$      t57 = nu*t17
81*59599516SKenneth E. Jansenc$$$      t60 = t8*y
82*59599516SKenneth E. Jansenc$$$      t65 = t2**2
83*59599516SKenneth E. Jansenc$$$      t66 = t65*t5
84*59599516SKenneth E. Jansenc$$$      t73 = t22*t5
85*59599516SKenneth E. Jansenc$$$      t77 = t2*t1*t5
86*59599516SKenneth E. Jansenc$$$      t80 = -196.D0*t43*t19-96.D0*t50*t9-122.D0*t43*t24-32.D0*t34*t24-4.
87*59599516SKenneth E. Jansenc$$$     #D0*t57*y+36.D0*t12*t60+192.D0*t50*t18+98.D0*t66*t8+24.D0*t12*t18+2
88*59599516SKenneth E. Jansenc$$$     #8.D0*t66*t9-24.D0*t73*t60-336.D0*t77*t60
89*59599516SKenneth E. Jansenc$$$      t106 = -336.D0*t50*t8-12.D0*t12*t9+8.D0*t73*t9-140.D0*t6*t8+28.D0*
90*59599516SKenneth E. Jansenc$$$     #t73*t8+12.D0*t15*t17+12.D0*t43*t17+288.D0*t50*t60+112.D0*t77*t9-22
91*59599516SKenneth E. Jansenc$$$     #4.D0*t77*t18-56.D0*t66*t18+120.D0*t6*t60
92*59599516SKenneth E. Jansenc$$$      t131 = -8.D0*t57*t18-20.D0*t6*t7+56.D0*t77*t7-42.D0*t12*t8-24.D0*t
93*59599516SKenneth E. Jansenc$$$     #23*t17-48.D0*t50*t7+4.D0*t73*t7-16.D0*t73*t18+392.D0*t77*t8+14.D0*
94*59599516SKenneth E. Jansenc$$$     #t66*t7-84.D0*t66*t60+12.D0*t57*t7
95*59599516SKenneth E. Jansenc$$$      fx = t46+t80+t106+t131
96*59599516SKenneth E. Jansenc$$$
97*59599516SKenneth E. Jansenc$$$
98*59599516SKenneth E. Jansenc$$$      t1 = x**2
99*59599516SKenneth E. Jansenc$$$      t2 = t1**2
100*59599516SKenneth E. Jansenc$$$      t3 = nu*t2
101*59599516SKenneth E. Jansenc$$$      t5 = dexp(7.D0*x)
102*59599516SKenneth E. Jansenc$$$      t6 = t5*y
103*59599516SKenneth E. Jansenc$$$      t9 = y**2
104*59599516SKenneth E. Jansenc$$$      t10 = t9**2
105*59599516SKenneth E. Jansenc$$$      t11 = nu*t10
106*59599516SKenneth E. Jansenc$$$      t18 = t1*x
107*59599516SKenneth E. Jansenc$$$      t22 = nu*x
108*59599516SKenneth E. Jansenc$$$      t25 = nu*t18
109*59599516SKenneth E. Jansenc$$$      t32 = dexp(14.D0*x)
110*59599516SKenneth E. Jansenc$$$      t33 = t2*t32
111*59599516SKenneth E. Jansenc$$$      t39 = t2*t1*t32
112*59599516SKenneth E. Jansenc$$$      t40 = t10*t9
113*59599516SKenneth E. Jansenc$$$      t43 = t1*t32
114*59599516SKenneth E. Jansenc$$$      t46 = t9*y
115*59599516SKenneth E. Jansenc$$$      t47 = t10*t46
116*59599516SKenneth E. Jansenc$$$      t50 = -84.D0*t3*t6+343.D0*t11*t5*t2+66.D0*t11*t5*x-98.D0*t11*t5*t1
117*59599516SKenneth E. Jansenc$$$     #8-24.D0*t22*t6+120.D0*t25*t6-287.D0*t11*t5*t1-140.D0*t33*t10+14.D0
118*59599516SKenneth E. Jansenc$$$     #*t3*t5-56.D0*t39*t40-28.D0*t43*t40+8.D0*t43*t47
119*59599516SKenneth E. Jansenc$$$      t52 = t2*x*t32
120*59599516SKenneth E. Jansenc$$$      t55 = t18*t32
121*59599516SKenneth E. Jansenc$$$      t62 = t10*y
122*59599516SKenneth E. Jansenc$$$      t69 = nu*t5
123*59599516SKenneth E. Jansenc$$$      t80 = 120.D0*t52*t10-32.D0*t55*t47+56.D0*t33*t47+30.D0*t11*t5-216.
124*59599516SKenneth E. Jansenc$$$     #D0*t52*t62-144.D0*t55*t62+72.D0*t39*t62-60.D0*t69*t46+30.D0*t69*t9
125*59599516SKenneth E. Jansenc$$$     #-20.D0*t25*t5-20.D0*t43*t10+36.D0*t43*t62
126*59599516SKenneth E. Jansenc$$$      t86 = nu*t1
127*59599516SKenneth E. Jansenc$$$      t89 = t5*t9
128*59599516SKenneth E. Jansenc$$$      t92 = t5*t46
129*59599516SKenneth E. Jansenc$$$      t109 = 112.D0*t55*t40+80.D0*t55*t10-12.D0*t86*t6-275.D0*t86*t89+57
130*59599516SKenneth E. Jansenc$$$     #4.D0*t86*t92-218.D0*t25*t89+196.D0*t25*t92+90.D0*t22*t89-132.D0*t2
131*59599516SKenneth E. Jansenc$$$     #2*t92+427.D0*t3*t89+8.D0*t39*t46+4.D0*t43*t46
132*59599516SKenneth E. Jansenc$$$      t134 = 16.D0*t39*t47+28.D0*t33*t46-48.D0*t52*t47+252.D0*t33*t62+16
133*59599516SKenneth E. Jansenc$$$     #8.D0*t52*t40-24.D0*t52*t46-686.D0*t3*t92+2.D0*t86*t5-40.D0*t39*t10
134*59599516SKenneth E. Jansenc$$$     #-196.D0*t33*t40-16.D0*t55*t46+4.D0*t22*t5
135*59599516SKenneth E. Jansenc$$$      fy = t50+t80+t109+t134
136*59599516SKenneth E. Jansenc$$$
137*59599516SKenneth E. Jansenc$$$
138*59599516SKenneth E. Jansenc$$$      t3 = om3*x
139*59599516SKenneth E. Jansenc$$$      t5 = dexp(7.D0*x)
140*59599516SKenneth E. Jansenc$$$      t7 = x**2
141*59599516SKenneth E. Jansenc$$$      t12 = y**2
142*59599516SKenneth E. Jansenc$$$      t15 = (-1.D0+y)**2
143*59599516SKenneth E. Jansenc$$$      fxa = 2.D0*om2*w+2.D0*t3*t5*(2.D0+x-10.D0*t7+7.D0*t7*x)*t12*t15+om
144*59599516SKenneth E. Jansenc$$$     #2*(om1*y-om2*x)-om3*(t3-om1*z)
145*59599516SKenneth E. Jansenc$$$
146*59599516SKenneth E. Jansenc$$$      t1 = x**2
147*59599516SKenneth E. Jansenc$$$      t4 = (-1.D0+x)**2
148*59599516SKenneth E. Jansenc$$$      t7 = dexp(7.D0*x)
149*59599516SKenneth E. Jansenc$$$      t10 = y**2
150*59599516SKenneth E. Jansenc$$$      fya = 4.D0*om3*t1*t4*t7*y*(1.D0-3.D0*y+2.D0*t10)-2.D0*om1*w+om3*(o
151*59599516SKenneth E. Jansenc$$$     #m2*z-om3*y)-om1*(om1*y-om2*x)
152*59599516SKenneth E. Jansenc$$$
153*59599516SKenneth E. Jansenc$$$
154*59599516SKenneth E. Jansenc$$$      t3 = dexp(7.D0*x)
155*59599516SKenneth E. Jansenc$$$      t5 = x**2
156*59599516SKenneth E. Jansenc$$$      t10 = y**2
157*59599516SKenneth E. Jansenc$$$      t13 = (-1.D0+y)**2
158*59599516SKenneth E. Jansenc$$$      t19 = (-1.D0+x)**2
159*59599516SKenneth E. Jansenc$$$      fza = -2.D0*om1*x*t3*(2.D0+x-10.D0*t5+7.D0*t5*x)*t10*t13-4.D0*om2*
160*59599516SKenneth E. Jansenc$$$     #t5*t19*t3*y*(1.D0-3.D0*y+2.D0*t10)+om1*(om3*x-om1*z)-om2*(om2*z-om
161*59599516SKenneth E. Jansenc$$$     #3*y)
162*59599516SKenneth E. Jansenc$$$
163*59599516SKenneth E. Jansenc$$$      src(iel,1) =  fx + fxa
164*59599516SKenneth E. Jansenc$$$      src(iel,2) =  fy + fya
165*59599516SKenneth E. Jansenc$$$      src(iel,3) =       fza
166*59599516SKenneth E. Jansenc$$$      enddo
167*59599516SKenneth E. Jansen!Analytic lid case
168*59599516SKenneth E. Jansen
169*59599516SKenneth E. Jansen      do iel = 1, npro
170*59599516SKenneth E. Jansen            x = xx(iel,1)
171*59599516SKenneth E. Jansen            y = xx(iel,2)
172*59599516SKenneth E. Jansen            z = xx(iel,3)
173*59599516SKenneth E. Jansen            Re = 1.0/datmat(1,2,1)
174*59599516SKenneth E. Jansenc
175*59599516SKenneth E. Jansenc  The following are MAPLE generated forcing functions for
176*59599516SKenneth E. Jansenc  a Lid Driven cavity flow with an analytic solution
177*59599516SKenneth E. Jansenc
178*59599516SKenneth E. Jansen            t2 = x**2
179*59599516SKenneth E. Jansen            t3 = t2**2
180*59599516SKenneth E. Jansen            t4 = t3*x
181*59599516SKenneth E. Jansen            t7 = t2*x
182*59599516SKenneth E. Jansen            t8 = 8*t7
183*59599516SKenneth E. Jansen            t13 = y**2
184*59599516SKenneth E. Jansen            t20 = t13**2
185*59599516SKenneth E. Jansen            t21 = t20-t13
186*59599516SKenneth E. Jansen            t26 = t3**2
187*59599516SKenneth E. Jansen            t30 = t3*t2
188*59599516SKenneth E. Jansen            t37 = t13*y
189*59599516SKenneth E. Jansen            t54 = -8/Re*(24.E0/5.E0*t4-12*t3+t8+2*(4*t7-6*t2+2*x)*(12
190*59599516SKenneth E. Jansen     &           *t13-2)+(24*x-12)*t21)-64*(t26/2-2*t3*t7+3*t30-2*t4+t3
191*59599516SKenneth E. Jansen     &           /2)*(-24*t20*y+8*t37-4*y)+64*t21*(4*t37-2*y)*(-4*t30+12
192*59599516SKenneth E. Jansen     &           *t4-14*t3+t8-2*t2)
193*59599516SKenneth E. Jansen
194*59599516SKenneth E. Jansen            src(iel,1) = 0.0
195*59599516SKenneth E. Jansen            src(iel,2) = t54
196*59599516SKenneth E. Jansen            src(iel,3) = 0.0
197*59599516SKenneth E. Jansen         enddo
198*59599516SKenneth E. Jansen
199*59599516SKenneth E. Jansen      endif
200*59599516SKenneth E. Jansen      return
201*59599516SKenneth E. Jansen      end
202*59599516SKenneth E. Jansen
203*59599516SKenneth E. Jansen
204*59599516SKenneth E. Jansen!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
205*59599516SKenneth E. Jansen!******************************************************************
206*59599516SKenneth E. Jansen!-------------------------------------------------------------------
207*59599516SKenneth E. Jansen
208*59599516SKenneth E. Jansen
209*59599516SKenneth E. Jansen
210*59599516SKenneth E. Jansen      subroutine e3sourceSclr ( Sclr,      Sdot,      gradS,
211*59599516SKenneth E. Jansen     &                          dwl,       shape_funct,     shg,
212*59599516SKenneth E. Jansen     &                          yl,        dxidx,     rmu,
213*59599516SKenneth E. Jansen     &                          u1,        u2,        u3,   xl,
214*59599516SKenneth E. Jansen     &                          srcR,      srcL,      uMod,
215*59599516SKenneth E. Jansen     &                          srcRat)
216*59599516SKenneth E. Jansen
217*59599516SKenneth E. Jansen
218*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
219*59599516SKenneth E. Jansenc
220*59599516SKenneth E. Jansenc     calculate the coefficients for the Spalart-Allmaras
221*59599516SKenneth E. Jansenc     turbulence model and its jacobian
222*59599516SKenneth E. Jansenc
223*59599516SKenneth E. Jansenc     input:
224*59599516SKenneth E. Jansenc             Sclr  :   The scalar that is being transported in this pde.
225*59599516SKenneth E. Jansenc             gradS :   spatial gradient of Sclr
226*59599516SKenneth E. Jansenc             rmu   :   kinematic viscosity
227*59599516SKenneth E. Jansenc
228*59599516SKenneth E. Jansenc     output:
229*59599516SKenneth E. Jansenc             rmu   :   diffusion for eddy viscosity equation
230*59599516SKenneth E. Jansenc             srcR  :   residual terms for
231*59599516SKenneth E. Jansenc                         (srcR * turb)
232*59599516SKenneth E. Jansenc             srcL :   jacobian of srcR
233*59599516SKenneth E. Jansenc
234*59599516SKenneth E. Jansenc-----------------------------------------------------------------------
235*59599516SKenneth E. Jansen      use     turbSA
236*59599516SKenneth E. Jansen      use turbKE
237*59599516SKenneth E. Jansen      include "common.h"
238*59599516SKenneth E. Jansenc coming in
239*59599516SKenneth E. Jansen      real*8  Sclr(npro),          Sdot(npro),
240*59599516SKenneth E. Jansen     &        gradS(npro,nsd),     dwl(npro,nenl),
241*59599516SKenneth E. Jansen     &        shape_funct(npro,nshl),    shg(npro,nshl,nsd),
242*59599516SKenneth E. Jansen     &        yl(npro,nenl,ndof),  dxidx(npro,nsd,nsd),
243*59599516SKenneth E. Jansen     &        rmu(npro),           u1(npro),
244*59599516SKenneth E. Jansen     &        u2(npro),            u3(npro),
245*59599516SKenneth E. Jansen     &        xl(npro,nenl,nsd)
246*59599516SKenneth E. Jansenc going out
247*59599516SKenneth E. Jansen      real*8  srcR(npro),          srcL(npro),
248*59599516SKenneth E. Jansen     &        uMod(npro,nsd)
249*59599516SKenneth E. Jansenc used locally
250*59599516SKenneth E. Jansen
251*59599516SKenneth E. Jansen      real*8  gradv(npro,nsd,nsd), absVort(npro),
252*59599516SKenneth E. Jansen     &        dwall(npro)
253*59599516SKenneth E. Jansen      real*8  chi,    chiP3,   fv1,   fv2,     st,    r,
254*59599516SKenneth E. Jansen     &        g,      fw,      s,     viscInv, k2d2Inv,
255*59599516SKenneth E. Jansen     &        dP2Inv, sixth,   tmp,   tmp1,    p,     dp,
256*59599516SKenneth E. Jansen     &        fv1p,   fv2p,    stp,         gp,      fwp,   rp,
257*59599516SKenneth E. Jansen     &        chiP2,  mytmp(npro),         sign_levelset(npro),
258*59599516SKenneth E. Jansen     &        sclr_ls(npro)
259*59599516SKenneth E. Jansenc    Kay-Epsilon
260*59599516SKenneth E. Jansen      real*8 Ry(npro), Rt(npro), RtP2(npro), f2(npro), f1(npro),
261*59599516SKenneth E. Jansen     &       kay(npro), epsilon(npro), fmu(npro), fmui(npro),
262*59599516SKenneth E. Jansen     &       srcjac(npro,4), srcRat(npro)
263*59599516SKenneth E. Jansen      integer e
264*59599516SKenneth E. Jansenc----------------------------------------------------------------------
265*59599516SKenneth E. Jansenc
266*59599516SKenneth E. Jansenc           --             --              --           --
267*59599516SKenneth E. Jansenc           |      cb2      |           1  |             |
268*59599516SKenneth E. Jansenc   phi,  + |u  - -----phi, | phi,  - -----|(nu+phi)phi, |
269*59599516SKenneth E. Jansenc       t   | i   sigma    i|     i   sigma|            i|
270*59599516SKenneth E. Jansenc           --             --              --           --,
271*59599516SKenneth E. Jansenc                                                          i
272*59599516SKenneth E. Jansenc
273*59599516SKenneth E. Jansenc                  ~              phi^2
274*59599516SKenneth E. Jansenc            - cb1*S*phi + cw1*fw*-----
275*59599516SKenneth E. Jansenc                                  d^2
276*59599516SKenneth E. Jansenc
277*59599516SKenneth E. Jansenc----------------------------------------------------------------------
278*59599516SKenneth E. Jansenc
279*59599516SKenneth E. Jansen      if(iRANS.eq.-1) then    ! spalart almaras model
280*59599516SKenneth E. Jansenc
281*59599516SKenneth E. Jansenc.... compute the global gradient of u
282*59599516SKenneth E. Jansenc
283*59599516SKenneth E. Jansen         gradV = zero
284*59599516SKenneth E. Jansen         do n = 1, nshl
285*59599516SKenneth E. Jansenc
286*59599516SKenneth E. Jansenc          du_i/dx_j
287*59599516SKenneth E. Jansenc
288*59599516SKenneth E. Jansenc           i j   indices match array where V is the velocity (u in our notes)
289*59599516SKenneth E. Jansen            gradV(:,1,1) = gradV(:,1,1) + shg(:,n,1) * yl(:,n,2)
290*59599516SKenneth E. Jansen            gradV(:,2,1) = gradV(:,2,1) + shg(:,n,1) * yl(:,n,3)
291*59599516SKenneth E. Jansen            gradV(:,3,1) = gradV(:,3,1) + shg(:,n,1) * yl(:,n,4)
292*59599516SKenneth E. Jansenc
293*59599516SKenneth E. Jansen            gradV(:,1,2) = gradV(:,1,2) + shg(:,n,2) * yl(:,n,2)
294*59599516SKenneth E. Jansen            gradV(:,2,2) = gradV(:,2,2) + shg(:,n,2) * yl(:,n,3)
295*59599516SKenneth E. Jansen            gradV(:,3,2) = gradV(:,3,2) + shg(:,n,2) * yl(:,n,4)
296*59599516SKenneth E. Jansenc
297*59599516SKenneth E. Jansen            gradV(:,1,3) = gradV(:,1,3) + shg(:,n,3) * yl(:,n,2)
298*59599516SKenneth E. Jansen            gradV(:,2,3) = gradV(:,2,3) + shg(:,n,3) * yl(:,n,3)
299*59599516SKenneth E. Jansen            gradV(:,3,3) = gradV(:,3,3) + shg(:,n,3) * yl(:,n,4)
300*59599516SKenneth E. Jansenc                                             a j     u   a i
301*59599516SKenneth E. Jansenc from our notes where we had N_{a,j} = dN_a/dx_j  note that i is off by one because p was first in yl vector
302*59599516SKenneth E. Jansenc
303*59599516SKenneth E. Jansen         enddo
304*59599516SKenneth E. Jansenc
305*59599516SKenneth E. Jansenc.... magnitude of vorticity
306*59599516SKenneth E. Jansenc
307*59599516SKenneth E. Jansen         absVort  = sqrt( (gradV(:,2,3) - gradV(:,3,2)) ** 2
308*59599516SKenneth E. Jansen     &                  + (gradV(:,3,1) - gradV(:,1,3)) ** 2
309*59599516SKenneth E. Jansen     &		        + (gradV(:,1,2) - gradV(:,2,1)) ** 2 )
310*59599516SKenneth E. Jansen         dwall = zero
311*59599516SKenneth E. Jansen         do n = 1, nenl
312*59599516SKenneth E. Jansen            dwall     = dwall + shape_funct(:,n) * dwl(:,n)
313*59599516SKenneth E. Jansen         enddo
314*59599516SKenneth E. Jansen
315*59599516SKenneth E. Jansen         sixth   = 1.0/6.0
316*59599516SKenneth E. Jansenc
317*59599516SKenneth E. Jansenc.... compute source and its jacobian
318*59599516SKenneth E. Jansenc
319*59599516SKenneth E. Jansen         do e = 1, npro
320*59599516SKenneth E. Jansen            SclrNN= max(Sclr(e),zero)
321*59599516SKenneth E. Jansen
322*59599516SKenneth E. Jansen            if(iles.lt.0) then     ! for DES
323*59599516SKenneth E. Jansen               dx=maxval(xl(e,:,1))-minval(xl(e,:,1))
324*59599516SKenneth E. Jansen               dy=maxval(xl(e,:,2))-minval(xl(e,:,2))
325*59599516SKenneth E. Jansen               dz=maxval(xl(e,:,3))-minval(xl(e,:,3))
326*59599516SKenneth E. Jansen               dmax=max(dx,max(dy,dz))
327*59599516SKenneth E. Jansen               dmax=0.65*dmax
328*59599516SKenneth E. Jansenc
329*59599516SKenneth E. Jansenc Next 5 lines are DDES
330*59599516SKenneth E. Jansenc
331*59599516SKenneth E. Jansen         qfac    = sqrt(gradV(e,1,1)**2+gradV(e,1,2)**2+gradV(e,1,3)**2
332*59599516SKenneth E. Jansen     &                 +gradV(e,2,1)**2+gradV(e,2,2)**2+gradV(e,2,3)**2
333*59599516SKenneth E. Jansen     &                 +gradV(e,3,1)**2+gradV(e,3,2)**2+gradV(e,3,3)**2)
334*59599516SKenneth E. Jansen	       rd=SclrNN*saKappaP2Inv/(dwall(e)**2*qfact)
335*59599516SKenneth E. Jansen               fd=one-tanh((8.0000000000000000d0*rd)**3)
336*59599516SKenneth E. Jansen               dwall(e)=dwall(e)-fd*max(zero,dwall(e)-dmax)
337*59599516SKenneth E. Jansenc DES97               dwall(e)=min(dmax,dwall(e))
338*59599516SKenneth E. Jansen            endif
339*59599516SKenneth E. Jansen
340*59599516SKenneth E. Jansen            if (dwall(e) .ne. 0) dP2Inv = 1.0 / dwall(e)**2
341*59599516SKenneth E. Jansen
342*59599516SKenneth E. Jansen            k2d2Inv = dP2Inv * saKappaP2Inv
343*59599516SKenneth E. Jansen
344*59599516SKenneth E. Jansen            viscInv = 1.0/datmat(1,2,1)
345*59599516SKenneth E. Jansen            chi     = SclrNN * viscInv !skipping the kiy kie for now
346*59599516SKenneth E. Jansen            chiP2   = chi * chi
347*59599516SKenneth E. Jansen            chiP3   = chi * chiP2
348*59599516SKenneth E. Jansen
349*59599516SKenneth E. Jansen            ep=zero
350*59599516SKenneth E. Jansen            if(Sclr(e).gt.zero) ep=one
351*59599516SKenneth E. Jansen
352*59599516SKenneth E. Jansen            tmp     = 1.0 / ( chiP3 + saCv1P3 )
353*59599516SKenneth E. Jansen            fv1     = chiP3 * tmp
354*59599516SKenneth E. Jansen            fv1p    = 3.0 * viscInv * saCv1P3 * chiP2 * tmp**2
355*59599516SKenneth E. Jansen
356*59599516SKenneth E. Jansen            tmp     = 1.0 / (1.0 + chi * fv1)
357*59599516SKenneth E. Jansen            fv2     = 1.0 - chi * tmp
358*59599516SKenneth E. Jansen            fv2p    = (chiP2 * fv1p - viscInv) * tmp**2
359*59599516SKenneth E. Jansen
360*59599516SKenneth E. Jansen            s       = absVort(e) !prd(e,2)
361*59599516SKenneth E. Jansen            tmp     = SclrNN * k2d2Inv !eOkdP2
362*59599516SKenneth E. Jansen            st      = s + fv2 * tmp
363*59599516SKenneth E. Jansen            stp     = ep*k2d2Inv * fv2 + tmp*fv2p
364*59599516SKenneth E. Jansen
365*59599516SKenneth E. Jansen            r       = zero
366*59599516SKenneth E. Jansen            rp      = zero
367*59599516SKenneth E. Jansen            if (st .gt. epsM ) then
368*59599516SKenneth E. Jansen                r = tmp / st
369*59599516SKenneth E. Jansen                r       = min( max(r, -8.0d0), 8.0d0)
370*59599516SKenneth E. Jansen                rp = k2d2Inv / st**2 * (ep*st - SclrNN*stp)
371*59599516SKenneth E. Jansen            endif
372*59599516SKenneth E. Jansen            rP5     = r * (r * r) ** 2
373*59599516SKenneth E. Jansen            tmp     = 1.0 + saCw2 * (rP5 - 1.0)
374*59599516SKenneth E. Jansen            g       = r * tmp
375*59599516SKenneth E. Jansen            gp      = rp * (tmp + 5.0 * saCw2 * rP5)
376*59599516SKenneth E. Jansen
377*59599516SKenneth E. Jansen            gP6     = (g * g * g) ** 2
378*59599516SKenneth E. Jansen            tmp     = 1.0 / (gP6 + saCw3P6)
379*59599516SKenneth E. Jansen            tmp1    = ( (1.0 + saCw3P6) * tmp ) ** sixth
380*59599516SKenneth E. Jansen            fw      = g * tmp1
381*59599516SKenneth E. Jansen            fwp     = gp * tmp1 * saCw3P6 * tmp
382*59599516SKenneth E. Jansen            if(abs(r).gt.7.0d0) fwp=zero
383*59599516SKenneth E. Jansen
384*59599516SKenneth E. Jansen            tmp1     = saCw1 * dP2Inv*SclrNN
385*59599516SKenneth E. Jansen            prat     = ep*saCb1*st  !prodRat
386*59599516SKenneth E. Jansen            srcRat(e)= prat - ep*fw * tmp1
387*59599516SKenneth E. Jansen
388*59599516SKenneth E. Jansen            tmp     = zero
389*59599516SKenneth E. Jansen            if(prat .lt. zero) tmp=prat
390*59599516SKenneth E. Jansen            if(stp .lt. zero) tmp=tmp+SaCb1*stp*SclrNN
391*59599516SKenneth E. Jansen            if(fw   .gt. zero) tmp=tmp-two*tmp1*fw
392*59599516SKenneth E. Jansen            if(fwp  .gt. zero) tmp=tmp-tmp1*SclrNN*fwp
393*59599516SKenneth E. Jansen
394*59599516SKenneth E. Jansen            srcL(e)  = -1.0*tmp
395*59599516SKenneth E. Jansen            srcR(e)   = srcRat(e)*SclrNN
396*59599516SKenneth E. Jansen         enddo
397*59599516SKenneth E. Jansenc
398*59599516SKenneth E. Jansenc.... One source term has the form (beta_i)*(phi,_i).  Since
399*59599516SKenneth E. Jansenc     the convective term has (u_i)*(phi,_i), it is useful to treat
400*59599516SKenneth E. Jansenc     beta_i as a "correction" to the velocity.  In calculating the
401*59599516SKenneth E. Jansenc     stabilization terms, the new "modified" velocity (u_i-beta_i) is
402*59599516SKenneth E. Jansenc     then used in place of the pure velocity for stabilization terms,
403*59599516SKenneth E. Jansenc     and the source term sneaks into the RHS and LHS.  Here, the term
404*59599516SKenneth E. Jansenc     is given by beta_i=(cb_2)*(phi,_i)/(sigma)
405*59599516SKenneth E. Jansenc
406*59599516SKenneth E. Jansen
407*59599516SKenneth E. Jansen         tmp = saCb2 * saSigmaInv
408*59599516SKenneth E. Jansen         uMod(:,1) = u1(:) - tmp * gradS(:,1)
409*59599516SKenneth E. Jansen         uMod(:,2) = u2(:) - tmp * gradS(:,2)
410*59599516SKenneth E. Jansen         uMod(:,3) = u3(:) - tmp * gradS(:,3)
411*59599516SKenneth E. Jansen
412*59599516SKenneth E. Jansen      elseif(iRANS.eq.-2) then    ! K-epsilon
413*59599516SKenneth E. Jansenc.... compute the global gradient of u
414*59599516SKenneth E. Jansenc
415*59599516SKenneth E. Jansen         gradV = zero
416*59599516SKenneth E. Jansen         do n = 1, nshl
417*59599516SKenneth E. Jansenc
418*59599516SKenneth E. Jansenc          du_i/dx_j
419*59599516SKenneth E. Jansenc
420*59599516SKenneth E. Jansenc           i j   indices match array where V is the velocity (u in our notes)
421*59599516SKenneth E. Jansen            gradV(:,1,1) = gradV(:,1,1) + shg(:,n,1) * yl(:,n,2)
422*59599516SKenneth E. Jansen            gradV(:,2,1) = gradV(:,2,1) + shg(:,n,1) * yl(:,n,3)
423*59599516SKenneth E. Jansen            gradV(:,3,1) = gradV(:,3,1) + shg(:,n,1) * yl(:,n,4)
424*59599516SKenneth E. Jansenc
425*59599516SKenneth E. Jansen            gradV(:,1,2) = gradV(:,1,2) + shg(:,n,2) * yl(:,n,2)
426*59599516SKenneth E. Jansen            gradV(:,2,2) = gradV(:,2,2) + shg(:,n,2) * yl(:,n,3)
427*59599516SKenneth E. Jansen            gradV(:,3,2) = gradV(:,3,2) + shg(:,n,2) * yl(:,n,4)
428*59599516SKenneth E. Jansenc
429*59599516SKenneth E. Jansen            gradV(:,1,3) = gradV(:,1,3) + shg(:,n,3) * yl(:,n,2)
430*59599516SKenneth E. Jansen            gradV(:,2,3) = gradV(:,2,3) + shg(:,n,3) * yl(:,n,3)
431*59599516SKenneth E. Jansen            gradV(:,3,3) = gradV(:,3,3) + shg(:,n,3) * yl(:,n,4)
432*59599516SKenneth E. Jansenc                                             a j     u   a i
433*59599516SKenneth E. Jansenc from our notes where we had N_{a,j} = dN_a/dx_j  note that i is off by one because p was first in yl vector
434*59599516SKenneth E. Jansenc
435*59599516SKenneth E. Jansen         enddo
436*59599516SKenneth E. Jansen
437*59599516SKenneth E. Jansen         dwall = zero
438*59599516SKenneth E. Jansen         do n = 1, nenl
439*59599516SKenneth E. Jansen            dwall     = dwall + shape_funct(:,n) * dwl(:,n)
440*59599516SKenneth E. Jansen         enddo
441*59599516SKenneth E. Jansen
442*59599516SKenneth E. Jansen         kay(:)=zero
443*59599516SKenneth E. Jansen         epsilon(:)=zero
444*59599516SKenneth E. Jansen         do ii=1,npro
445*59599516SKenneth E. Jansen            do jj = 1, nshl
446*59599516SKenneth E. Jansen               kay(ii) =  kay(ii) + shape_funct(ii,jj) * yl(ii,jj,6)
447*59599516SKenneth E. Jansen               epsilon(ii) =  epsilon(ii)
448*59599516SKenneth E. Jansen     &              + shape_funct(ii,jj) * yl(ii,jj,7)
449*59599516SKenneth E. Jansen            enddo
450*59599516SKenneth E. Jansen         enddo
451*59599516SKenneth E. Jansenc        no source term in the LHS yet
452*59599516SKenneth E. Jansen         srcL(:)=zero
453*59599516SKenneth E. Jansen
454*59599516SKenneth E. Jansen         call elm3keps   (kay,     epsilon,
455*59599516SKenneth E. Jansen     &                    dwall,   gradV,
456*59599516SKenneth E. Jansen     &                    srcRat,  srcR,     srcJac)
457*59599516SKenneth E. Jansen         if(isclr.eq.1) srcL = srcJac(:,1)
458*59599516SKenneth E. Jansen         if(isclr.eq.2) srcL = srcJac(:,4)
459*59599516SKenneth E. Jansen         iadvdiff=0 ! scalar advection-diffusion flag
460*59599516SKenneth E. Jansen         if(iadvdiff.eq.1)then
461*59599516SKenneth E. Jansen            srcL(:)=zero
462*59599516SKenneth E. Jansen            srcR(:)=zero
463*59599516SKenneth E. Jansen            srcRat(:)=zero
464*59599516SKenneth E. Jansen         endif
465*59599516SKenneth E. Jansenc
466*59599516SKenneth E. Jansenc.... No source terms with the form (beta_i)*(phi,_i) for K or E
467*59599516SKenneth E. Jansenc
468*59599516SKenneth E. Jansen         uMod(:,1) = u1(:) - zero
469*59599516SKenneth E. Jansen         uMod(:,2) = u2(:) - zero
470*59599516SKenneth E. Jansen         uMod(:,3) = u3(:) - zero
471*59599516SKenneth E. Jansen
472*59599516SKenneth E. Jansen
473*59599516SKenneth E. Jansen      elseif (iLSet.ne.0) then
474*59599516SKenneth E. Jansen         if (isclr.eq.1)  then
475*59599516SKenneth E. Jansen
476*59599516SKenneth E. Jansen            srcR = zero
477*59599516SKenneth E. Jansen            srcL = zero
478*59599516SKenneth E. Jansen
479*59599516SKenneth E. Jansen         elseif (isclr.eq.2) then !we are redistancing level-sets
480*59599516SKenneth E. Jansen
481*59599516SKenneth E. JansenCAD   If Sclr(:,1).gt.zero, result of sign_term function 1
482*59599516SKenneth E. JansenCAD   If Sclr(:,1).eq.zero, result of sign_term function 0
483*59599516SKenneth E. JansenCAD   If Sclr(:,1).lt.zero, result of sign_term function -1
484*59599516SKenneth E. Jansen
485*59599516SKenneth E. Jansen            sclr_ls = zero      !zero out temp variable
486*59599516SKenneth E. Jansen
487*59599516SKenneth E. Jansen            do ii=1,npro
488*59599516SKenneth E. Jansen
489*59599516SKenneth E. Jansen               do jj = 1, nshl  ! first find the value of levelset at point ii
490*59599516SKenneth E. Jansen
491*59599516SKenneth E. Jansen                  sclr_ls(ii) =  sclr_ls(ii)
492*59599516SKenneth E. Jansen     &                        + shape_funct(ii,jj) * yl(ii,jj,6)
493*59599516SKenneth E. Jansen
494*59599516SKenneth E. Jansen               enddo
495*59599516SKenneth E. Jansen
496*59599516SKenneth E. Jansen               if (sclr_ls(ii) .gt. epsilon_ls)then
497*59599516SKenneth E. Jansen
498*59599516SKenneth E. Jansen                  sign_levelset(ii) = one
499*59599516SKenneth E. Jansen
500*59599516SKenneth E. Jansen               elseif  (abs(sclr_ls(ii)) .le. epsilon_ls)then
501*59599516SKenneth E. Jansen
502*59599516SKenneth E. Jansen                  sign_levelset(ii) = sclr_ls(ii)/epsilon_ls
503*59599516SKenneth E. Jansen     &                  + sin(pi*sclr_ls(ii)/epsilon_ls)/pi
504*59599516SKenneth E. Jansen
505*59599516SKenneth E. Jansen               elseif (sclr_ls(ii) .lt. epsilon_ls) then
506*59599516SKenneth E. Jansen
507*59599516SKenneth E. Jansen                  sign_levelset(ii) = - one
508*59599516SKenneth E. Jansen
509*59599516SKenneth E. Jansen               endif
510*59599516SKenneth E. Jansen
511*59599516SKenneth E. Jansen               srcR(ii) = sign_levelset(ii)
512*59599516SKenneth E. Jansen
513*59599516SKenneth E. Jansen            enddo
514*59599516SKenneth E. Jansen
515*59599516SKenneth E. Jansen            srcL =  zero
516*59599516SKenneth E. Jansen
517*59599516SKenneth E. Jansencad   The redistancing equation can be written in the following form
518*59599516SKenneth E. Jansencad
519*59599516SKenneth E. Jansencad   d_{,t} + sign(phi)*( d_{,i}/|d_{,i}| )* d_{,i} = sign(phi)
520*59599516SKenneth E. Jansencad
521*59599516SKenneth E. Jansencad   This is rewritten in the form
522*59599516SKenneth E. Jansencad
523*59599516SKenneth E. Jansencad   d_{,t} + u * d_{,i} = sign(phi)
524*59599516SKenneth E. Jansencad
525*59599516SKenneth E. Jansen
526*59599516SKenneth E. JansenCAD   For the redistancing equation the "velocity" term is calculated as
527*59599516SKenneth E. JansenCAD   follows
528*59599516SKenneth E. Jansen
529*59599516SKenneth E. Jansen
530*59599516SKenneth E. Jansen
531*59599516SKenneth E. Jansen            mytmp = sign_levelset / sqrt( gradS(:,1)*gradS(:,1)
532*59599516SKenneth E. Jansen     &                            + gradS(:,2)*gradS(:,2)
533*59599516SKenneth E. Jansen     &                            + gradS(:,3)*gradS(:,3))
534*59599516SKenneth E. Jansen
535*59599516SKenneth E. Jansen            uMod(:,1) = mytmp * gradS(:,1) ! These are for the LHS
536*59599516SKenneth E. Jansen            uMod(:,2) = mytmp * gradS(:,2)
537*59599516SKenneth E. Jansen            uMod(:,3) = mytmp * gradS(:,3)
538*59599516SKenneth E. Jansen
539*59599516SKenneth E. Jansen            u1 = umod(:,1)      ! These are for the RHS
540*59599516SKenneth E. Jansen            u2 = umod(:,2)
541*59599516SKenneth E. Jansen            u3 = umod(:,3)
542*59599516SKenneth E. Jansen
543*59599516SKenneth E. Jansen
544*59599516SKenneth E. Jansen         endif    ! close of scalar 2 of level set
545*59599516SKenneth E. Jansen
546*59599516SKenneth E. Jansen      else        ! NOT turbulence and NOT level set so this is a simple
547*59599516SKenneth E. Jansen                  ! scalar. If your scalar equation has a source term
548*59599516SKenneth E. Jansen                  ! then add your own like the above but leave an unforced case
549*59599516SKenneth E. Jansen                  ! as an option like you see here
550*59599516SKenneth E. Jansen
551*59599516SKenneth E. Jansen         srcR = zero
552*59599516SKenneth E. Jansen         srcL = zero
553*59599516SKenneth E. Jansen      endif
554*59599516SKenneth E. Jansen
555*59599516SKenneth E. Jansen      return
556*59599516SKenneth E. Jansen      end
557*59599516SKenneth E. Jansen
558*59599516SKenneth E. Jansen
559*59599516SKenneth E. Jansen
560