c c c ===================================================== subroutine rpn2meu(ixy,maxm,meqn,mwaves,mbc,mx,ql,qr,maux, & auxl,auxr,wave,s,amdq,apdq) c ===================================================== c c # Solve Riemann problems for the 2D two-component Euler equations c # using HLLC. Use flux difference splitting formulation for full c # compatibility to Wave Propagation Method. 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 c # This data is along a slice in the x-direction if ixy=1 c # or the y-direction if ixy=2. 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 routines, this routine is called with ql = qr 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 wave(1-mbc:maxm+mbc, meqn, mwaves) dimension s(1-mbc:maxm+mbc, mwaves) dimension ql(1-mbc:maxm+mbc, meqn) dimension qr(1-mbc:maxm+mbc, meqn) dimension auxl(1-mbc:maxm+mbc, maux) dimension auxr(1-mbc:maxm+mbc, maux) dimension apdq(1-mbc:maxm+mbc, meqn) dimension amdq(1-mbc:maxm+mbc, meqn) c c local arrays -- common block comroe is passed to rpt2eu c ------------ parameter (maxm2 = 10005) !# assumes at most 10000x10000 grid with mbc=5 parameter (minm2 = -4) !# assumes at most mbc=5 dimension qls(4), qrs(4) logical roespeed c common /comroe/ u2v2(minm2:maxm2), & u(minm2:maxm2),v(minm2:maxm2),enth(minm2:maxm2), & a(minm2:maxm2),g1a2(minm2:maxm2),euv(minm2:maxm2), & p(minm2:maxm2) c data roespeed /.false./ !# use Roe average for wave speed estimation c c # Riemann solver returns flux differences c ------------ common /rpnflx/ mrpnflx mrpnflx = 0 c if (minm2.gt.1-mbc .or. maxm2 .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 to the orthogonal c # momentum: c if (ixy.eq.1) then mu = 2 mv = 3 else mu = 3 mv = 2 endif c c # note that notation for u and v reflects assumption that the c # Riemann problems are in the x-direction with u in the normal c # direciton and v in the orthogonal direcion, but with the above c # definitions of mu and mv the routine also works with ixy=2 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 rpt2eu to do the transverse wave splitting. c do 10 i = 2-mbc, mx+mbc if (qr(i-1,1).le.0.d0.or.ql(i,1).le.0.d0) then write (6,*) 'Unrecoverable error in density',i write (6,*) qr(i-1,1),ql(i,1) stop endif c rl = qr(i-1,1) ul = qr(i-1,mu)/rl vl = qr(i-1,mv)/rl El = qr(i-1,4) gammal1 = 1.d0/qr(i-1,5) gammal = gammal1 + 1.d0 pinfl = qr(i-1,6)*gammal1/gammal pl = (El - 0.5d0*(ul**2+vl**2)*rl - qr(i-1,6))/qr(i-1,5) if (pl+pinfl.le.0.d0.or.gammal.le.0.d0) then write (6,*) 'Unrecoverable error in speed of sound l',i write (6,*) pl,pinfl,pl+pinfl,gammal stop endif al = dsqrt(gammal*(pl+pinfl)/rl) c rr = ql(i ,1) ur = ql(i ,mu)/rr vr = ql(i ,mv)/rr Er = ql(i ,4) gammar1 = 1.d0/ql(i ,5) gammar = gammar1 + 1.d0 pinfr = ql(i ,6)*gammar1/gammar pr = (Er - 0.5d0*(ur**2+vr**2)*rr - ql(i ,6))/ql(i ,5) if (pr+pinfr.le.0.d0.or.gammar.le.0.d0) then write (6,*) 'Unrecoverable error in speed of sound r',i write (6,*) pr,pinfr,pr+pinfr,gammar stop endif ar = dsqrt(gammar*(pr+pinfr)/rr) c rhsqrtl = dsqrt(qr(i-1,1)) rhsqrtr = dsqrt(ql(i,1)) rhsq2 = rhsqrtl + rhsqrtr gamma1 = rhsq2 / ( qr(i-1,5)*rhsqrtl + ql(i,5)*rhsqrtr ) xjota = ( pl*qr(i-1,5)*rhsqrtl + pr*ql(i,5)*rhsqrtr ) / rhsq2 p(i) = xjota*gamma1 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 u2v2(i) = u(i)**2 + v(i)**2 enth(i) = (((qr(i-1,4)+pl)/rhsqrtl & + (ql(i,4)+pr)/rhsqrtr)) / rhsq2 a2 = gamma1*(enth(i) - .5d0*u2v2(i)) a(i) = dsqrt(a2) if (a2.le.0.d0) a(i) = dmax1(al,ar) g1a2(i) = gamma1 / a2 euv(i) = enth(i) - u2v2(i) c if (roespeed) then sl = u(i)-a(i) sr = u(i)+a(i) else sl = dmin1(ul-al,ur-ar) sr = dmax1(ul+al,ur+ar) endif ss = (pr-pl+rl*ul*(sl-ul)-rr*ur*(sr-ur))/ & (rl*(sl-ul)-rr*(sr-ur)) c qrs(1) = rr*(sr-ur)/(sr-ss) qrs(mu) = qrs(1)*ss qrs(mv) = qrs(1)*vr qrs(4) = qrs(1)*(Er/rr+ & (ss-ur)*(ss+pr/(rr*(sr-ur)))) c qls(1) = rl*(sl-ul)/(sl-ss) qls(mu) = qls(1)*ss qls(mv) = qls(1)*vl qls(4) = qls(1)*(El/rl+ & (ss-ul)*(ss+pl/(rl*(sl-ul)))) c do m=1,4 wave(i,m,1) = qls(m) - qr(i-1,m) wave(i,m,2) = qrs(m) - qls(m) wave(i,m,3) = ql(i,m) - qrs(m) enddo do m=5,6 wave(i,m,1) = 0.d0 wave(i,m,2) = ql(i,m) - qr(i-1,m) wave(i,m,3) = 0.d0 enddo c s(i,1) = sl s(i,2) = ss s(i,3) = sr c do m=1,meqn amdq(i,m) = 0.d0 apdq(i,m) = 0.d0 do 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 enddo enddo 10 continue return end c