xref: /phasta/phSolver/compressible/e3source.f (revision 513954ef803c300cddba2bb96b4a5dac0b93489a)
159599516SKenneth E. Jansen        subroutine e3source  (ri,        rmi,          rlyi,
259599516SKenneth E. Jansen     &                        rho,       u1,           u2,
359599516SKenneth E. Jansen     &                        u3,        pres,         sforce,
459599516SKenneth E. Jansen     &                        dui,       dxidx,        ytargetl,
559599516SKenneth E. Jansen     &                        xl,        shpfun,       bcool)
659599516SKenneth E. Jansenc
759599516SKenneth E. Jansenc----------------------------------------------------------------------
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc This routine calculates the contribution of the bodyforce and surface
1059599516SKenneth E. Jansenc tension force operator to the RHS vector and LHS tangent matrix. The
1159599516SKenneth E. Jansenc temporary results are put in ri.
1259599516SKenneth E. Jansenc
1359599516SKenneth E. Jansenc  u1    (npro)               : x1-velocity component
1459599516SKenneth E. Jansenc  u2    (npro)               : x2-velocity component
1559599516SKenneth E. Jansenc  u3    (npro)               : x3-velocity component
1659599516SKenneth E. Jansenc  ri     (npro,nflow*(nsd+1)) : partial residual
1759599516SKenneth E. Jansenc  rmi    (npro,nflow*(nsd+1)) : partial modified residual
1859599516SKenneth E. Jansenc  rLyi  (npro,nflow)          : least-squares residual vector
1959599516SKenneth E. Jansenc  shape (npro,nshl)          : element shape functions
2059599516SKenneth E. Jansenc  g1yti  (npro)                : grad-Sclr in direction 1 at intpt
2159599516SKenneth E. Jansenc  g2yti  (npro)                : grad-Sclr in direction 2 at intpt
2259599516SKenneth E. Jansenc  g3yti  (npro)                : grad-Sclr in direction 3 at intpt
2359599516SKenneth E. Jansenc
2459599516SKenneth E. Jansen      use turbSA
2559599516SKenneth E. Jansen      use specialBC
2659599516SKenneth E. Jansen      include "common.h"
2759599516SKenneth E. Jansenc
2859599516SKenneth E. Jansen      dimension ri(npro,nflow*(nsd+1)),     rmi(npro,nflow*(nsd+1)),
2959599516SKenneth E. Jansen     &            u1(npro),                  u2(npro),
3059599516SKenneth E. Jansen     &            u3(npro),                  rho(npro),
3159599516SKenneth E. Jansen     &            pres(npro),
3259599516SKenneth E. Jansen     &            rLyi(npro,nflow),          sforce(npro,3),
3359599516SKenneth E. Jansen     &            shpfun(npro,nshl),
3459599516SKenneth E. Jansen     &            xl(npro,nenl,3),           xx(npro,3)
3559599516SKenneth E. Jansen
3659599516SKenneth E. Jansen      real*8 ytargeti(npro,nflow), ytargetl(npro,nshl,nflow)
3759599516SKenneth E. Jansen
3859599516SKenneth E. Jansen      real*8 src(npro,nflow), bcool(npro),
3959599516SKenneth E. Jansen     &       dui(npro,nflow), duitarg(npro,nflow),
4059599516SKenneth E. Jansen     &       dxidx( npro, nsd, nsd), xfind( npro ), delta(npro), rat
4159599516SKenneth E. Jansenc
4259599516SKenneth E. Jansenc......contribution of body force
4359599516SKenneth E. Jansenc
4459599516SKenneth E. Jansen      bcool=zero
4559599516SKenneth E. Jansen      src=zero
4659599516SKenneth E. Jansenc
4759599516SKenneth E. Jansen
4859599516SKenneth E. Jansen
4959599516SKenneth E. Jansen      if(matflg(5,1).eq.1) then ! usual case
5059599516SKenneth E. Jansen         src(:,1) = zero
5159599516SKenneth E. Jansen         src(:,2) = rho(:) * datmat(1,5,1)
5259599516SKenneth E. Jansen         src(:,3) = rho(:) * datmat(2,5,1)
5359599516SKenneth E. Jansen         src(:,4) = rho(:) * datmat(3,5,1)
5459599516SKenneth E. Jansen         src(:,5) = u1*src(:,2) + u2*src(:,3) + u3*src(:,4)
5559599516SKenneth E. Jansen      else if(matflg(5,1).eq.3) then ! user supplied white noise
5659599516SKenneth E. Jansen
5759599516SKenneth E. Jansen         xsor = 18
5859599516SKenneth E. Jansenc            ampl = spamp(lstep+1)
5959599516SKenneth E. Jansenc            rat = Delt(1)/0.1
6059599516SKenneth E. Jansen         ampl = 0.002*exp(-(0.1248222*(lstep)-2.9957323)**2)
6159599516SKenneth E. Jansenc            if((myrank.eq.zero).and.(intp.eq.ngauss)) write(*,*) ampl
6259599516SKenneth E. Jansen         delta(:) = 0.5*sqrt(dxidx(:,1,1)*dxidx(:,1,1) ! 1/dx
6359599516SKenneth E. Jansen     .            +dxidx(:,2,1)*dxidx(:,2,1)
6459599516SKenneth E. Jansen     .            +dxidx(:,3,1)*dxidx(:,3,1))
6559599516SKenneth E. Jansen         do i=1,npro
6659599516SKenneth E. Jansen            xfind(i) = (xsor-minval(xl(i,:,1)))
6759599516SKenneth E. Jansen     &               *(maxval(xl(i,:,1))-xsor)
6859599516SKenneth E. Jansen         enddo
6959599516SKenneth E. Jansen
7059599516SKenneth E. Jansen         where ( xfind .ge. 0. )
7159599516SKenneth E. Jansen            src(:,2) = rho(:) * ampl * delta
7259599516SKenneth E. Jansenc             scaling by element size is removed not to mess up
7359599516SKenneth E. Jansenc             refinement  studies
7459599516SKenneth E. Jansenc              src(:,2) = rho(:) * ampl
7559599516SKenneth E. Jansen            src(:,5) = u1*src(:,2)
7659599516SKenneth E. Jansen         endwhere
7759599516SKenneth E. Jansen
7859599516SKenneth E. Jansen      else if(matflg(5,1).ge.4) then ! cool case (sponge outside of a
7959599516SKenneth E. Jansen                                     ! revolved box defined from zinSponge to
8059599516SKenneth E. Jansen                                     ! zoutSponge in axial extent and 0
8159599516SKenneth E. Jansen                                     ! to radSponge in radial extent for
8259599516SKenneth E. Jansen                                     ! all theta)
8359599516SKenneth E. Jansen
8459599516SKenneth E. Jansenc           determine coordinates of quadrature pt
8559599516SKenneth E. Jansen            xx=zero
8659599516SKenneth E. Jansen            do n  = 1,nenl
8759599516SKenneth E. Jansen               xx(:,1) = xx(:,1)  + shpfun(:,n) * xl(:,n,1)
8859599516SKenneth E. Jansen               xx(:,2) = xx(:,2)  + shpfun(:,n) * xl(:,n,2)
8959599516SKenneth E. Jansen               xx(:,3) = xx(:,3)  + shpfun(:,n) * xl(:,n,3)
9059599516SKenneth E. Jansen            enddo
9159599516SKenneth E. Jansen            ytargeti=zero
9259599516SKenneth E. Jansen            do j=1,nflow
9359599516SKenneth E. Jansen               do n=1,nshl
9459599516SKenneth E. Jansen                  ytargeti(:,j) = ytargeti(:,j)
9559599516SKenneth E. Jansen     &                          + shpfun(:,n)*ytargetl(:,n,j)
9659599516SKenneth E. Jansen               enddo
9759599516SKenneth E. Jansen            enddo
9859599516SKenneth E. Jansen
9959599516SKenneth E. Jansen
10059599516SKenneth E. Jansenc            we=3.0*29./682.
10159599516SKenneth E. Jansen            rsteep=3.0
10259599516SKenneth E. Jansen            src=zero
10359599516SKenneth E. Jansen            radsts=radSponge*radSponge
10459599516SKenneth E. Jansen            CoefRatioI2O = grthISponge/grthOSponge
10559599516SKenneth E. Jansen            do id=1,npro
10659599516SKenneth E. Jansen               radsqr=xx(id,2)**2+xx(id,1)**2
10759599516SKenneth E. Jansen               if(xx(id,3).lt. zinSponge) then  ! map this into big outflow
10859599516SKenneth E. Jansen                                                ! sponge to keep logic
10959599516SKenneth E. Jansen                                                ! below simple
11059599516SKenneth E. Jansen
11159599516SKenneth E. Jansen                  xx(id,3)=(zinSponge-xx(id,3))*CoefRatioI2O
11259599516SKenneth E. Jansen     &                    + zoutSponge
11359599516SKenneth E. Jansen !
11459599516SKenneth E. Jansen !    CoefRatioI2O is the ratio of the inlet quadratic coefficient to the
11559599516SKenneth E. Jansen !    outlet quadratic coeficient (basically how much faster sponge
11659599516SKenneth E. Jansen !    coefficient grows in inlet region relative to outlet region)
11759599516SKenneth E. Jansen !
11859599516SKenneth E. Jansen               endif
11959599516SKenneth E. Jansen               if((xx(id,3).gt.zoutSponge).or.(radsqr.gt.radsts))  then
12059599516SKenneth E. Jansen                  rad=sqrt(radsqr)
12159599516SKenneth E. Jansen                  radc=max(rad,radSponge)
12259599516SKenneth E. Jansen                  zval=max(xx(id,3),zoutSponge)
12359599516SKenneth E. Jansen                  bcool(id)=grthOSponge*((zval-zoutSponge)**2
12459599516SKenneth E. Jansen     &                                   +(radc-radSponge)**2)
12559599516SKenneth E. Jansen                  bcool(id)=min(bcool(id),betamax)
12659599516SKenneth E. Jansenc     Determine the resulting density and energies
12759599516SKenneth E. Jansen               den   = ytargeti(id,1) / (Rgas * ytargeti(id,5))
12859599516SKenneth E. Jansen               ei    = ytargeti(id,5) * ( Rgas / gamma1 )
12959599516SKenneth E. Jansen               rk    = pt5 * ( ytargeti(id,2)**2+ytargeti(id,3)**2
13059599516SKenneth E. Jansen     &                                         +ytargeti(id,4)**2 )
13159599516SKenneth E. Jansenc     Determine the resulting conservation variables
13259599516SKenneth E. Jansen               duitarg(id,1) = den
13359599516SKenneth E. Jansen               duitarg(id,2) = den * ytargeti(id,2)
13459599516SKenneth E. Jansen               duitarg(id,3) = den * ytargeti(id,3)
13559599516SKenneth E. Jansen               duitarg(id,4) = den * ytargeti(id,4)
13659599516SKenneth E. Jansen               duitarg(id,5) = den * (ei + rk)
13759599516SKenneth E. Jansenc     Apply the sponge
13859599516SKenneth E. Jansen               if(spongeContinuity.eq.1)
13959599516SKenneth E. Jansen     &           src(id,1) = -bcool(id)*(dui(id,1) - duitarg(id,1))
14059599516SKenneth E. Jansen               if(spongeMomentum1.eq.1)
14159599516SKenneth E. Jansen     &           src(id,2) = -bcool(id)*(dui(id,2) - duitarg(id,2))
14259599516SKenneth E. Jansen               if(spongeMomentum2.eq.1)
14359599516SKenneth E. Jansen     &           src(id,3) = -bcool(id)*(dui(id,3) - duitarg(id,3))
14459599516SKenneth E. Jansen               if(spongeMomentum3.eq.1)
14559599516SKenneth E. Jansen     &           src(id,4) = -bcool(id)*(dui(id,4) - duitarg(id,4))
14659599516SKenneth E. Jansen               if(spongeEnergy.eq.1)
14759599516SKenneth E. Jansen     &           src(id,5) = -bcool(id)*(dui(id,5) - duitarg(id,5))
14859599516SKenneth E. Jansen            endif
14959599516SKenneth E. Jansen         enddo
15059599516SKenneth E. Jansen      else
15159599516SKenneth E. Jansen         if(isurf .ne. 1) then
15259599516SKenneth E. Jansen            write(*,*) 'only vector (1) and cooling (4) implemented'
15359599516SKenneth E. Jansen            stop
15459599516SKenneth E. Jansen         endif
15559599516SKenneth E. Jansen      endif
15659599516SKenneth E. Jansen
15759599516SKenneth E. Jansen      if (isurf .eq. 1) then    ! add the surface tension force
15859599516SKenneth E. Jansen         src(:,2) = src(:,2) +  rho(:)*sforce(:,1)
15959599516SKenneth E. Jansen         src(:,3) = src(:,3) +  rho(:)*sforce(:,2)
16059599516SKenneth E. Jansen         src(:,4) = src(:,4) +  rho(:)*sforce(:,3)
16159599516SKenneth E. Jansen         src(:,5) = src(:,5) + (u1*sforce(:,1)+u2*sforce(:,2)
16259599516SKenneth E. Jansen     &                       + u3*sforce(:,3))*rho(:)
16359599516SKenneth E. Jansen      endif
16459599516SKenneth E. Jansen
16559599516SKenneth E. Jansenc
16659599516SKenneth E. Jansenc==========================>>  IRES = 1 or 3  <<=======================
16759599516SKenneth E. Jansenc
16859599516SKenneth E. Jansen      if (ivart.gt.1) then
16959599516SKenneth E. Jansen         rLyi(:,1) = rLyi(:,1) - src(:,1)
17059599516SKenneth E. Jansen         rLyi(:,2) = rLyi(:,2) - src(:,2)
17159599516SKenneth E. Jansen         rLyi(:,3) = rLyi(:,3) - src(:,3)
17259599516SKenneth E. Jansen         rLyi(:,4) = rLyi(:,4) - src(:,4)
17359599516SKenneth E. Jansen         rLyi(:,5) = rLyi(:,5) - src(:,5)
17459599516SKenneth E. Jansen      endif
17559599516SKenneth E. Jansen
17659599516SKenneth E. Jansen      if ((ires .eq. 1) .or. (ires .eq. 3)) then ! we need ri built
17759599516SKenneth E. Jansen         ri (:,16) = ri (:,16) -  src(:,1)
17859599516SKenneth E. Jansen         ri (:,17) = ri (:,17) -  src(:,2)
17959599516SKenneth E. Jansen         ri (:,18) = ri (:,18) -  src(:,3)
18059599516SKenneth E. Jansen         ri (:,19) = ri (:,19) -  src(:,4)
18159599516SKenneth E. Jansen         ri (:,20) = ri (:,20) -  src(:,5)
18259599516SKenneth E. Jansen
18359599516SKenneth E. Jansen      endif
18459599516SKenneth E. Jansen
18559599516SKenneth E. Jansen      if ((ires.eq.2) .or. (ires.eq.3)) then ! we need rmi built
18659599516SKenneth E. Jansen         rmi (:,16) = rmi (:,16) -  src(:,1)
18759599516SKenneth E. Jansen         rmi (:,17) = rmi (:,17) -  src(:,2)
18859599516SKenneth E. Jansen         rmi (:,18) = rmi (:,18) -  src(:,3)
18959599516SKenneth E. Jansen         rmi (:,19) = rmi (:,19) -  src(:,4)
19059599516SKenneth E. Jansen         rmi (:,20) = rmi (:,20) -  src(:,5)
19159599516SKenneth E. Jansen      endif
19259599516SKenneth E. Jansenc
19359599516SKenneth E. Jansen      return
19459599516SKenneth E. Jansen      end
19559599516SKenneth E. Jansenc
19659599516SKenneth E. Jansenc
19759599516SKenneth E. Jansenc
19859599516SKenneth E. Jansen      subroutine e3sourceSclr(Sclr,    rho,    rmu,
199*513954efSKenneth E. Jansen     &                          dist2w,  vort,   gVnrm, con,
20059599516SKenneth E. Jansen     &                          g1yti,   g2yti,  g3yti,
20159599516SKenneth E. Jansen     &                          rti,     rLyti,  srcp,
20259599516SKenneth E. Jansen     &                          ycl,      shape,  u1,
20359599516SKenneth E. Jansen     &                          u2,      u3,      xl, elDwl)
20459599516SKenneth E. Jansenc
20559599516SKenneth E. Jansenc---------------------------------------------------------------------
20659599516SKenneth E. Jansenc
20759599516SKenneth E. Jansenc  This routine calculates the source term indicated in the Spalart-
20859599516SKenneth E. Jansenc  Allmaras eddy viscosity model.  After term is stored in rti(:,4),
20959599516SKenneth E. Jansenc  for later use by e3wmltSclr, and in rLyti(:) for later use by e3lsSclr.
21059599516SKenneth E. Jansenc
21159599516SKenneth E. Jansenc input:
21259599516SKenneth E. Jansenc  Sclr   (npro)              : working turbulence variable
21359599516SKenneth E. Jansenc  rho    (npro)              : density at intpt
21459599516SKenneth E. Jansenc  rmu    (npro)              : molecular viscosity
21559599516SKenneth E. Jansenc  dist2w (npro)              : distance from intpt to the nearest wall
21659599516SKenneth E. Jansenc  vort   (npro)              : magnitude of the vorticity
217*513954efSKenneth E. Jansenc  gVnrm  (npro)              : magnitude of the velocity gradient
21859599516SKenneth E. Jansenc  con    (npro)              : conductivity
21959599516SKenneth E. Jansenc  g1yti  (npro)              : grad-Sclr in direction 1
22059599516SKenneth E. Jansenc  g2yti  (npro)              : grad-Sclr in direction 2
22159599516SKenneth E. Jansenc  g3yti  (npro)              : grad-Sclr in direction 3
22259599516SKenneth E. Jansenc
22359599516SKenneth E. Jansenc output:
22459599516SKenneth E. Jansenc  rti    (npro,4)            : components of residual at intpt
22559599516SKenneth E. Jansenc  rLyti  (npro)              : GLS stabilization
22659599516SKenneth E. Jansenc
22759599516SKenneth E. Jansenc---------------------------------------------------------------------
22859599516SKenneth E. Jansenc
22959599516SKenneth E. Jansen      use turbSA
23059599516SKenneth E. Jansen      include "common.h"
23159599516SKenneth E. Jansenc
23259599516SKenneth E. Jansen      dimension Sclr   (npro),        ycl(npro,nshl,ndof),
23359599516SKenneth E. Jansen     &          dist2w (npro),          shape(npro,nshl),
234*513954efSKenneth E. Jansen     &          vort   (npro), gVnrm(npro), rho   (npro),
23559599516SKenneth E. Jansen     &          rmu    (npro),             con    (npro),
23659599516SKenneth E. Jansen     &          g1yti  (npro),             g2yti  (npro),
23759599516SKenneth E. Jansen     &          g3yti  (npro),             u1     (npro),
23859599516SKenneth E. Jansen     &          u2     (npro),             u3     (npro)
23959599516SKenneth E. Jansenc
24059599516SKenneth E. Jansen      dimension rti    (npro,4),           rLyti  (npro)
24159599516SKenneth E. Jansenc
24259599516SKenneth E. Jansen      dimension ft1    (npro),
24359599516SKenneth E. Jansenc unfix later -- pieces used in acusim:
24459599516SKenneth E. Jansen     &          srcrat (npro),             vdgn   (npro),
24559599516SKenneth E. Jansenc    &          term1  (npro),             term2  (npro),
24659599516SKenneth E. Jansenc    &          term3  (npro),
24759599516SKenneth E. Jansenc
24859599516SKenneth E. Jansen     &          chi    (npro),             fv1    (npro),
24959599516SKenneth E. Jansen     &          fv2    (npro),             Stilde (npro),
25059599516SKenneth E. Jansen     &          r      (npro),             g      (npro),
25159599516SKenneth E. Jansen     &          fw     (npro),             ft2    (npro),
25259599516SKenneth E. Jansen     &          fv1p   (npro),             fv2p   (npro),
25359599516SKenneth E. Jansen     &          stp    (npro),             rp     (npro),
25459599516SKenneth E. Jansen     &          gp     (npro),             fwp    (npro),
25559599516SKenneth E. Jansen     &          bf     (npro),             srcp   (npro),
25659599516SKenneth E. Jansen     &          gp6    (npro),             tmp    (npro),
25759599516SKenneth E. Jansen     &          tmp1   (npro),             fwog   (npro)
25859599516SKenneth E. Jansen      real*8 elDwl(npro) ! local quadrature point DES dvar
25959599516SKenneth E. Jansen      real*8 sclrm(npro) ! modified for non-negativity
260*513954efSKenneth E. Jansen      real*8 saCb1Scale(npro)  !Hack to change the production term and BL thickness
261*513954efSKenneth E. Jansen      real*8 xl_xbar(npro)     !Hack to store mean x location of element.
26259599516SKenneth E. Jansenc... for levelset
26359599516SKenneth E. Jansen      real*8  sign_levelset(npro), sclr_ls(npro), mytmp(npro),
26459599516SKenneth E. Jansen     &        xl(npro,nenl,nsd)
26559599516SKenneth E. Jansen
26659599516SKenneth E. Jansenc
26759599516SKenneth E. Jansen      if(iRANS.lt.0) then    ! spalart almaras model
26859599516SKenneth E. Jansen        sclrm=max(rmu/100.0,Sclr)
26959599516SKenneth E. Jansen        if(iles.lt.0) then
27059599516SKenneth E. Jansen          do i=1,npro
27159599516SKenneth E. Jansen            dx=maxval(xl(i,:,1))-minval(xl(i,:,1))
27259599516SKenneth E. Jansen            dy=maxval(xl(i,:,2))-minval(xl(i,:,2))
27359599516SKenneth E. Jansen            dz=maxval(xl(i,:,3))-minval(xl(i,:,3))
27459599516SKenneth E. Jansen            dmax=max(dx,max(dy,dz))
27559599516SKenneth E. Jansen            dmax=0.65d0*dmax
276*513954efSKenneth E. Jansen            if( iles.eq.-1) then !original DES97
27759599516SKenneth E. Jansen               dist2w(i)=min(dmax,dist2w(i))
278*513954efSKenneth E. Jansen            elseif(iles.eq.-2) then ! DDES
279*513954efSKenneth E. Jansen               rd=sclrm(i)*saKappaP2Inv/(dist2w(i)**2*gVnrm(i)+1.0d-12)
280*513954efSKenneth E. Jansen               fd=one-tanh((8.0000000000000000d0*rd)**3)
281*513954efSKenneth E. Jansen               dist2w(i)=dist2w(i)-fd*max(zero,dist2w(i)-dmax)
282*513954efSKenneth E. Jansen            endif
28359599516SKenneth E. Jansen          enddo
28459599516SKenneth E. Jansen        endif
28559599516SKenneth E. Jansen
28659599516SKenneth E. Jansen        elDwl(:)=elDwl(:)+dist2w(:)
28759599516SKenneth E. Jansenc
28859599516SKenneth E. Jansenc  determine chi
28959599516SKenneth E. Jansen        chi = rho*sclrm/rmu
29059599516SKenneth E. Jansenc  determine f_v1
29159599516SKenneth E. Jansen        fv1 = chi**3/(chi**3+saCv1**3)
29259599516SKenneth E. Jansenc  determine f_v2
29359599516SKenneth E. Jansen        fv2 = one - chi/(one+chi*fv1)
29459599516SKenneth E. Jansenc  determine Stilde
29559599516SKenneth E. Jansen        Stilde = vort + sclrm*fv2/(saKappa*dist2w)**2!unfix
29659599516SKenneth E. Jansenc  determine r
29759599516SKenneth E. Jansen        where(Stilde(:).ne.zero)
29859599516SKenneth E. Jansen           r(:) = sclrm(:)/Stilde(:)/(saKappa*dist2w(:))**2
29959599516SKenneth E. Jansen        elsewhere
30059599516SKenneth E. Jansen           r(:) = 1.0d32
30159599516SKenneth E. Jansen        endwhere
30259599516SKenneth E. Jansenc  determine g
30359599516SKenneth E. Jansen        saCw3l=saCw3
30459599516SKenneth E. Jansen        g = r + saCw2*(r**6-r)
30559599516SKenneth E. Jansen        sixth = 1.0/6.0
30659599516SKenneth E. Jansenc            gp      = rp * (tmp + 5 * saCw2 * rP5)
30759599516SKenneth E. Jansenc
30859599516SKenneth E. Jansenc            gP6     = (g * g * g) ** 2
30959599516SKenneth E. Jansenc            tmp     = 1 / (gP6 + (saCw3*saCw3*saCw3)**2)
31059599516SKenneth E. Jansenc            tmp1    = ( (1 + (saCw3*saCw3*saCw3)**2) * tmp ) ** sixth
31159599516SKenneth E. Jansenc            fw      = g * tmp1
31259599516SKenneth E. Jansenc            fwp     = gp * tmp1 * (saCw3*saCw3*saCw3)**2 * tmp
31359599516SKenneth E. Jansenc  determine f_w and f_w/g
31459599516SKenneth E. Jansen        fwog = ((one+saCw3**6)/(g**6+saCw3**6))**sixth
31559599516SKenneth E. Jansen        fw   = g*fwog
31659599516SKenneth E. Jansenc  determine f_t2
31759599516SKenneth E. Jansenc        ft2 = ct3*exp(-ct4*chi**2)
31859599516SKenneth E. Jansen        ft2 = zero
31959599516SKenneth E. Jansen
32059599516SKenneth E. Jansenc        srcrat=saCb1*(one-ft2)*Stilde*sclrm
32159599516SKenneth E. Jansenc     &      -(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w)**2
32259599516SKenneth E. Jansenc        srcrat=srcrat/sclrm
32359599516SKenneth E. Jansen
324*513954efSKenneth E. Jansen!----------------------------------------------------------------------------
325*513954efSKenneth E. Jansen!HACK: lower the EV production rate within a region to decrease BL thickness.
326*513954efSKenneth E. Jansen! Appear NM was not finished yet        if(scrScaleEnable) then
327*513954efSKenneth E. Jansen        if(one.eq.zero) then
328*513954efSKenneth E. Jansen          do i = 1,nenl !average the x-locations
329*513954efSKenneth E. Jansen            xl_xbar(:) = xl_xbar(:) + xl(:,i,1)
330*513954efSKenneth E. Jansen          enddo
331*513954efSKenneth E. Jansen          xl_xbar = xl_xbar/nenl
332*513954efSKenneth E. Jansen
333*513954efSKenneth E. Jansen          saCb1Scale = one
334*513954efSKenneth E. Jansen          where(xl_xbar < saCb1alterXmin .and. xl_xbar > saCb1alterXmax)
335*513954efSKenneth E. Jansen            saCb1Scale(:) = seCb1alter
336*513954efSKenneth E. Jansen          endwhere
337*513954efSKenneth E. Jansen
338*513954efSKenneth E. Jansen          srcrat = saCb1Scale*saCb1*(one-ft2)*Stilde
339*513954efSKenneth E. Jansen     &         -(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w/dist2w)
340*513954efSKenneth E. Jansen        else
34159599516SKenneth E. Jansen          srcrat=saCb1*(one-ft2)*Stilde
34259599516SKenneth E. Jansen     &         -(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w/dist2w)
343*513954efSKenneth E. Jansen        endif
344*513954efSKenneth E. Jansen
345*513954efSKenneth E. Jansen!Original:
346*513954efSKenneth E. Jansen!        srcrat=saCb1*(one-ft2)*Stilde
347*513954efSKenneth E. Jansen!     &       -(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w/dist2w)
348*513954efSKenneth E. Jansen!End Hack
349*513954efSKenneth E. Jansen!----------------------------------------------------------------------------
35059599516SKenneth E. Jansen
35159599516SKenneth E. Jansenc
35259599516SKenneth E. Jansenc        term1=saCb1*(one-ft2)*Stilde*sclrm
35359599516SKenneth E. Jansenc        term2=saCb2*saSigmaInv*(g1yti**2+g2yti**2+g3yti**2)
35459599516SKenneth E. Jansenc        term3=-(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w)**2
35559599516SKenneth E. Jansenc determine d()/d(sclrm)
35659599516SKenneth E. Jansen        fv1p = 3*(saCv1**3)*(chi**2)*rho
35759599516SKenneth E. Jansen          fv1p = fv1p/(rmu*(chi**3+saCv1**3)**2)
35859599516SKenneth E. Jansen        fv2p = (chi**2)*fv1p-(one/rmu)
35959599516SKenneth E. Jansen          fv2p = fv2p/(one+chi*fv1)**2
36059599516SKenneth E. Jansen        stp = fv2 + sclrm*fv2p
36159599516SKenneth E. Jansen          stp = stp/(saKappa*dist2w)**2
36259599516SKenneth E. Jansen        where(Stilde(:).ne.zero)
36359599516SKenneth E. Jansen             rp(:) = Stilde(:) - sclrm(:)*stp(:)
36459599516SKenneth E. Jansen             rp(:) = rp(:)/(saKappa*dist2w(:)*Stilde(:))**2
36559599516SKenneth E. Jansen        elsewhere
36659599516SKenneth E. Jansen             rp(:) = 1.0d32
36759599516SKenneth E. Jansen        endwhere
36859599516SKenneth E. Jansen        gp = one+saCw2*(6*r**5 - one)
36959599516SKenneth E. Jansen          gp = gp*rp
37059599516SKenneth E. Jansen        fwp = (saCw3**6)*fwog
37159599516SKenneth E. Jansen          fwp = fwp*gp/(g**6+saCw3**6)
37259599516SKenneth E. Jansenc  determine source term
37359599516SKenneth E. Jansen        bf = saCb2*saSigmaInv*(g1yti**2+g2yti**2+g3yti**2)
37459599516SKenneth E. Jansen     &      +saCb1*(one-ft2)*Stilde*sclrm
37559599516SKenneth E. Jansen     &      -(saCw1*fw - saCb1*ft2/saKappa**2)*(sclrm/dist2w)**2
37659599516SKenneth E. Jansen        bf = bf * rho
37759599516SKenneth E. Jansenc determine d(source)/d(sclrm)
37859599516SKenneth E. Jansen        srcp = rho*saCb1*(sclrm*stp+Stilde)
37959599516SKenneth E. Jansen     &        -rho*saCw1*(fwp*sclrm**2 + 2*sclrm*fw)/dist2w**2
38059599516SKenneth E. Jansen        do i=1, npro
38159599516SKenneth E. Jansen          if(srcp(i).le.zero .and. srcp(i).le.srcrat(i)) then
38259599516SKenneth E. Jansen            srcp(i)=srcp(i)
38359599516SKenneth E. Jansen          else if(srcrat(i).lt.zero) then
38459599516SKenneth E. Jansen            srcp(i)=srcrat(i)
38559599516SKenneth E. Jansen          else
38659599516SKenneth E. Jansen            srcp(i)=zero
38759599516SKenneth E. Jansen          endif
38859599516SKenneth E. Jansen        enddo
38959599516SKenneth E. Jansenc
39059599516SKenneth E. Jansenc==========================>>  IRES = 1 or 3  <<=======================
39159599516SKenneth E. Jansenc
39259599516SKenneth E. Jansenc          if ((ires .eq. 1) .or. (ires .eq. 3)) then
39359599516SKenneth E. Jansen             rti (:,4) = rti (:,4) -  bf(:)
39459599516SKenneth E. Jansenc          endif                 !ires
39559599516SKenneth E. Jansen
39659599516SKenneth E. Jansenc          rmti (:,4) = rmti (:,4) - bf(:)
39759599516SKenneth E. Jansen          rLyti(:) = rLyti(:) - bf(:)
39859599516SKenneth E. Jansenc
39959599516SKenneth E. Jansen       elseif (iLSet.ne.0) then
40059599516SKenneth E. Jansen          if (isclr.eq.1)  then
40159599516SKenneth E. Jansen             srcp = zero
40259599516SKenneth E. Jansen
40359599516SKenneth E. Jansen          elseif (isclr.eq.2) then !we are redistancing level-sets
40459599516SKenneth E. Jansen
40559599516SKenneth E. Jansen             sclr_ls = zero     !zero out temp variable
40659599516SKenneth E. Jansen
40759599516SKenneth E. Jansen             do ii=1,npro
40859599516SKenneth E. Jansen
40959599516SKenneth E. Jansen                do jj = 1, nshl ! first find the value of levelset at point ii
41059599516SKenneth E. Jansen
41159599516SKenneth E. Jansen                   sclr_ls(ii) =  sclr_ls(ii)
41259599516SKenneth E. Jansen     &                  + shape(ii,jj) * ycl(ii,jj,6)
41359599516SKenneth E. Jansen
41459599516SKenneth E. Jansen                enddo
41559599516SKenneth E. Jansen
41659599516SKenneth E. Jansen                if (sclr_ls(ii) .lt. - epsilon_ls)then
41759599516SKenneth E. Jansen
41859599516SKenneth E. Jansen                   sign_levelset(ii) = - one
41959599516SKenneth E. Jansen
42059599516SKenneth E. Jansen                elseif  (abs(sclr_ls(ii)) .le. epsilon_ls)then
42159599516SKenneth E. Jansenc     sign_levelset(ii) = zero
42259599516SKenneth E. Jansenc
42359599516SKenneth E. Jansen                   sign_levelset(ii) =sclr_ls(ii)/epsilon_ls
42459599516SKenneth E. Jansen     &                  + sin(pi*sclr_ls(ii)/epsilon_ls)/pi
42559599516SKenneth E. Jansen
42659599516SKenneth E. Jansen
42759599516SKenneth E. Jansen                elseif (sclr_ls(ii) .gt. epsilon_ls) then
42859599516SKenneth E. Jansen
42959599516SKenneth E. Jansen                   sign_levelset(ii) = one
43059599516SKenneth E. Jansen
43159599516SKenneth E. Jansen                endif
43259599516SKenneth E. Jansen                srcp(ii) = sign_levelset(ii)
43359599516SKenneth E. Jansen
43459599516SKenneth E. Jansen             enddo
43559599516SKenneth E. Jansenc
43659599516SKenneth E. Jansenc     ad   The redistancing equation can be written in the following form
43759599516SKenneth E. Jansenc     ad
43859599516SKenneth E. Jansenc     ad   d_{,t} + sign(phi)*( d_{,i}/|d_{,i}| )* d_{,i} = sign(phi)
43959599516SKenneth E. Jansenc     ad
44059599516SKenneth E. Jansenc     ad   This is rewritten in the form
44159599516SKenneth E. Jansenc     ad
44259599516SKenneth E. Jansenc     ad   d_{,t} + u * d_{,i} = sign(phi)
44359599516SKenneth E. Jansenc     ad
44459599516SKenneth E. Jansen
44559599516SKenneth E. Jansenc$$$  CAD   For the redistancing equation the "pseudo velocity" term is
44659599516SKenneth E. Jansenc$$$  CAD   calculated as follows
44759599516SKenneth E. Jansen
44859599516SKenneth E. Jansen
44959599516SKenneth E. Jansen
45059599516SKenneth E. Jansen             mytmp = srcp(:) / sqrt( g1yti(:) * g1yti(:)
45159599516SKenneth E. Jansen     &            + g2yti(:) * g2yti(:)
45259599516SKenneth E. Jansen     &            + g3yti(:) * g3yti(:) )
45359599516SKenneth E. Jansen
45459599516SKenneth E. Jansen             u1 = mytmp(:) * g1yti(:)
45559599516SKenneth E. Jansen             u2 = mytmp(:) * g2yti(:)
45659599516SKenneth E. Jansen             u3 = mytmp(:) * g3yti(:)
45759599516SKenneth E. Jansenc
45859599516SKenneth E. Jansenc==========================>>  IRES = 1 or 3  <<=======================
45959599516SKenneth E. Jansenc
46059599516SKenneth E. Jansenc     if ((ires .eq. 1) .or. (ires .eq. 3)) then
46159599516SKenneth E. Jansen             rti (:,4) = rti (:,4) - srcp(:)
46259599516SKenneth E. Jansenc     endif                 !ires
46359599516SKenneth E. Jansen
46459599516SKenneth E. Jansenc     rmti (:,4) = rmti (:,4) -  srcp(:)
46559599516SKenneth E. Jansen             rLyti(:) = rLyti(:) - srcp(:)
46659599516SKenneth E. Jansenc
46759599516SKenneth E. Jansen          endif                 ! close of scalar 2 of level set
46859599516SKenneth E. Jansen
46959599516SKenneth E. Jansen       else    ! NOT turbulence and NOT level set so this is a simple
47059599516SKenneth E. Jansen               ! scalar. If your scalar equation has a source term
47159599516SKenneth E. Jansen               ! then add your own like the above but leave an unforced case
47259599516SKenneth E. Jansen               ! as an option like you see here
47359599516SKenneth E. Jansen
47459599516SKenneth E. Jansen          srcp = zero
47559599516SKenneth E. Jansen       endif
47659599516SKenneth E. Jansen
47759599516SKenneth E. Jansen
47859599516SKenneth E. Jansenc
47959599516SKenneth E. Jansenc.... Return and end
48059599516SKenneth E. Jansenc
48159599516SKenneth E. Jansen        return
48259599516SKenneth E. Jansen        end
48359599516SKenneth E. Jansen
484