vtf-logo

src/3d/equations/euler/rprhok/rpn3eurhokefix.f

c
c     =========================================================
      subroutine rpn3eurhok(ixyz,maxm,meqn,mwaves,mbc,mx,ql,qr,maux,
     &                      auxl,auxr,wave,s,amdq,apdq)
c     =========================================================
c
c     # solve Riemann problems for the thermally perfect 1D multi-component 
c     # Euler equations using a Roe-type approximate Riemann solver.  
c     # Scheme is blended with HLL for robustness and uses a multi-dimensional 
c     # entropy correction to prevent the carbuncle phenomenon. 
c
c     # On input, ql contains the state vector at the left edge of each cell
c     #           qr contains the state vector at the right edge of each cell
c     # This data is along a slice in the x-direction if ixyz=1
c     #                               the y-direction if ixyz=2.
c     #                               the z-direction if ixyz=3.
c
c     # On output, wave contains the waves, s the speeds, 
c     # amdq and apdq the positive and negative flux.
c
c     # Note that the i'th Riemann problem has left state qr(i-1,:)
c     #                                    and right state ql(i,:)
c     # From the basic clawpack routine step1, rp is called with ql = qr = q.
c
c     # Copyright (C) 2002 Ralf Deiterding
c     # Brandenburgische Universitaet Cottbus
c
      implicit double precision (a-h,o-z)
c
      include "ck.i"
c
      dimension   ql(1-mbc:maxm+mbc, meqn)
      dimension   qr(1-mbc:maxm+mbc, meqn)
      dimension    s(1-mbc:maxm+mbc, mwaves)
      dimension wave(1-mbc:maxm+mbc, meqn, mwaves)
      dimension amdq(1-mbc:maxm+mbc, meqn)
      dimension apdq(1-mbc:maxm+mbc, meqn)
      dimension auxl(1-mbc:maxm+mbc, maux, 3)
      dimension auxr(1-mbc:maxm+mbc, maux, 3)
c
c
c     local arrays -- common block comroe is passed to rpt3eurhok
c     ------------
      parameter (maxmrp = 1005) !# assumes atmost max(mx,my,mz) = 1000 with mbc=5
      parameter (minmrp = -4)   !# assumes at most mbc=5
      common /comroe/ u(minmrp:maxmrp), v(minmrp:maxmrp), 
     &     w(minmrp:maxmrp), u2v2w2(minmrp:maxmrp), 
     &     enth(minmrp:maxmrp), a(minmrp:maxmrp), 
     &     g1a2(minmrp:maxmrp), dpY(minmrp:maxmrp), 
     &     Y(LeNsp,minmrp:maxmrp), pk(LeNsp,-1:maxmrp) 
      logical efix, pfix, hll, roe, hllfix
c
c     define local arrays
c
      dimension delta(LeNsp+4)
      dimension rkl(LeNsp), rkr(LeNsp)
      dimension hkl(LeNsp), hkr(LeNsp)
      dimension fl(minmrp:maxmrp,LeNsp+5), fr(minmrp:maxmrp,LeNsp+5)
c
      data efix /.true./   !# use entropy fix
      data pfix /.true./   !# use Larrouturou's positivity fix for species
      data hll  /.true./   !# use HLL instead of Roe solver, if unphysical values occur
      data roe  /.true./   !# turn off Roe solver when debugging HLL
c
c     # Riemann solver returns fluxes
c     ------------
      common /rpnflx/ mrpnflx
      mrpnflx = 0
c
      if (minmrp.gt.1-mbc .or. maxmrp .lt. maxm+mbc) then
         write(6,*) 'need to increase maxm2 in rpA'
         stop
      endif
c
c     # set mu to point to  the component of the system that corresponds
c     # to momentum in the direction of this slice, mv and mw to the 
c     # orthogonal momentums:
c
      if(ixyz .eq. 1)then
         mu = Nsp+1
         mv = Nsp+2
         mw = Nsp+3
      else if(ixyz .eq. 2)then
         mu = Nsp+2
         mv = Nsp+3
         mw = Nsp+1
      else
         mu = Nsp+3
         mv = Nsp+1
         mw = Nsp+2
      endif
      mE = Nsp+4
      mT = Nsp+5
c
c     # note that notation for u,v, and w reflects assumption that the 
c     # Riemann problems are in the x-direction with u in the normal
c     # direction and v and w in the orthogonal directions, but with the 
c     # above definitions of mu, mv, and mw the routine also works with 
c     # ixyz=2 and ixyz = 3
c     # and returns, for example, f0 as the Godunov flux g0 for the
c     # Riemann problems u_t + g(u)_y = 0 in the y-direction.
c
c     # compute the Roe-averaged variables needed in the Roe solver.
c     # These are stored in the common block comroe since they are
c     # later used in routine rpt2eurhok to do the transverse wave splitting.
c
      do 20 i=2-mbc,mx+mbc
         rhol = 0.d0
         rhor = 0.d0
         do k = 1, Nsp
            rkl(k) = qr(i-1,k)
            rkr(k) = ql(i  ,k)
            rhol = rhol + rkl(k)
            rhor = rhor + rkr(k)
         enddo
         if( rhol.le.1.d-10 ) then
            write(6,*) 'negative total density, left', rhol
            stop
         endif
         if( rhor.le.1.d-10 ) then
            write(6,*) 'negative total density, right', rhor
            stop
         endif
c
c        # compute left/right rho/W and rho*Cp
c     
         rhoWl = 0.d0
         rhoWr = 0.d0
         do k = 1, Nsp
            rhoWl = rhoWl + rkl(k)/Wk(k)
            rhoWr = rhoWr + rkr(k)/Wk(k)
         enddo
c
c        # left/right Temperatures already calculated
c
         rhoel = qr(i-1,mE)-0.5d0*
     &          (qr(i-1,mu)**2+qr(i-1,mv)**2+qr(i-1,mw)**2)/rhol
         call SolveTrhok(qr(i-1,mT),rhoel,rhoWl,rkl,Nsp,ifail) 
         rhoer = ql(i  ,mE)-0.5d0*
     &          (ql(i  ,mu)**2+ql(i  ,mv)**2+ql(i  ,mw)**2)/rhor
         call SolveTrhok(ql(i  ,mT),rhoer,rhoWr,rkr,Nsp,ifail) 
c
         Tl = qr(i-1,mT)
         Tr = ql(i  ,mT)
         pl = rhoWl*RU*Tl 
         pr = rhoWr*RU*Tr
c
c       # compute quantities for rho-average
c
         rhsqrtl = dsqrt(rhol)  
         rhsqrtr = dsqrt(rhor)
         rhsq2 = rhsqrtl + rhsqrtr
c
c        # find rho-averaged specific velocity and enthalpy
c
         u(i) = (qr(i-1,mu)/rhsqrtl + ql(i,mu)/rhsqrtr) / rhsq2
         v(i) = (qr(i-1,mv)/rhsqrtl + ql(i,mv)/rhsqrtr) / rhsq2
         w(i) = (qr(i-1,mw)/rhsqrtl + ql(i,mw)/rhsqrtr) / rhsq2
         u2v2w2(i) = u(i)**2 + v(i)**2 + w(i)**2
         enth(i) = (((qr(i-1,mE)+pl)/rhsqrtl
     &             + (ql(i  ,mE)+pr)/rhsqrtr)) / rhsq2  
c
c        # compute rho-averages for T, cp, and W
c
         T  = (Tl * rhsqrtl + Tr * rhsqrtr) / rhsq2
         Wm = rhsq2 / (rhoWl/rhsqrtl + rhoWr/rhsqrtr) 
c        
c        # evaluate left/right entropies and mean cp
c
         call tabintp( Tl, hkl, hms, Nsp )
         call tabintp( Tr, hkr, hms, Nsp )
         do k = 1, Nsp
            Y(k,i) = (rkl(k)/rhsqrtl + rkr(k)/rhsqrtr) / rhsq2
         enddo
c         
         Cp = Cpmix( Tl, Tr, hkl, hkr, Y(1,i) )
         gamma1 = RU / ( Wm*Cp - RU )
         gamma  = gamma1 + 1.d0
c
c        # find rho-averaged specific enthalpies,
c        # compute rho-averaged mass fractions and
c        # compute partial pressure derivatives
c
         tmp = gamma * RU * T / gamma1 
         ht = 0.d0
         do k = 1, Nsp
            hk     = (hkl(k)*rhsqrtl + hkr(k)*rhsqrtr) / rhsq2
            pk(k,i) = 0.5d0*u2v2w2(i) - hk + tmp / Wk(k)
         enddo
c
c        # compute speed of sound
c
         dpY(i) = 0.d0
         do k = 1, Nsp
            dpY(i) = dpY(i) + pk(k,i) * Y(k,i)
         enddo
         a2 = dpY(i) + enth(i)-u2v2w2(i)
         g1a2(i) = 1.d0 / a2
         a(i) = dsqrt(gamma1*a2) 
c     
   20 continue
c
c
      do 30 i=2-mbc,mx+mbc
c
c        # find a1 thru a5, the coefficients of the Nsp+4 eigenvectors:
c
         dpdr = 0.d0
         drho = 0.d0
         do k = 1, Nsp
            delta(k) = ql(i,k) - qr(i-1,k)
            drho = drho + delta(k)
            dpdr = dpdr + pk(k,i) * delta(k)
         enddo
         delta(mu) = ql(i,mu) - qr(i-1,mu)
         delta(mv) = ql(i,mv) - qr(i-1,mv)
         delta(mw) = ql(i,mw) - qr(i-1,mw)
         delta(mE) = ql(i,mE) - qr(i-1,mE)
c
         a2 = g1a2(i)*(dpdr - ( u(i)*delta(mu) + v(i)*delta(mv) + 
     &        w(i)*delta(mw) ) + delta(mE) )
         a3 = delta(mv) - v(i)*drho
         a4 = delta(mw) - w(i)*drho
         a5 = 0.5d0*( a2 - ( u(i)*drho - delta(mu) )/a(i) )
         a1 = a2 - a5 
c
c        # Compute the waves.
c        # Note that the 1+k-waves, for 1 .le. k .le. Nsp travel at
c        # the same speed and are lumped together in wave(.,.,2).
c        # The 3-wave is then stored in wave(.,.,3).
c
         do k = 1, Nsp
c           # 1-wave
            wave(i,k,1) = a1*Y(k,i)
c           # 2-wave
            wave(i,k,2) = delta(k) - Y(k,i)*a2
c           # 3-wave
            wave(i,k,3) = a5*Y(k,i)
         enddo
 
c        # 1-wave
         wave(i,mu,1) = a1*(u(i) - a(i))
         wave(i,mv,1) = a1*v(i)
         wave(i,mw,1) = a1*w(i)
         wave(i,mE,1) = a1*(enth(i) - u(i)*a(i))
         wave(i,mT,1) = 0.d0
         s(i,1) = u(i)-a(i)
c
c        # 2-wave
         wave(i,mu,2) = (drho - a2)*u(i)
         wave(i,mv,2) = (drho - a2)*v(i) + a3
         wave(i,mw,2) = (drho - a2)*w(i) + a4
         wave(i,mE,2) = (drho - a2)*u2v2w2(i) 
     &        - dpdr + dpY(i)*a2         + a3*v(i) + a4*w(i)
         wave(i,mT,2) = 0.d0
         s(i,2) = u(i)
c
c        # 3-wave
         wave(i,mu,3) = a5*(u(i) + a(i))
         wave(i,mv,3) = a5*v(i)
         wave(i,mw,3) = a5*w(i)
         wave(i,mE,3) = a5*(enth(i) + u(i)*a(i))
         wave(i,mT,3) = 0.d0
         s(i,3) = u(i)+a(i)
c                  
   30 continue
c
c     # compute flux differences as
c     #  (+/-)
c     # A     (Ur-Ul) = 0.5*( f(Ur)-f(Ul) +/- |A|(Ur-Ul) )
c     --------------------------
c
      call flx3(ixyz,maxm,meqn,mbc,mx,qr,maux,auxr,apdq)
      call flx3(ixyz,maxm,meqn,mbc,mx,ql,maux,auxl,amdq)
c
      do 35 i = 1-mbc, mx+mbc
         do 35 m=1,meqn
            fl(i,m) = amdq(i,m)
            fr(i,m) = apdq(i,m)
 35   continue      
c
      if (roe) then
         do 40 i = 2-mbc, mx+mbc
            do 40 m=1,meqn
               amdq(i,m) = 0.5d0*(fl(i,m)-fr(i-1,m))
 40      continue
c     
         do 50 i = 2-mbc, mx+mbc
            do 50 m=1,meqn
               sw = 0.d0
               do 60 mws=1,mwaves
                  sl = dabs(s(i,mws))
                  if (efix.and.dabs(s(i,mws)).lt.auxl(i,ixyz,2)) 
     &               sl = s(i,mws)**2/(2.d0*auxl(i,ixyz,2))+
     &                    0.5d0*auxl(i,ixyz,2)
                  sw = sw + sl*wave(i,m,mws)
 60            continue
               amdq(i,m) = amdq(i,m) - 0.5d0*sw
               apdq(i,m) = amdq(i,m) + sw
 50      continue
      endif
c
      if (hll) then
         do 55 i = 2-mbc, mx+mbc
c     # set this to hllfix = .true. when debugging HLL
            hllfix = .false.
            if (.not.roe) hllfix = .true.
c     
            rhol   = 0.d0
            rhoWl  = 0.d0
            do k = 1, Nsp
               rkl(k) = qr(i-1,k) + wave(i,k,1)
               rhol   = rhol   + rkl(k)
               rhoWl  = rhoWl  + rkl(k)/Wk(k)
            enddo
            rhou  = qr(i-1,mu) + wave(i,mu,1)
            rhov  = qr(i-1,mv) + wave(i,mv,1)
            rhow  = qr(i-1,mw) + wave(i,mw,1)
            rhoEl = qr(i-1,mE) + wave(i,mE,1)
            Tl    = qr(i-1,mT)
            rhoel = rhoEl - 0.5d0*(rhou**2+rhov**2+rhow**2)/rhol
            call SolveTrhok( Tl, rhoel, rhoWl, rkl, Nsp, ifail)
            rhoCpl = avgtabip( Tl, rkl, cpk, Nsp )
            gammal = RU / ( rhoCpl/rhoWl - RU ) + 1.d0
            pl = rhoWl*RU*Tl
            al = dsqrt(gammal*pl/rhol)
            if (rhol.le.0.d0.or.pl.le.0.d0) hllfix = .true.
c     
            rhor   = 0.d0
            rhoWr  = 0.d0
            do k = 1, Nsp
               rkr(k) = ql(i  ,k) - wave(i,k,3)
               rhor   = rhor   + rkr(k)
               rhoWr  = rhoWr  + rkr(k)/Wk(k)
            enddo
            rhou  = ql(i  ,mu) - wave(i,mu,3)
            rhov  = ql(i  ,mv) - wave(i,mv,3)
            rhow  = ql(i  ,mw) - wave(i,mw,3)
            rhoEr = ql(i  ,mE) - wave(i,mE,3)
            Tr    = ql(i  ,mT)
            rhoer = rhoEr - 0.5d0*(rhou**2+rhov**2+rhow**2)/rhor
            call SolveTrhok( Tr, rhoer, rhoWr, rkr, Nsp, ifail)
            rhoCpr = avgtabip( Tr, rkr, cpk, Nsp )
            gammar = RU / ( rhoCpr/rhoWr - RU ) + 1.d0
            pr = rhoWr*RU*Tr
            ar = dsqrt(gammar*pr/rhor)
            if (rhor.le.0.d0.or.pr.le.0.d0) hllfix = .true.
c     
            if (hllfix) then
c               if (roe) write (6,*) 'Switching to HLL in',i
c     
               rhol  = 0.d0
               rhoWl = 0.d0
               do k = 1, Nsp
                  rkl(k) = qr(i-1,k) 
                  rhol   = rhol  + qr(i-1,k)
                  rhoWl  = rhoWl + qr(i-1,k)/Wk(k)
               enddo
               ul = qr(i-1,mu)/rhol
               Tl = qr(i-1,mT)
               pl = rhoWl*RU*Tl
               rhoCpl = avgtabip( Tl, rkl, cpk, Nsp )
               gammal = RU / ( rhoCpl/rhoWl - RU ) + 1.d0
               al = dsqrt(gammal*pl/rhol)
c
               rhor  = 0.d0
               rhoWr = 0.d0
               do k = 1, Nsp
                  rkr(k) = ql(i  ,k) 
                  rhor   = rhor  + ql(i  ,k)
                  rhoWr  = rhoWr + ql(i  ,k)/Wk(k)
               enddo
               ur = ql(i  ,mu)/rhor
               Tr = ql(i  ,mT)
               pr = rhoWr*RU*Tr
               rhoCpr = avgtabip( Tr, rkr, cpk, Nsp )
               gammar = RU / ( rhoCpr/rhoWr - RU ) + 1.d0
               ar = dsqrt(gammar*pr/rhor)
c
               sl = dmin1(ul-al,ur-ar)
               sr = dmax1(ul+al,ur+ar)
c
               do m=1,meqn
                  if (sl.ge.0.d0) fg = fr(i-1,m)
                  if (sr.le.0.d0) fg = fl(i,m)
                  if (sl.lt.0.d0.and.sr.gt.0.d0) 
     &               fg = (sr*fr(i-1,m) - sl*fl(i,m) + 
     &                     sl*sr*(ql(i,m)-qr(i-1,m)))/ (sr-sl)
                  amdq(i,m) =   fg-fr(i-1,m)
                  apdq(i,m) = -(fg-fl(i  ,m))
               enddo
               amdq(i,mT) = 0.d0
               s(i,1) = sl
               s(i,2) = 0.d0
               s(i,3) = sr
            endif     
 55      continue
      endif
c
      if (pfix) then
         do 70 i=2-mbc,mx+mbc
            amdr = 0.d0
            apdr = 0.d0
            rhol = 0.d0
            rhor = 0.d0
            do k = 1, Nsp
               amdr = amdr + amdq(i,k)
               apdr = apdr + apdq(i,k)
               rhol = rhol + qr(i-1,k)
               rhor = rhor + ql(i  ,k)
            enddo
            do 70 k=1,Nsp
               if (qr(i-1,mu)+amdr.gt.0.d0) then
                  Z = qr(i-1,k)/rhol
               else
                  Z = ql(i  ,k)/rhor
               endif
               amdq(i,k) = Z*amdr + (Z-qr(i-1,k)/rhol)*qr(i-1,mu)
               apdq(i,k) = Z*apdr - (Z-ql(i  ,k)/rhor)*ql(i  ,mu)               
 70      continue  
      endif  
c
      return
      end
c
c
c =========================================================
      double precision function Cpmix( Tl, Tr, hl, hr, Y )
c =========================================================
      implicit double precision(a-h,o-z)
      include  "ck.i"
c
      dimension Y(*)
      dimension hl(*), hr(*)
      data Tol /1.d-6/
c
      if( dabs(Tr-Tl).gt.Tol ) then
         Cp = 0.d0
         do k = 1, Nsp
            Cp = Cp + (hr(k)-hl(k)) * Y(k) 
         enddo
         Cp = Cp / (Tr-Tl)
      else
         T = 0.5d0*(Tr+Tl)
         Cp = avgtabip( T, Y, cpk, Nsp )
      endif
      Cpmix = Cp
c
      return
      end

<