xref: /phasta/phSolver/compressible/getdiff.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine getDiff (T,      cp,     rho,    ycl,
2*59599516SKenneth E. Jansen     &                      rmu,    rlm,    rlm2mu, con, shp,
3*59599516SKenneth E. Jansen     &                      xmudmi, xl)
4*59599516SKenneth E. Jansen
5*59599516SKenneth E. Jansenc----------------------------------------------------------------------
6*59599516SKenneth E. Jansenc
7*59599516SKenneth E. Jansenc This routine calculates the fluid material properties.
8*59599516SKenneth E. Jansenc
9*59599516SKenneth E. Jansenc input:
10*59599516SKenneth E. Jansenc  T      (npro)          : temperature
11*59599516SKenneth E. Jansenc  cp     (npro)          : specific heat at constant pressure
12*59599516SKenneth E. Jansenc **************************************************************
13*59599516SKenneth E. Jansenc  rho    (npro)          : density
14*59599516SKenneth E. Jansenc  ycl    (npro,nshl,ndof): Y variables
15*59599516SKenneth E. Jansenc  shp    (npro,nshl)     : element shape-functions
16*59599516SKenneth E. Jansenc *************************************************************
17*59599516SKenneth E. Jansenc output:
18*59599516SKenneth E. Jansenc  rmu    (npro)        : Mu
19*59599516SKenneth E. Jansenc  rlm    (npro)        : Lambda
20*59599516SKenneth E. Jansenc  rlm2mu (npro)        : Lambda + 2 Mu
21*59599516SKenneth E. Jansenc  con    (npro)        : Conductivity
22*59599516SKenneth E. Jansenc
23*59599516SKenneth E. Jansenc Note: material type flags
24*59599516SKenneth E. Jansenc         matflg(2):
25*59599516SKenneth E. Jansenc          eq. 0, constant viscosity
26*59599516SKenneth E. Jansenc          eq. 1, generalized Sutherland viscosity
27*59599516SKenneth E. Jansenc         matflg(3):
28*59599516SKenneth E. Jansenc          eq. 0, Stokes approximation
29*59599516SKenneth E. Jansenc          eq. 1, shear proportional bulk viscosity
30*59599516SKenneth E. Jansenc
31*59599516SKenneth E. Jansenc
32*59599516SKenneth E. Jansenc Farzin Shakib, Winter 1987.
33*59599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
34*59599516SKenneth E. Jansenc----------------------------------------------------------------------
35*59599516SKenneth E. Jansenc
36*59599516SKenneth E. Jansen      use turbSA
37*59599516SKenneth E. Jansen      use pointer_data
38*59599516SKenneth E. Jansen      include "common.h"
39*59599516SKenneth E. Jansenc
40*59599516SKenneth E. Jansen      dimension T(npro),                   cp(npro),
41*59599516SKenneth E. Jansen     &     rho(npro),                 Sclr(npro),
42*59599516SKenneth E. Jansen     &     rmu(npro),                 rlm(npro),
43*59599516SKenneth E. Jansen     &     rlm2mu(npro),              con(npro),
44*59599516SKenneth E. Jansen     &     ycl(npro,nshl,ndof),        shp(npro,nshl),
45*59599516SKenneth E. Jansen     &     xmudmi(npro,ngauss),        xl(npro,nenl,nsd),
46*59599516SKenneth E. Jansen     &     xx(npro)
47*59599516SKenneth E. Jansenc
48*59599516SKenneth E. Jansen      dimension xmut(npro)
49*59599516SKenneth E. Jansen      real*8 prop_blend(npro),test_it(npro)
50*59599516SKenneth E. Jansen
51*59599516SKenneth E. Jansen      integer n, e
52*59599516SKenneth E. Jansen      integer wallmask(nshl)
53*59599516SKenneth E. Jansen      real*8  xki, xki3, fv1, evisc
54*59599516SKenneth E. Jansenc
55*59599516SKenneth E. Jansenc
56*59599516SKenneth E. Jansenc.... constant viscosity
57*59599516SKenneth E. Jansenc
58*59599516SKenneth E. Jansen      if (matflg(2,1) .eq. 0) then
59*59599516SKenneth E. Jansenc
60*59599516SKenneth E. Jansen         if (iLSet .ne. 0)then  !two fluid properties used in this model
61*59599516SKenneth E. Jansen            Sclr = zero
62*59599516SKenneth E. Jansen            isc=abs(iRANS)+6
63*59599516SKenneth E. Jansen            do n = 1, nshl
64*59599516SKenneth E. Jansen               Sclr = Sclr + shp(:,n) * ycl(:,n,isc)
65*59599516SKenneth E. Jansen            enddo
66*59599516SKenneth E. Jansen            test_it = 0.5*(one + Sclr/epsilon_ls +
67*59599516SKenneth E. Jansen     &           (sin(pi*Sclr/epsilon_ls))/pi )
68*59599516SKenneth E. Jansen
69*59599516SKenneth E. Jansen            prop_blend = max( min(test_it(:), one ), zero  )
70*59599516SKenneth E. Jansen            rmu = datmat(1,2,2) + (datmat(1,2,1)-datmat(1,2,2))
71*59599516SKenneth E. Jansen     &           *prop_blend
72*59599516SKenneth E. Jansen         elseif(irampViscOutlet.eq.1)then ! increase viscosity near outlet
73*59599516SKenneth E. Jansenc.............ramp rmu near outlet (for a NGC geometry)
74*59599516SKenneth E. Jansen            xx=zero
75*59599516SKenneth E. Jansen            do n=1,nenl
76*59599516SKenneth E. Jansen               xx(:)=xx(:) + shp(:,n) * xl(:,n,1)
77*59599516SKenneth E. Jansen            enddo
78*59599516SKenneth E. Jansen            fmax=10.0
79*59599516SKenneth E. Jansen            where(xx(:).le. 0.42) !healfway btwn AIP and exit
80*59599516SKenneth E. Jansen               rmu(:)=datmat(1,2,1)
81*59599516SKenneth E. Jansen            elsewhere(xx(:).ge. 0.75) !2/3 of the way to the exit
82*59599516SKenneth E. Jansen               rmu(:)=fmax*datmat(1,2,1)
83*59599516SKenneth E. Jansen            elsewhere
84*59599516SKenneth E. Jansen               rmu(:)= datmat(1,2,1)*(
85*59599516SKenneth E. Jansen     &          (55.65294821-55.65294821*fmax)*xx(:)*xx(:)*xx(:)
86*59599516SKenneth E. Jansen     &          +(-97.67092412+97.67092412*fmax)*xx(:)*xx(:)
87*59599516SKenneth E. Jansen     &          +(52.59203606-52.59203606*fmax)*xx(:)
88*59599516SKenneth E. Jansen     &          -7.982719760+8.982719760*fmax)
89*59599516SKenneth E. Jansen            endwhere
90*59599516SKenneth E. Jansen         else ! constant viscosity
91*59599516SKenneth E. Jansen            rmu = datmat(1,2,1)
92*59599516SKenneth E. Jansen         endif
93*59599516SKenneth E. Jansenc
94*59599516SKenneth E. Jansen      else
95*59599516SKenneth E. Jansenc
96*59599516SKenneth E. Jansenc.... generalized Sutherland viscosity
97*59599516SKenneth E. Jansenc
98*59599516SKenneth E. Jansen         rmu = datmat(1,2,1) * (T/datmat(2,2,1))*sqrt(T/datmat(2,2,1))
99*59599516SKenneth E. Jansen     &        * ( datmat(2,2,1) + datmat(3,2,1) ) / (T + datmat(3,2,1))
100*59599516SKenneth E. Jansenc
101*59599516SKenneth E. Jansen      endif
102*59599516SKenneth E. Jansenc
103*59599516SKenneth E. Jansenc.... calculate the second viscosity coefficient
104*59599516SKenneth E. Jansenc
105*59599516SKenneth E. Jansen      if (matflg(3,1) .eq. 0) then
106*59599516SKenneth E. Jansen         rlm = -pt66 * rmu
107*59599516SKenneth E. Jansen      else
108*59599516SKenneth E. Jansen         rlm = (datmat(1,3,1) - pt66) * rmu
109*59599516SKenneth E. Jansen      endif
110*59599516SKenneth E. Jansenc
111*59599516SKenneth E. Jansenc.... calculate the remaining quantities
112*59599516SKenneth E. Jansenc
113*59599516SKenneth E. Jansen      con    = rmu * cp / pr
114*59599516SKenneth E. Jansenc
115*59599516SKenneth E. Jansenc-------------Eddy Viscosity Calculation-----------------
116*59599516SKenneth E. Jansenc
117*59599516SKenneth E. Jansenc.... dynamic model
118*59599516SKenneth E. Jansenc
119*59599516SKenneth E. Jansen      if (iLES .gt. 0. and. iRANS.eq.0) then  ! simple LES
120*59599516SKenneth E. Jansen         xmut = xmudmi(:,intp)
121*59599516SKenneth E. Jansen      else if (iRANS .eq. 0 .and. iLES.eq.0 ) then !DNS
122*59599516SKenneth E. Jansen         xmut = zero
123*59599516SKenneth E. Jansen      else if (iRANS .lt. 0) then ! calculate RANS viscosity
124*59599516SKenneth E. Jansenc
125*59599516SKenneth E. Jansenc.... RANS
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansen         do e = 1, npro
128*59599516SKenneth E. Jansen            wallmask = 0
129*59599516SKenneth E. Jansen            if(itwmod.eq.-2) then ! effective viscosity
130*59599516SKenneth E. Jansenc mark the wall nodes for this element, if there are any
131*59599516SKenneth E. Jansen               do n = 1, nshl
132*59599516SKenneth E. Jansenc
133*59599516SKenneth E. Jansenc  note that we are using ycl here so that means that these
134*59599516SKenneth E. Jansenc  terms are not perturbed for MFG difference and therefore
135*59599516SKenneth E. Jansenc  NOT in the LHS.  As they only give the evisc near the wall
136*59599516SKenneth E. Jansenc  I doubt this is a problem.
137*59599516SKenneth E. Jansenc
138*59599516SKenneth E. Jansen                  u1=ycl(e,n,2)
139*59599516SKenneth E. Jansen                  u2=ycl(e,n,3)
140*59599516SKenneth E. Jansen                  u3=ycl(e,n,4)
141*59599516SKenneth E. Jansen                  if((u1.eq.zero).and.(u2.eq.zero).and.(u3.eq.zero))
142*59599516SKenneth E. Jansen     &                 then
143*59599516SKenneth E. Jansen                     wallmask(n)=1
144*59599516SKenneth E. Jansen                  endif
145*59599516SKenneth E. Jansen               enddo
146*59599516SKenneth E. Jansen            endif
147*59599516SKenneth E. Jansenc
148*59599516SKenneth E. Jansen            if( any(wallmask.eq.1) ) then
149*59599516SKenneth E. Jansenc if there are wall nodes for this elt in an effective-viscosity wall
150*59599516SKenneth E. Jansenc modeled case,then eddy viscosity has been stored at the wall nodes
151*59599516SKenneth E. Jansenc in place of the spalart-allmaras variable; the eddy viscosity for
152*59599516SKenneth E. Jansenc the whole element is taken to be the avg of wall values
153*59599516SKenneth E. Jansen               evisc = zero
154*59599516SKenneth E. Jansen               nwnode=0
155*59599516SKenneth E. Jansen               do n = 1, nshl
156*59599516SKenneth E. Jansen                  if(wallmask(n).eq.1) then
157*59599516SKenneth E. Jansen                     evisc = evisc + ycl(e,n,6)
158*59599516SKenneth E. Jansen                     nwnode = nwnode + 1
159*59599516SKenneth E. Jansen                  endif
160*59599516SKenneth E. Jansen               enddo
161*59599516SKenneth E. Jansen               evisc = evisc/nwnode
162*59599516SKenneth E. Jansen               xmut(e)= abs(evisc)
163*59599516SKenneth E. Jansenc this is what we would use instead of the above if we were allowing
164*59599516SKenneth E. Jansenc the eddy viscosity to vary through the element based on non-wall nodes
165*59599516SKenneth E. Jansenc$$$               evisc = zero
166*59599516SKenneth E. Jansenc$$$               Turb = zero
167*59599516SKenneth E. Jansenc$$$               do n = 1, nshl
168*59599516SKenneth E. Jansenc$$$                  if(wallmask(n).eq.1) then
169*59599516SKenneth E. Jansenc$$$                     evisc = evisc + shape(e,n) * ycl(e,n,6)
170*59599516SKenneth E. Jansenc$$$                  else
171*59599516SKenneth E. Jansenc$$$                     Turb = Turb + shape(e,n) * ycl(e,n,6)
172*59599516SKenneth E. Jansenc$$$                  endif
173*59599516SKenneth E. Jansenc$$$               enddo
174*59599516SKenneth E. Jansenc$$$               xki    = abs(Turb)/rmu(e)
175*59599516SKenneth E. Jansenc$$$               xki3   = xki * xki * xki
176*59599516SKenneth E. Jansenc$$$               fv1    = xki3 / (xki3 + saCv1P3)
177*59599516SKenneth E. Jansenc$$$               rmu(e) = rmu(e) + fv1*abs(Turb)
178*59599516SKenneth E. Jansenc$$$               rmu(e) = rmu(e) + abs(evisc)
179*59599516SKenneth E. Jansen            else
180*59599516SKenneth E. Jansenc else one of the following is the case:
181*59599516SKenneth E. Jansenc   using effective-viscosity, but no wall nodes on this elt
182*59599516SKenneth E. Jansenc   using slip-velocity
183*59599516SKenneth E. Jansenc   using no model; walls are resolved
184*59599516SKenneth E. Jansenc in all of these cases, eddy viscosity is calculated normally
185*59599516SKenneth E. Jansen               savar = zero
186*59599516SKenneth E. Jansen               do n = 1, nshl
187*59599516SKenneth E. Jansen                  savar = savar + shp(e,n) * ycl(e,n,6)
188*59599516SKenneth E. Jansen               enddo
189*59599516SKenneth E. Jansen               xki    = abs(savar)/rmu(e)
190*59599516SKenneth E. Jansen               xki3   = xki * xki * xki
191*59599516SKenneth E. Jansen               fv1    = xki3 / (xki3 + saCv1P3)
192*59599516SKenneth E. Jansen               xmut(e) = fv1*abs(savar)
193*59599516SKenneth E. Jansen            endif
194*59599516SKenneth E. Jansen         enddo                  ! end loop over elts
195*59599516SKenneth E. Jansen
196*59599516SKenneth E. Jansen         if (iLES.gt.0) then    ! this is DES so we have to blend in
197*59599516SKenneth E. Jansen                                ! xmudmi based on max edge length of
198*59599516SKenneth E. Jansen                                ! element
199*59599516SKenneth E. Jansen            call EviscDES (xl,xmut,xmudmi)
200*59599516SKenneth E. Jansen         endif
201*59599516SKenneth E. Jansen      endif                     ! check for LES or RANS
202*59599516SKenneth E. Jansen
203*59599516SKenneth E. Jansen      rlm    = rlm - pt66*xmuT
204*59599516SKenneth E. Jansen      rmu    = rmu + xmuT
205*59599516SKenneth E. Jansen      rlm2mu = rlm + two * rmu
206*59599516SKenneth E. Jansen      con    = con + xmuT*cp/pr
207*59599516SKenneth E. Jansenc
208*59599516SKenneth E. Jansenc.... return
209*59599516SKenneth E. Jansenc
210*59599516SKenneth E. Jansen      return
211*59599516SKenneth E. Jansen      end
212*59599516SKenneth E. Jansenc
213*59599516SKenneth E. Jansenc
214*59599516SKenneth E. Jansenc
215*59599516SKenneth E. Jansen      subroutine getDiffSclr (T,      cp,   rmu,   rlm,
216*59599516SKenneth E. Jansen     &     rlm2mu, con, rho, Sclr)
217*59599516SKenneth E. Jansenc
218*59599516SKenneth E. Jansenc----------------------------------------------------------------------
219*59599516SKenneth E. Jansenc
220*59599516SKenneth E. Jansenc This routine calculates the fluid material properties.
221*59599516SKenneth E. Jansenc
222*59599516SKenneth E. Jansenc input:
223*59599516SKenneth E. Jansenc  T      (npro)        : temperature
224*59599516SKenneth E. Jansenc  cp     (npro)        : specific heat at constant pressure
225*59599516SKenneth E. Jansenc
226*59599516SKenneth E. Jansenc output:
227*59599516SKenneth E. Jansenc  rmu    (npro)        : Mu
228*59599516SKenneth E. Jansenc  rlm    (npro)        : Lambda
229*59599516SKenneth E. Jansenc  rlm2mu (npro)        : Lambda + 2 Mu
230*59599516SKenneth E. Jansenc  con    (npro)        : Conductivity
231*59599516SKenneth E. Jansenc
232*59599516SKenneth E. Jansenc Note: material type flags
233*59599516SKenneth E. Jansenc         matflg(2):
234*59599516SKenneth E. Jansenc          eq. 0, constant viscosity
235*59599516SKenneth E. Jansenc          eq. 1, generalized Sutherland viscosity
236*59599516SKenneth E. Jansenc         matflg(3):
237*59599516SKenneth E. Jansenc          eq. 0, Stokes approximation
238*59599516SKenneth E. Jansenc          eq. 1, shear proportional bulk viscosity
239*59599516SKenneth E. Jansenc
240*59599516SKenneth E. Jansenc
241*59599516SKenneth E. Jansenc Farzin Shakib, Winter 1987.
242*59599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
243*59599516SKenneth E. Jansenc----------------------------------------------------------------------
244*59599516SKenneth E. Jansenc
245*59599516SKenneth E. Jansen        include "common.h"
246*59599516SKenneth E. Jansenc
247*59599516SKenneth E. Jansen        dimension T(npro),                   cp(npro),
248*59599516SKenneth E. Jansen     &            rmu(npro),                 rlm(npro),
249*59599516SKenneth E. Jansen     &            rlm2mu(npro),              con(npro),
250*59599516SKenneth E. Jansen     &            rho(npro),                 Sclr(npro)
251*59599516SKenneth E. Jansen
252*59599516SKenneth E. Jansen
253*59599516SKenneth E. Jansen
254*59599516SKenneth E. Jansenc
255*59599516SKenneth E. Jansenc
256*59599516SKenneth E. Jansenc.... constant viscosity
257*59599516SKenneth E. Jansenc
258*59599516SKenneth E. Jansen        if (matflg(2,1) .eq. 0) then
259*59599516SKenneth E. Jansenc
260*59599516SKenneth E. Jansen          rmu = datmat(1,2,1)
261*59599516SKenneth E. Jansenc
262*59599516SKenneth E. Jansen        else
263*59599516SKenneth E. Jansenc
264*59599516SKenneth E. Jansenc.... generalized Sutherland viscosity
265*59599516SKenneth E. Jansenc
266*59599516SKenneth E. Jansen          rmu = datmat(1,2,1) * (T/datmat(2,2,1))*sqrt(T/datmat(2,2,1))
267*59599516SKenneth E. Jansen     &        * ( datmat(2,2,1) + datmat(3,2,1) ) / (T + datmat(3,2,1))
268*59599516SKenneth E. Jansenc
269*59599516SKenneth E. Jansen        endif
270*59599516SKenneth E. Jansenc
271*59599516SKenneth E. Jansen*************************check****************************
272*59599516SKenneth E. Jansenc        if (iRANS(1).lt.zero) then
273*59599516SKenneth E. Jansenc           rmu = saSigmaInv*rho*((rmu/rho)+Sclr)
274*59599516SKenneth E. Jansenc        endif
275*59599516SKenneth E. Jansenc This augmentation of viscosity is performed in e3viscsclr
276*59599516SKenneth E. Jansenc The Spalart -Allmaras model will need molecular viscosity
277*59599516SKenneth E. Jansenc  in subsequent calculations.
278*59599516SKenneth E. Jansenc.... calculate the second viscosity coefficient
279*59599516SKenneth E. Jansenc
280*59599516SKenneth E. Jansen        if (matflg(3,1) .eq. 0) then
281*59599516SKenneth E. Jansenc
282*59599516SKenneth E. Jansen          rlm = -pt66 * rmu
283*59599516SKenneth E. Jansenc
284*59599516SKenneth E. Jansen        else
285*59599516SKenneth E. Jansenc
286*59599516SKenneth E. Jansen          rlm = (datmat(1,3,1) - pt66) * rmu
287*59599516SKenneth E. Jansenc
288*59599516SKenneth E. Jansen        endif
289*59599516SKenneth E. Jansenc
290*59599516SKenneth E. Jansenc.... calculate the remaining quantities
291*59599516SKenneth E. Jansenc
292*59599516SKenneth E. Jansen
293*59599516SKenneth E. Jansen
294*59599516SKenneth E. Jansen
295*59599516SKenneth E. Jansen        rlm2mu = rlm + two * rmu
296*59599516SKenneth E. Jansen        con    = rmu * cp / pr
297*59599516SKenneth E. Jansen
298*59599516SKenneth E. Jansen
299*59599516SKenneth E. Jansen
300*59599516SKenneth E. Jansenc
301*59599516SKenneth E. Jansenc.... return
302*59599516SKenneth E. Jansenc
303*59599516SKenneth E. Jansen        return
304*59599516SKenneth E. Jansen        end
305*59599516SKenneth E. Jansen
306*59599516SKenneth E. Jansen      subroutine EviscDES(xl,xmut,xmudmi)
307*59599516SKenneth E. Jansen
308*59599516SKenneth E. Jansen      include "common.h"
309*59599516SKenneth E. Jansen      real*8 xmut(npro),xl(npro,nenl,nsd),xmudmi(npro,ngauss)
310*59599516SKenneth E. Jansen
311*59599516SKenneth E. Jansen
312*59599516SKenneth E. Jansen      do i=1,npro
313*59599516SKenneth E. Jansen         dx=maxval(xl(i,:,1))-minval(xl(i,:,1))
314*59599516SKenneth E. Jansen         dy=maxval(xl(i,:,2))-minval(xl(i,:,2))
315*59599516SKenneth E. Jansen         dz=maxval(xl(i,:,3))-minval(xl(i,:,3))
316*59599516SKenneth E. Jansen         emax=max(dx,max(dy,dz))
317*59599516SKenneth E. Jansen         if(emax.lt.eles) then  ! pure les
318*59599516SKenneth E. Jansen            xmut(i)=xmudmi(i,intp)
319*59599516SKenneth E. Jansen         else if(emax.lt.two*eles) then ! blend
320*59599516SKenneth E. Jansen            xi=(emax-eles)/(eles)
321*59599516SKenneth E. Jansen            xmut(i)=xi*xmut(i)+(one-xi)*xmudmi(1,intp)
322*59599516SKenneth E. Jansen         endif                  ! leave at RANS value as edge is twice pure les
323*59599516SKenneth E. Jansen      enddo
324*59599516SKenneth E. Jansen
325*59599516SKenneth E. Jansen      return
326*59599516SKenneth E. Jansen      end
327