• VTF
  • FSI
  • AMROC
  • SFC
  • Motion
  • STLIB
  • Main Page
  • src/1d/equations/euler/rpm/rp1eum.f

    c
    c
    c     =========================================================
          subroutine rp1eum(maxmx,meqn,mwaves,mbc,mx,ql,qr,maux,
         &                  auxl,auxr,wave,s,amdq,apdq)
    c     =========================================================
    c
    c     # solve Riemann problems for the 1D two-component 
    c     # Euler equations using Roe's approximate Riemann solver.  
    c     
    c     # Keh-Ming Shyue "An efficient shock-capturing algorithm for
    c     # compressible multicomponent problems", J. Comput. Phys., Vol. 142, 
    c     # pp 208-242, 1998
    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     # On output, wave contains the waves, 
    c     #            s the speeds, 
    c     #            amdq the  left-going flux difference  A^- \Delta q
    c     #            apdq the right-going flux difference  A^+ \Delta q
    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 routine step1, rp is called with ql = qr = q.
    c     
    c     # Copyright (C) 2002 Ralf Deiterding
    c     # Brandenburgische Universitaet Cottbus
    c
    c     # Copyright (C) 2003-2007 California Institute of Technology
    c     # Ralf Deiterding, ralf@cacr.caltech.edu
    c
          implicit double precision (a-h,o-z)
          dimension   ql(1-mbc:maxmx+mbc, meqn)
          dimension   qr(1-mbc:maxmx+mbc, meqn)
          dimension    s(1-mbc:maxmx+mbc, mwaves)
          dimension wave(1-mbc:maxmx+mbc, meqn, mwaves)
          dimension amdq(1-mbc:maxmx+mbc, meqn)
          dimension apdq(1-mbc:maxmx+mbc, meqn)
    c
    c     # local storage
    c     ---------------
          parameter (max2 = 100002)  !# assumes at most 100000 grid points with mbc=2
          dimension delta(5)
          dimension u(-1:max2),enth(-1:max2),a(-1:max2),
         &     g1a2(-1:max2),euv(-1:max2),p(-1:max2) 
          logical efix
    c
          data efix /.true./    !# use entropy fix for transonic rarefactions
    c
    c     # Riemann solver returns flux differences
    c     ------------
          common /rpnflx/ mrpnflx
          mrpnflx = 0
    c
    c     # Compute Roe-averaged quantities:
    c
          do 10 i = 2-mbc, mx+mbc
             rhsqrtl = dsqrt(qr(i-1,1))
             rhsqrtr = dsqrt(ql(i,1))
             pl = (qr(i-1,3) - 0.5d0*(qr(i-1,2)**2)/qr(i-1,1) 
         &        - qr(i-1,5) ) / qr(i-1,4)
             pr = (ql(i,3) - 0.5d0*(ql(i,2)**2)/ql(i,1) 
         &        - ql(i,5) ) / ql(i,4)
             rhsq2 = rhsqrtl + rhsqrtr
             gamma1 = rhsq2 / ( qr(i-1,4)*rhsqrtl + ql(i,4)*rhsqrtr ) 
             xjota = ( pl*qr(i-1,4)*rhsqrtl + pr*ql(i,4)*rhsqrtr ) / rhsq2
             p(i) = xjota*gamma1                      
             u(i) = (qr(i-1,2)/rhsqrtl + ql(i,2)/rhsqrtr) / rhsq2
             enth(i) = (((qr(i-1,3)+pl)/rhsqrtl
         &             + (ql(i,3)+pr)/rhsqrtr)) / rhsq2                      
             a2 = gamma1*(enth(i) - .5d0*u(i)**2)
             a(i) = dsqrt(a2)
             g1a2(i) = gamma1 / a2
             euv(i) = enth(i) - u(i)**2
     10   continue
    c
    c
    c     # now split the jump in q at each interface into waves
    c
    c     # find a1 thru a5, the coefficients of the 5 eigenvectors:
          do 20 i = 2-mbc, mx+mbc
             do n = 1, 5
                delta(n) = ql(i,n) - qr(i-1,n)
             enddo
             a2 = g1a2(i) * (euv(i)*delta(1) + u(i)*delta(2) - delta(3)
         &      + p(i)*delta(4) + delta(5))
             a3 = (delta(2) + (a(i)-u(i))*delta(1) - a(i)*a2) / (2.d0*a(i))
             a1 = delta(1) - a2 - a3
             a4 = delta(4)
             a5 = delta(5)
    c
    c        # Compute the waves.
    c        # Note that the 2-wave as well as the 4-wave and 5-wave
    c        # travel at the same speed and are lumped together in wave(.,.,2).
    c        # The 3-wave is then stored in wave(.,.,3).
    c
             wave(i,1,1) = a1
             wave(i,2,1) = a1*(u(i)-a(i))
             wave(i,3,1) = a1*(enth(i) - u(i)*a(i))
             wave(i,4,1) = 0.d0
             wave(i,5,1) = 0.d0 
             s(i,1) = u(i)-a(i)
    c
             wave(i,1,2) = a2
             wave(i,2,2) = a2*u(i)
             wave(i,3,2) = a2*0.5d0*u(i)**2 + a4*p(i) + a5
             wave(i,4,2) =                    a4
             wave(i,5,2) =                              a5
             s(i,2) = u(i)
    c
             wave(i,1,3) = a3
             wave(i,2,3) = a3*(u(i)+a(i))
             wave(i,3,3) = a3*(enth(i)+u(i)*a(i))
             wave(i,4,3) = 0.d0
             wave(i,5,3) = 0.d0 
             s(i,3) = u(i)+a(i)
       20 continue
    c
    c
    c     # compute flux differences amdq and apdq.
    c     ---------------------------------------
    c
          if (efix) go to 110
    c
    c     # no entropy fix
    c     ----------------
    c
    c     # amdq = SUM s*wave   over left-going waves
    c     # apdq = SUM s*wave   over right-going waves
    c
          do 100 m=1,5
             do 100 i=2-mbc, mx+mbc
                amdq(i,m) = 0.d0
                apdq(i,m) = 0.d0
                do 90 mw=1,mwaves
                   if (s(i,mw) .lt. 0.d0) then
                      amdq(i,m) = amdq(i,m) + s(i,mw)*wave(i,m,mw)
                   else
                      apdq(i,m) = apdq(i,m) + s(i,mw)*wave(i,m,mw)
                   endif
       90       continue
      100 continue
          go to 900     
    c
    c-----------------------------------------------------
    c
      110 continue
    c
    c     # With entropy fix
    c     ------------------
    c
    c     # compute flux differences amdq and apdq.
    c     # First compute amdq as sum of s*wave for left going waves.
    c     # Incorporate entropy fix by adding a modified fraction of wave
    c     # if s should change sign.
    c
          do 200 i = 2-mbc, mx+mbc
    c
    c        # check 1-wave:
    c        ---------------
    c     
             rhoim1 = qr(i-1,1)
             pim1 = (qr(i-1,3) - 0.5d0*(qr(i-1,2)**2)/qr(i-1,1) 
         &        - qr(i-1,5) ) / qr(i-1,4)
             gamma1 = 1.d0/qr(i-1,4)
             gamma = gamma1 + 1.d0
             pinf = qr(i-1,5)*gamma1/gamma
             cim1 = dsqrt(gamma*(pim1+pinf)/rhoim1)
             s0 = qr(i-1,2)/rhoim1 - cim1 !# u-c in left state (cell i-1)
             
    c        # check for fully supersonic case:
             if (s0.ge.0.d0 .and. s(i,1).gt.0.d0)  then
    c           # everything is right-going
                do 60 m=1,5
                   amdq(i,m) = 0.d0
     60         continue
                go to 200 
             endif
    c
             rho1 = qr(i-1,1) + wave(i,1,1)
             rhou1 = qr(i-1,2) + wave(i,2,1)
             en1 = qr(i-1,3) + wave(i,3,1)
             p1 = (en1-0.5d0*(rhou1**2)/rho1-qr(i-1,5))/qr(i-1,4)
             c1 = dsqrt(gamma*(p1+pinf)/rho1)
             s1 = rhou1/rho1 - c1   !# u-c to right of 1-wave
             if (s0.lt.0.d0 .and. s1.gt.0.d0) then
    c           # transonic rarefaction in the 1-wave
                sfract = s0 * (s1-s(i,1)) / (s1-s0)
             else if (s(i,1) .lt. 0.d0) then
    c           # 1-wave is leftgoing
                sfract = s(i,1)
             else
    c           # 1-wave is rightgoing
                sfract = 0.d0       !# this shouldn't happen since s0 < 0
             endif
             do 120 m=1,5
                amdq(i,m) = sfract*wave(i,m,1)
     120     continue
    c
    c        # check 2-wave:
    c        ---------------
    c
             if (s(i,2) .ge. 0.d0) go to 200 !# 2- and 3- waves are rightgoing
             do 140 m=1,5
                amdq(i,m) = amdq(i,m) + s(i,2)*wave(i,m,2)
     140     continue
    c
    c        # check 3-wave:
    c        ---------------
    c
             rhoi = ql(i,1)
             pi = (ql(i,3) - 0.5d0*(ql(i,2)**2)/ql(i,1) 
         &        - ql(i,5) ) / ql(i,4)
             gamma1 = 1.d0/ql(i,4)
             gamma = gamma1 + 1.d0
             pinf = ql(i,5)*gamma1/gamma
             ci = dsqrt(gamma*(pi+pinf)/rhoi)
             s3 = ql(i,2)/rhoi + ci     !# u+c in right state  (cell i)
    c
             rho2 = ql(i,1) - wave(i,1,3)
             rhou2 = ql(i,2) - wave(i,2,3)
             en2 = ql(i,3) - wave(i,3,3)
             p2 = (en2-0.5d0*(rhou2**2)/rho2-ql(i,5))/ql(i,4)
             c2 = dsqrt(gamma*(p2+pinf)/rho2)
             s2 = rhou2/rho2 + c2   !# u+c to left of 3-wave
             if (s2 .lt. 0.d0 .and. s3.gt.0.d0) then
    c           # transonic rarefaction in the 3-wave
                sfract = s2 * (s3-s(i,3)) / (s3-s2)
             else if (s(i,3) .lt. 0.d0) then
    c           # 3-wave is leftgoing
                sfract = s(i,3)
             else 
    c           # 3-wave is rightgoing
                go to 200
             endif
    c
             do 160 m=1,5
                amdq(i,m) = amdq(i,m) + sfract*wave(i,m,3)
     160     continue
     200  continue
    c
    c     # compute the rightgoing flux differences:
    c     # df = SUM s*wave   is the total flux difference and apdq = df - amdq
    c
          do 220 m=1,5
             do 220 i = 2-mbc, mx+mbc
                df = 0.d0
                do 210 mw=1,mwaves
                   df = df + s(i,mw)*wave(i,m,mw)
     210        continue
                apdq(i,m) = df - amdq(i,m)
     220  continue
    c
     900  continue
          return
          end
    
    

<