Scrambler
1
|
00001 !######################################################################### 00002 ! 00003 ! Copyright (C) 2003-2012 Department of Physics and Astronomy, 00004 ! University of Rochester, 00005 ! Rochester, NY 00006 ! 00007 ! riemann_solvers.f90 is part of AstroBEAR. 00008 ! 00009 ! AstroBEAR is free software: you can redistribute it and/or modify 00010 ! it under the terms of the GNU General Public License as published by 00011 ! the Free Software Foundation, either version 3 of the License, or 00012 ! (at your option) any later version. 00013 ! 00014 ! AstroBEAR is distributed in the hope that it will be useful, 00015 ! but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 ! GNU General Public License for more details. 00018 ! 00019 ! You should have received a copy of the GNU General Public License 00020 ! along with AstroBEAR. If not, see <http://www.gnu.org/licenses/>. 00021 ! 00022 !######################################################################### 00025 00029 00030 00033 MODULE RiemannSolvers 00034 USE GlobalDeclarations 00035 USE HyperbolicDeclarations 00036 USE PhysicsDeclarations 00037 USE EOS 00038 IMPLICIT NONE 00039 INTEGER, PARAMETER, PUBLIC :: iHLL=3, iHLLC=2, iHLLD=6, iExactRS=1, iHLLC_ISO=4, iHLLD_ISO=8, iHLL_MHD=7 00040 CONTAINS 00041 00047 FUNCTION calc_flux(wl, wr, f, lambda_max) 00048 REAL(KIND=qPREC) :: calc_flux 00049 REAL(KIND=qPrec), DIMENSION(:) :: f, wl, wr 00050 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00051 INTEGER :: i 00052 IF (iSolver == 0) THEN !select default solver 00053 IF (lMHD) THEN 00054 IF (iEOS == EOS_ISOTHERMAL) THEN 00055 iSolver=iHLLD_ISO 00056 ELSE 00057 iSolver = iHLLD 00058 END IF 00059 ELSE 00060 IF (iEOS == EOS_ISOTHERMAL) THEN 00061 iSolver = iHLLC_ISO 00062 ELSE 00063 iSolver = iHLLC 00064 END IF 00065 END IF 00066 END IF 00067 IF (present(lambda_max)) THEN 00068 SELECT CASE (iSolver) 00069 CASE (iHLL) 00070 calc_flux=HLL(wl,wr,f,lambda_max) 00071 CASE(iHLL_MHD) 00072 calc_flux=HLL_MHD(wl,wr,f,lambda_max) 00073 CASE (iHLLC) 00074 calc_flux=HLLC(wl,wr,f,lambda_max) 00075 CASE (iHLLC_ISO) 00076 calc_flux=HLLC_ISO(wl,wr,f,lambda_max) 00077 CASE (iHLLD) 00078 calc_flux=HLLD(wl,wr,f,lambda_max) 00079 CASE (iHLLD_ISO) 00080 calc_flux=HLLD_ISO(wl,wr,f,lambda_max) 00081 CASE (iExactRS) 00082 calc_flux=ExactRS(wl,wr,f,lambda_max) 00083 END SELECT 00084 ELSE 00085 SELECT CASE (iSolver) 00086 CASE (iHLL) 00087 calc_flux=HLL(wl,wr,f) 00088 CASE(iHLL_MHD) 00089 calc_flux=HLL_MHD(wl,wr,f) 00090 CASE (iHLLC) 00091 calc_flux=HLLC(wl,wr,f) 00092 CASE (iHLLC_ISO) 00093 calc_flux=HLLC_ISO(wl,wr,f) 00094 CASE (iHLLD) 00095 calc_flux=HLLD(wl,wr,f) 00096 CASE (iHLLD_ISO) 00097 calc_flux=HLLD_ISO(wl,wr,f) 00098 CASE (iExactRS) 00099 calc_flux=ExactRS(wl,wr,f) 00100 END SELECT 00101 END IF 00102 DO i=1,size(f) 00103 IF (ISNAN(f(i))) THEN 00104 IF (.NOT. lRequestRestart) THEN 00105 write(*,*) 'processor', MPI_ID, 'requesting restart due to nan in flux' 00106 write(*,*) 'wl=', wl 00107 write(*,*) 'wr=', wr 00108 END IF 00109 lRequestRestart=.true. 00110 00111 ! write(*,*) 'found nan flux' 00112 ! write(*,'(15E13.2)') f 00113 ! write(*,'(15E13.2)') wl 00114 ! write(*,'(15E13.2)') wr 00115 ! STOP 00116 END IF 00117 00118 END DO 00119 END FUNCTION CALC_FLUX 00120 00126 FUNCTION HLLC_ISO(wl,wr,f, lambda_max) 00127 !wl and wr = (/rho,vx,vy,vz,Bx,By,Bz/) 00128 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00129 REAL(KIND=qPrec), DIMENSION(:) :: f, wl, wr 00130 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: ul, ur, usl, usr, usc, fl, fr 00131 REAL(KIND=qPrec) :: SL, SR, DM, SM, MM, pTr, pTl,ml, mr, f2l, f2r, sqrtdl, sqrtdr, um 00132 INTEGER :: n 00133 REAL(KIND=qPREC) :: HLLC_ISO 00134 n=SIZE(wl) 00135 ALLOCATE(ul(n), ur(n), usl(n), usr(n), fl(n), fr(n)) 00136 00137 !Einfeldt version 00138 sqrtdl=SQRT(wl(1)) 00139 sqrtdr=SQRT(wr(1)) 00140 um=(sqrtdl*wl(2)+sqrtdr*wr(2))/(sqrtdl+sqrtdr) 00141 SL = MIN(wl(2), um) - Iso_Speed 00142 SR = MAX(wr(2), um) + Iso_Speed 00143 00144 !Davis Version 00145 ! SL=min(wl(2),wr(2))-Iso_Speed 00146 ! SR=max(wl(2),wr(2))+Iso_Speed 00147 00148 00149 ! H Viscosity protection 00150 ! Sanders 9a 00151 ! IF (present(lambda_max)) THEN 00152 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00153 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00154 ! END IF 00155 00156 ! Sanders 9b 00157 IF (present(lambda_max)) THEN 00158 SL=SIGN(abs(SL)+lambda_max, SL) 00159 SR=SIGN(abs(SR)+lambda_max, SR) 00160 END IF 00161 00162 IF (SL >= 0d0) THEN 00163 CALL F_prim_Iso(wl, f) 00164 ELSE IF (SR <= 0d0) THEN 00165 CALL F_prim_Iso(wr, f) 00166 ELSE !Solve for common star variables DM, SM, and f(1:2) 00167 pTr=wr(1)*Iso_Speed2 00168 pTl=wl(1)*Iso_Speed2 00169 mr=wr(1)*wr(2) 00170 ml=wl(1)*wl(2) 00171 f2r=mr*wr(2)+pTr 00172 f2l=ml*wl(2)+pTl 00173 DM = ((SR*wr(1)-SL*wl(1))-(mr-ml))/(SR-SL) 00174 MM = ((SR*mr-SL*ml)-(f2r-f2l))/(SR-SL) 00175 f(1)=((SR*ml-SL*mr)+SR*SL*(wr(1)-wl(1)))/(SR-SL) 00176 f(2)=((SR*f2l-SL*f2r)+SR*SL*(mr-ml))/(SR-SL) 00177 SM=f(1)/DM !advective choice for middle speed (instead of MM/DM) 00178 00179 IF (SM >= 0) THEN !Need to find sleft values 00180 usl(1:2)=(/ DM, & 00181 MM /) 00182 usl(3:m_high)=DM*wl(3:m_high) 00183 usl(m_high+1:n)=DM*wl(m_high+1:n)/wl(1) 00184 00185 CALL F_prim_Iso(wl,fl) 00186 CALL rprim_to_cons_Iso(wl, ul) 00187 fl = fl+SL*(usl-ul) 00188 f=fl 00189 END IF 00190 IF (SM <= 0) THEN 00191 usr(1:2)=(/ DM, & 00192 MM /) 00193 00194 usr(3:m_high)=DM*wr(3:m_high) 00195 usr(m_high+1:n)=DM*wr(m_high+1:n)/wr(1) 00196 00197 CALL F_prim_Iso(wr,fr) 00198 CALL rprim_to_cons_Iso(wr, ur) 00199 fr = fr+SR*(usr-ur) 00200 f=fr 00201 END IF 00202 IF (SM == 0d0) THEN 00203 f=.5d0*(fl+fr) 00204 END IF 00205 00206 END IF 00207 HLLC_ISO=max(ABS(SL),ABS(SR)) 00208 DEALLOCATE(ul, ur, usl, usr, fl, fr) 00209 END FUNCTION HLLC_ISO 00210 00211 00217 FUNCTION HLLD_ISO(wl,wr,f, lambda_max) 00218 !wl and wr = (/rho,vx,vy,vz,Bx,By,Bz/) 00219 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00220 REAL(KIND=qPrec), DIMENSION(:) :: f, wl, wr 00221 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: ul, ur, usl, usr, usc 00222 REAL(KIND=qPrec) :: SL, SR, DM, SM, AM, MM, SsL, SsR, CFL, CFR, sqrtdr, sqrtdl, 00223 pTr, pTl, Bx, Bx_2, temp, temp2, ml, mr, f2l, f2r, X, um, Bm(3), bx2, b2, a2b2, CFM, sbx, Xi 00224 INTEGER :: n, i 00225 REAL(KIND=qPREC) :: HLLD_ISO 00226 mr=0 00227 ml=0 00228 n=SIZE(wl) 00229 ALLOCATE(ul(n), ur(n), usl(n), usr(n), usc(n)) 00230 Bx = half*(wr(5)+wl(5)) 00231 CFL = fast_speed_Iso(wl) 00232 CFR = fast_speed_Iso(wr) 00233 Bx_2 = Bx**2 00234 IF (Bx == 0) THEN 00235 sBx=0d0 00236 ELSE 00237 sBx=SIGN(1d0,Bx) 00238 END IF 00239 00240 !Einfeldt version (more computationally expensive but more accurate) 00241 sqrtdl=SQRT(wl(1)) 00242 sqrtdr=SQRT(wr(1)) 00243 um=(sqrtdl*wl(2)+sqrtdr*wr(2))/(sqrtdl+sqrtdr) 00244 Bm=(sqrtdl*wl(5:7)+sqrtdr*wr(5:7))/(sqrtdl+sqrtdr) 00245 dm=sqrtdl*sqrtdr 00246 bx2=Bm(1)**2/dm; 00247 b2=DOT_PRODUCT(Bm(2:3),Bm(2:3))/dm+bx2 00248 a2b2 = (Iso_Speed2+b2)/2d0 00249 CFM = sqrt((a2b2+sqrt(a2b2**2-Iso_Speed2*bx2))) 00250 00251 SL = MIN(wl(2)-CFL, um-CFM) 00252 SR = MAX(wr(2)+CFR, um+CFM) 00253 00254 !Davis Version 00255 ! SL=min(wl(2)-CFL,wr(2)-CFR) 00256 ! SR=max(wl(2)+CFL,wr(2)+CFR) 00257 00258 00259 ! H Viscosity protection 00260 ! Sanders 9a 00261 ! IF (present(lambda_max)) THEN 00262 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00263 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00264 ! END IF 00265 00266 ! Sanders 9b 00267 IF (present(lambda_max)) THEN 00268 SL=SIGN(abs(SL)+lambda_max, SL) 00269 SR=SIGN(abs(SR)+lambda_max, SR) 00270 END IF 00271 00272 00273 IF (SL >= 0d0) THEN 00274 CALL F_prim_MHD_Iso(wl, f) 00275 ELSE IF (SR <= 0d0) THEN 00276 CALL F_prim_MHD_Iso(wr, f) 00277 ELSE !Solve for common star variables DM, SM, and f(1:2) 00278 pTr=wr(1)*Iso_Speed2+half*(Bx_2+DOT_PRODUCT(wr(6:7),wr(6:7))) 00279 pTl=wl(1)*Iso_Speed2+half*(Bx_2+DOT_PRODUCT(wl(6:7),wl(6:7))) 00280 mr=wr(1)*wr(2) 00281 ml=wl(1)*wl(2) 00282 f2r=mr*wr(2)+pTr-Bx_2 00283 f2l=ml*wl(2)+pTl-Bx_2 00284 DM = ((SR*wr(1)-SL*wl(1))-(mr-ml))/(SR-SL) 00285 MM = ((SR*mr-SL*ml)-(f2r-f2l))/(SR-SL) 00286 f(1)=((SR*ml-SL*mr)+SR*SL*(wr(1)-wl(1)))/(SR-SL) 00287 f(2)=((SR*f2l-SL*f2r)+SR*SL*(mr-ml))/(SR-SL) 00288 00289 SM=f(1)/DM !advective choice for middle speed (instead of MM/DM) 00290 AM=ABS(BX)/sqrt(DM) !star alfven speed 00291 SsL=SM-AM !Slow waves 00292 SsR=SM+AM !Slow waves 00293 IF (SsR >= 0) THEN !Need to find sleft values 00294 IF (miniscule(SL-SsL)) THEN 00295 temp2=0d0 00296 temp=1d0 00297 ELSE 00298 temp=1d0/((SL-SsL)*(SL-SsR)) 00299 temp2=Bx*(SM-wl(2))*temp 00300 temp=(wl(1)*(SL-wl(2))**2 - Bx_2)*temp/DM 00301 END IF 00302 00303 usl(1:7)=(/ DM, & 00304 MM, & 00305 DM*wl(3)-temp2*wl(6), & 00306 DM*wl(4)-temp2*wl(7), & 00307 Bx, & 00308 wl(6)*temp, & 00309 wl(7)*temp & 00310 /) 00311 usl(8:n) = DM*wl(8:n)/wl(1) 00312 END IF 00313 IF (SsL <= 0) THEN !Need to find sright values 00314 IF (miniscule(SR-SsR)) THEN 00315 temp2=0d0 00316 temp=1d0 00317 ELSE 00318 temp=1d0/((SR-SsR)*(SR-SsL)) 00319 temp2=Bx*(SM-wr(2))*temp 00320 temp=(wr(1)*(SR-wr(2))**2 - Bx_2)*temp/DM 00321 END IF 00322 00323 00324 usr(1:7)=(/ DM, & 00325 MM, & 00326 DM*wr(3)-temp2*wr(6), & 00327 DM*wr(4)-temp2*wr(7), & 00328 Bx, & 00329 wr(6)*temp, & 00330 wr(7)*temp & 00331 /) 00332 usr(8:n) = DM*wr(8:n)/wr(1) 00333 END IF 00334 00335 !Now need to find region 00336 IF (SsR >= 0 .AND. Ssl <= 0) THEN !in center region 00337 X=sbx*sqrt(DM) !sign(sqrt(DM),Bx) 00338 Xi=1d0/X 00339 IF (Bx==0) Xi=0 00340 usc(3:7) = (/ half*((usl(3)+usr(3))+X*(usr(6)-usl(6))), & 00341 half*((usl(4)+usr(4))+X*(usr(7)-usl(7))), & 00342 Bx, & 00343 half*((usr(6)+usl(6))+(usr(3)-usl(3))*Xi), & 00344 half*((usr(7)+usl(7))+(usr(4)-usl(4))*Xi) & 00345 /) 00346 00347 f(3:7) = (/ usc(3)*SM-Bx*usc(6) ,& 00348 usc(4)*SM-Bx*usc(7) ,& 00349 0d0 ,& 00350 usc(6)*SM-Bx*usc(3)/DM ,& 00351 usc(7)*SM-BX*usc(4)/DM & 00352 /) 00353 f(8:n)=((SR*wl(2)*wl(8:n)-SL*wr(2)*wr(8:n))+SR*SL*(wr(8:n)-wl(8:n)))/(SR-SL) 00354 ELSEIF (SsL >= 0) THEN !in sleft region 00355 CALL F_prim_MHD_Iso(wl,f) 00356 CALL rprim_to_cons_MHD_Iso(wl, ul) 00357 f = f+SL*(usl-ul) 00358 ELSE !in right region 00359 CALL F_prim_MHD_Iso(wr,f) 00360 CALL rprim_to_cons_MHD_Iso(wr, ur) 00361 f = f+SR*(usr-ur) 00362 END IF 00363 END IF 00364 f(5)=0d0 !kill normal magnetic field 00365 HLLD_ISO=max(ABS(SL),ABS(SR)) 00366 ! DO i=1,size(f) 00367 ! IF (ISNAN(f(i))) THEN 00368 ! write(*,*) 'HLLD_ISO returned a nan' 00369 ! write(*,'(A,10E25.16)') 'wl=', wl 00370 ! write(*,'(A,10E25.16)') 'wr=', wr 00371 ! write(*,'(A,10E25.16)') 'f=', f 00372 ! write(*,*) Sl, Ssl, SM, SsR, SR, ((SR*ml-SL*mr)+SR*SL*(wr(1)-wl(1)))/(SR-SL), DM, mr, ml 00373 ! write(*,*) '=============================' 00374 ! STOP 00375 ! END IF 00376 ! END DO 00377 DEALLOCATE(ul, ur, usl, usr, usc) 00378 END FUNCTION HLLD_ISO 00379 00380 00386 FUNCTION HLLD(wl,wr,f,lambda_max) 00387 REAL(KIND=qPREC) :: HLLD 00388 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00389 REAL(KIND=qPrec), DIMENSION(:) :: f, wl, wr 00390 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: ul, ur, rr, rl, rm, usl, usr, ussl, ussr, fl, fr 00391 REAL(KIND=qPrec) :: SL, SR, SM, SsL, SsR, sqrtdl,sqrtdr,sqrtdsl, sqrtdsr, CFL, CFR, CFM, 00392 pTr, pTl, pT, Bx, Bx_2, dsr, dsl, temp, Brat, vrat, tempr, templ, a2, aT2, b2, X, sumsqrtds, sBx 00393 REAL(KIND=qPrec), DIMENSION(3) :: vss, Bss 00394 REAL(KIND=qPrec), DIMENSION(2) :: db 00395 INTEGER :: n,i 00396 n=SIZE(wl) 00397 ALLOCATE(ul(n), ur(n), rr(n), rl(n), rm(n), usl(n), usr(n), ussl(n), ussr(n), fl(n), fr(n)) 00398 CALL rprim_to_cons_MHD(wl,ul) 00399 CALL rprim_to_cons_MHD(wr,ur) 00400 00401 Bx = half*(wr(6)+wl(6)) 00402 IF (Bx == 0) THEN 00403 sBx=0d0 00404 ELSE 00405 sBx=SIGN(1d0,Bx) 00406 END IF 00407 sqrtdl=SQRT(wl(1)) 00408 sqrtdr=SQRT(wr(1)) 00409 CFL = fast_speed(wl) 00410 CFR = fast_speed(wr) 00411 Bx_2 = Bx**2 00412 00413 !Einfeldt version (more computationally expensive but more accurate) 00414 CALL prim_to_Roe_MHD(wl,rl) 00415 CALL prim_to_Roe_MHD(wr,rr) 00416 sqrtdl=SQRT(wl(1)) 00417 sqrtdr=SQRT(wr(1)) 00418 rm(2:)=(sqrtdl*rl(2:)+sqrtdr*rr(2:))/(sqrtdl+sqrtdr) 00419 rm(1)=sqrtdl*sqrtdr 00420 CFM=fast_speed_Roe(rm) 00421 SL = MIN(wl(3)-CFL, rm(3)-CFM) 00422 SR = MAX(wr(3)+CFR, rm(3)+CFM) 00423 00424 !Davis Version 00425 ! SL=min(wl(3)-CFL,wr(3)-CFR) 00426 ! SR=max(wl(3)+CFL,wr(3)+CFR) 00427 00428 !protection to ensure proper behavior when flow is supersonic 00429 00430 00431 ! H Viscosity protection 00432 ! Sanders 9a 00433 ! IF (present(lambda_max)) THEN 00434 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00435 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00436 ! END IF 00437 00438 ! Sanders 9b 00439 ! write(*,'(A,4E25.16)') 'A', SL, SR, CFL, CFR 00440 ! write(*,'(A,4E25.16)') 'B', rm(3), CFM 00441 IF (present(lambda_max)) THEN 00442 SL=SIGN(abs(SL)+lambda_max, SL) 00443 SR=SIGN(abs(SR)+lambda_max, SR) 00444 ! write(*,'(A,4E25.16)') 'C', SL, SR, lambda_max 00445 END IF 00446 00447 IF (SL >= 0d0) THEN !Region 1 00448 CALL F_prim_MHD(wl, f) 00449 ELSE IF (SR <= 0d0) THEN !Region 6 00450 CALL F_prim_MHD(wr, f) 00451 ELSE !Regions 2-5 00452 pTr=wr(2)+half*DOT_PRODUCT(wr(6:8),wr(6:8)) 00453 pTl=wl(2)+half*DOT_PRODUCT(wl(6:8),wl(6:8)) 00454 tempr=(SR-wr(3))*wr(1) 00455 templ=(SL-wl(3))*wl(1) 00456 SM = ((tempr*wr(3)-templ*wl(3))-(pTr-pTl))/(tempr - templ) 00457 pT = ((tempr*pTl - templ*pTr) + tempr*templ*(wr(3)-wl(3)))/(tempr-templ) 00458 dsl=templ/(SL-SM) 00459 dsr=tempr/(SR-SM) 00460 sqrtdsl = SQRT(dsl) 00461 sqrtdsr = SQRT(dsr) 00462 SsL=SM-ABS(Bx)/sqrtdsl 00463 SsR=SM+ABS(Bx)/sqrtdsr 00464 IF (SsL <= 0.d0) THEN !Regions !3-5 00465 IF (tempr*(SR-SM) == Bx_2) THEN 00466 vrat = 0d0 00467 Brat = 1d0 00468 ELSE 00469 temp = 1.d0/(tempr*(SR-SM)-Bx_2) 00470 vrat = (SM-wr(3))*temp 00471 Brat = (tempr**2/wr(1) - Bx_2)*temp 00472 END IF 00473 usr(1:8) = (/dsr, & 00474 0.d0, & 00475 dsr*SM, & 00476 dsr*(wr(4)-Bx*wr(7)*vrat), & 00477 dsr*(wr(5)-Bx*wr(8)*vrat), & 00478 Bx, & 00479 wr(7)*Brat, & 00480 wr(8)*Brat /) 00481 usr(2) = ((SR-wr(3))*ur(2)-pTr*wr(3)+pT*SM+Bx* & 00482 (DOT_PRODUCT(wr(3:5),wr(6:8))-DOT_PRODUCT(usr(3:5),usr(6:8))/dsr))/(SR-SM) 00483 usr(9:n) = dsr*wr(9:n)/wr(1) 00484 END IF 00485 IF (SsR >= 0d0) THEN !Regions 2-4 00486 IF (templ*(SL-SM) == Bx_2) THEN 00487 vrat=0.0 00488 Brat=1d0 00489 ELSE 00490 temp = 1.d0/(templ*(SL-SM)-Bx_2) 00491 vrat = (SM-wl(3))*temp 00492 Brat = (templ**2/wl(1) - Bx_2)*temp 00493 END IF 00494 usl(1:8) = (/dsl, & 00495 0.d0, & 00496 dsl*SM, & 00497 dsl*(wl(4)-Bx*wl(7)*vrat), & 00498 dsl*(wl(5)-Bx*wl(8)*vrat), & 00499 Bx, & 00500 wl(7)*Brat, & 00501 wl(8)*Brat /) 00502 usl(2) = ((SL-wl(3))*ul(2)-pTl*wl(3)+pT*SM+Bx* & 00503 (DOT_PRODUCT(wl(3:5),wl(6:8))-DOT_PRODUCT(usl(3:5),usl(6:8))/dsl))/(SL-SM) 00504 usl(9:n) = dsl*wl(9:n)/wl(1) 00505 END IF 00506 IF (SsL <= 0.d0 .AND. SsR >= 0.d0) THEN !Regions 3-4 00507 sumsqrtds = sqrtdsl+sqrtdsr 00508 vss = (/SM, & 00509 (usl(4)/sqrtdsl+usr(4)/sqrtdsr + (usr(7)-usl(7))*sbx)/(sumsqrtds), & 00510 (usl(5)/sqrtdsl+usr(5)/sqrtdsr + (usr(8)-usl(8))*sbx)/(sumsqrtds)/) 00511 Bss = (/Bx, & 00512 (sqrtdsl*usr(7)+sqrtdsr*usl(7)+(sqrtdsl/sqrtdsr*usr(4)-sqrtdsr/sqrtdsl*usl(4))*sbx)/ & 00513 sumsqrtds, & 00514 (sqrtdsl*usr(8)+sqrtdsr*usl(8)+(sqrtdsl/sqrtdsr*usr(5)-sqrtdsr/sqrtdsl*usl(5))*sbx)/ & 00515 sumsqrtds/) 00516 END IF 00517 00518 IF (SM >= 0.d0) THEN !Regions 2-3 00519 CALL F_prim_MHD(wl,fl) 00520 IF (SsL >= 0.d0) THEN !Region 2 00521 fl = fl+SL*(usl-ul) 00522 ELSE !Region 3 00523 ussl(1:8) = (/usl(1), & 00524 usl(2)-sbx*(DOT_PRODUCT(usl(3:5),usl(6:8))/sqrtdsl-sqrtdsl*DOT_PRODUCT(vss,Bss)), & 00525 usl(1)*vss(1), & 00526 usl(1)*vss(2), & 00527 usl(1)*vss(3), & 00528 Bss(1), & 00529 Bss(2), & 00530 Bss(3)/) 00531 ussl(9:n)=usl(9:n) 00532 fl = fl+SsL*ussl - (SsL-SL)*usl-SL*ul 00533 END IF 00534 f=fl 00535 END IF 00536 IF (SM <= 0.d0) THEN !Regions 4-5 00537 CALL F_prim_MHD(wr,fr) 00538 IF (SsR <= 0.d0) THEN !Region 5 00539 fr = fr+SR*(usr-ur) 00540 ELSE !Region 4 00541 ussr(1:8) = (/usr(1), & 00542 usr(2)+sbx*(DOT_PRODUCT(usr(3:5),usr(6:8))/sqrtdsr-sqrtdsr*DOT_PRODUCT(vss,Bss)), & 00543 usr(1)*vss(1), & 00544 usr(1)*vss(2), & 00545 usr(1)*vss(3), & 00546 Bss(1), & 00547 Bss(2), & 00548 Bss(3)/) 00549 ussr(9:n)=usr(9:n) 00550 fr = fr+SsR*ussr - (SsR-SR)*usr-SR*ur 00551 END IF 00552 f=fr 00553 END IF 00554 IF (SM == 0.d0) THEN 00555 f=.5d0*(fl+fr) 00556 END IF 00557 END IF 00558 f(6)=0d0 !kill normal magnetic flux 00559 HLLD = MAX(ABS(SL),ABS(SR))!, ABS(wleft(2))+CL, ABS(wright(2))+CR) 00560 DEALLOCATE(ul, ur, rr, rl, rm, usl, usr, ussl, ussr, fl, fr) 00561 00562 END FUNCTION HLLD 00563 00569 FUNCTION HLLC(wl,wr,f, lambda_max) 00570 REAL(KIND=qPREC) :: HLLC 00571 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00572 REAL(KIND=qPrec), DIMENSION(:) :: wl, wr, f 00573 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: ul, ur, rl, rr, rm, usl, usr, fl, fr 00574 REAL(KIND=qPrec) :: SL, SR, SM, sqrtdl,sqrtdr, CL, CR, pT, tempr, templ, a, dsl, dsr 00575 INTEGER :: n, x(3) 00576 n = SIZE(f) 00577 00578 ALLOCATE(ur(n), ul(n), rl(n), rr(n), rm(n), usl(n), usr(n), fl(n), fr(n)) 00579 CL=sound_speed(wl) 00580 CR=sound_speed(wr) 00581 00582 !Einfeldt version 00583 CALL prim_to_Roe(wl,rl) 00584 CALL prim_to_Roe(wr,rr) 00585 sqrtdl=SQRT(wl(1)) 00586 sqrtdr=SQRT(wr(1)) 00587 rm(2:)=(sqrtdl*rl(2:)+sqrtdr*rr(2:))/(sqrtdl+sqrtdr) 00588 rm(1)=sqrtdl*sqrtdr 00589 a=sound_speed_Roe(rm) 00590 SL = MIN(wl(3)-CL, rm(3)-a) 00591 SR = MAX(wr(3)+CR, rm(3)+a) 00592 00593 !Davis Version 00594 ! SL=min(wl(3)-CL,wr(3)-CR) 00595 ! SR=max(wl(3)+CL,wr(3)+CR) 00596 00597 !protection to ensure proper behavior when flow is supersonic 00598 ! IF (present(lambda_max)) THEN 00599 ! IF (lambda_max <. SR .OR. lambda_max < -SL) write(*,*) SL, SR, lambda_max 00600 ! END IF 00601 00602 ! H Viscosity protection 00603 ! Sanders 9a 00604 ! IF (present(lambda_max)) THEN 00605 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00606 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00607 ! END IF 00608 00609 ! Sanders 9b 00610 IF (present(lambda_max)) THEN 00611 SL=SIGN(abs(SL)+lambda_max, SL) 00612 SR=SIGN(abs(SR)+lambda_max, SR) 00613 END IF 00614 00615 00616 tempr = (SR-wr(3))*wr(1) 00617 templ = (SL-wl(3))*wl(1) 00618 SM = ((tempr*wr(3)-templ*wl(3))-(wr(2)-wl(2)))/(tempr - templ) 00619 pT = ((tempr*wl(2)-templ*wr(2))+tempr*templ*(wr(3)-wl(3)))/ & 00620 (tempr-templ) 00621 IF (SL >= 0d0) THEN 00622 CALL F_prim(wl, f) 00623 ELSE IF (SR <= 0d0) THEN 00624 CALL F_prim(wr, f) 00625 ELSE 00626 IF (SM >= 0d0) THEN 00627 CALL F_prim(wl, fl) 00628 CALL rprim_to_cons(wl,ul) 00629 dsl = templ/(SL-SM) 00630 usl(1:3)=(/dsl,((SL-wl(3))*ul(2)-wl(2)*wl(3)+pT*SM)/(SL-SM),dsl*SM/) 00631 usl(4:m_high+1)=dsl*wl(4:m_high+1) 00632 usl(m_high+2:n)=dsl*wl(m_high+2:n)/wl(1) 00633 fl = fl+SL*(usl-ul) 00634 f=fl 00635 END IF 00636 IF (SM <= 0d0) THEN 00637 CALL F_prim(wr,fr) 00638 CALL rprim_to_cons(wr,ur) 00639 dsr=tempr/(SR-SM) 00640 usr(1:3)=(/dsr,((SR-wr(3))*ur(2)-wr(2)*wr(3)+pT*SM)/(SR-SM),dsr*SM/) 00641 usr(4:m_high+1)=dsr*wr(4:m_high+1) 00642 usr(m_high+2:n)=dsr*wr(m_high+2:n)/wr(1) 00643 fr= fr + SR*(usr-ur) 00644 f=fr 00645 END IF 00646 IF (SM == 0d0) THEN 00647 f=.5d0*(fr+fl) 00648 END IF 00649 END IF 00650 ! IF (MAX(ABS(SL), ABS(SR)) > 90000) THEN 00651 ! write(*,*) wl, wr, CL, CR, SL, SR 00652 ! END IF 00653 HLLC = max(ABS(SL),ABS(SR)) 00654 DEALLOCATE(ur, ul, rl, rr, rm, usl, usr, fl, fr) 00655 END FUNCTION HLLC 00656 00662 00663 FUNCTION HLL(wleft,wright,f, lambda_max) 00664 REAL(KIND=qPREC) :: HLL 00665 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00666 REAL(KIND=qPrec), DIMENSION(:) :: wleft, wright, f 00667 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: uleft, uright, fl, fr, rl, rr, rm 00668 REAL(KIND=qPrec) :: SL, SR, sqrtdl,sqrtdr, CL, CR, a 00669 INTEGER :: n 00670 n=size(f) 00671 ALLOCATE(uleft(n), uright(n), fl(n), fr(n), rl(n), rr(n), rm(n)) 00672 CALL rprim_to_cons(wleft,uleft) 00673 CALL rprim_to_cons(wright,uright) 00674 sqrtdl=SQRT(wleft(1)) 00675 sqrtdr=SQRT(wright(1)) 00676 CL=sound_speed(wleft) 00677 CR=sound_speed(wright) 00678 00679 !Einfeldt version 00680 CALL prim_to_Roe(wleft,rl) 00681 CALL prim_to_Roe(wright,rr) 00682 sqrtdl=SQRT(wleft(1)) 00683 sqrtdr=SQRT(wright(1)) 00684 rm(2:)=(sqrtdl*rl(2:)+sqrtdr*rr(2:))/(sqrtdl+sqrtdr) 00685 rm(1)=sqrtdl*sqrtdr 00686 a=sound_speed_Roe(rm) 00687 SL = MIN(wleft(3)-CL, rm(3)-a) 00688 SR = MAX(wright(3)+CR, rm(3)+a) 00689 00690 !Davis Version 00691 ! SL=min(wleft(3)-CL,wright(3)-CR) 00692 ! SR=max(wleft(3)+CL,wright(3)+CR) 00693 00694 !protection to ensure proper behavior when flow is supersonic 00695 00696 ! Sanders 9a 00697 ! IF (present(lambda_max)) THEN 00698 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00699 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00700 ! END IF 00701 00702 ! Sanders 9b 00703 IF (present(lambda_max)) THEN 00704 SL=SIGN(abs(SL)+lambda_max, SL) 00705 SR=SIGN(abs(SR)+lambda_max, SR) 00706 END IF 00707 00708 IF (SL >= 0d0) THEN 00709 CALL F_prim(wleft, f) 00710 ELSE IF (SR <= 0d0) THEN 00711 CALL F_prim(wright,f) 00712 ELSE 00713 CALL rprim_to_cons(wleft,uleft) 00714 CALL rprim_to_cons(wright,uright) 00715 CALL F_prim(wleft, fl) 00716 CALL F_prim(wright, fr) 00717 f = ((SR*fl-SL*fr)+SR*SL*(uright-uleft))/(SR-SL) 00718 END IF 00719 HLL = MAX(ABS(SL),ABS(SR)) 00720 DEALLOCATE(uleft, uright, fl, fr, rl, rr, rm) 00721 END FUNCTION HLL 00722 00728 00729 FUNCTION HLL_MHD(wleft,wright,f, lambda_max) 00730 REAL(KIND=qPREC) :: HLL_MHD 00731 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00732 REAL(KIND=qPrec), DIMENSION(:) :: wleft, wright, f 00733 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: uleft, uright, fl, fr, rl, rr, rm 00734 REAL(KIND=qPrec) :: SL, SR, sqrtdl,sqrtdr, CL, CR, a 00735 INTEGER :: n 00736 n=size(f) 00737 ALLOCATE(uleft(n), uright(n), fl(n), fr(n), rl(n), rr(n), rm(n)) 00738 CALL rprim_to_cons_MHD(wleft,uleft) 00739 CALL rprim_to_cons_MHD(wright,uright) 00740 sqrtdl=SQRT(wleft(1)) 00741 sqrtdr=SQRT(wright(1)) 00742 CL=fast_speed(wleft) 00743 CR=fast_speed(wright) 00744 00745 !Einfeldt version 00746 CALL prim_to_Roe_MHD(wleft,rl) 00747 CALL prim_to_Roe_MHD(wright,rr) 00748 sqrtdl=SQRT(wleft(1)) 00749 sqrtdr=SQRT(wright(1)) 00750 rm(2:)=(sqrtdl*rl(2:)+sqrtdr*rr(2:))/(sqrtdl+sqrtdr) 00751 rm(1)=sqrtdl*sqrtdr 00752 a=fast_speed_Roe(rm) 00753 SL = MIN(wleft(3)-CL, rm(3)-a) 00754 SR = MAX(wright(3)+CR, rm(3)+a) 00755 00756 !Davis Version 00757 ! SL=min(wleft(3)-CL,wright(3)-CR) 00758 ! SR=max(wleft(3)+CL,wright(3)+CR) 00759 00760 !protection to ensure proper behavior when flow is supersonic 00761 00762 ! Sanders 9a 00763 ! IF (present(lambda_max)) THEN 00764 ! SL=SIGN(max(abs(SL),lambda_max), SL) 00765 ! SR=SIGN(max(abs(SR),lambda_max), SR) 00766 ! END IF 00767 00768 ! Sanders 9b 00769 IF (present(lambda_max)) THEN 00770 SL=SIGN(abs(SL)+lambda_max, SL) 00771 SR=SIGN(abs(SR)+lambda_max, SR) 00772 END IF 00773 00774 IF (SL >= 0d0) THEN 00775 CALL F_prim_MHD(wleft, f) 00776 ELSE IF (SR <= 0d0) THEN 00777 CALL F_prim_MHD(wright,f) 00778 ELSE 00779 ! CALL rprim_to_cons(wleft,uleft) 00780 ! CALL rprim_to_cons(wright,uright) 00781 CALL F_prim_MHD(wleft, fl) 00782 CALL F_prim_MHD(wright, fr) 00783 f = ((SR*fl-SL*fr)+SR*SL*(uright-uleft))/(SR-SL) 00784 END IF 00785 HLL_MHD = MAX(ABS(SL),ABS(SR)) 00786 DEALLOCATE(uleft, uright, fl, fr, rl, rr, rm) 00787 END FUNCTION HLL_MHD 00788 00794 FUNCTION ExactRS(wl, wr, f, lambda_max) 00795 REAL(KIND=qPREC) :: ExactRS,max_speed 00796 REAL(KIND=qPREC), OPTIONAL :: lambda_max 00797 REAL(KIND=qPrec), DIMENSION(:) :: wl, wr, f 00798 REAL(KIND=qPrec), DIMENSION(:), ALLOCATABLE :: wstar 00799 ALLOCATE (wstar(SIZE(wl))) 00800 CALL riemann(wl,wr,wstar,max_speed) 00801 CALL f_prim(wstar,f) 00802 ExactRS=max_speed 00803 DEALLOCATE(wstar) 00804 END FUNCTION ExactRS 00805 00806 00807 00812 SUBROUTINE riemann(WLeft,WRight,WMiddle,max_speed) 00813 REAL(KIND=qPrec), DIMENSION(:) :: WLeft, WRight, WMiddle 00814 REAL(KIND=qPrec), DIMENSION(3) :: Left, Right, Middle 00815 REAL(KIND=qPrec) :: S, UM,max_speed 00816 S=0.0 00817 Left(1:3)=WLeft(1:3) 00818 Right(1:3)=WRight(1:3) 00819 CALL vacuum_solve(Left, Right, Middle, UM, S,max_speed) 00820 Wmiddle(1:3)=Middle(1:3) 00821 IF (SIZE(WLeft) >= 4) THEN 00822 IF (WMiddle(3) >= 0) THEN 00823 WMiddle(4:)=WLeft(4:) 00824 ELSE 00825 WMiddle(4:)=WRight(4:) 00826 END IF 00827 END IF 00828 END SUBROUTINE riemann 00829 00836 SUBROUTINE vacuum_solve(Left, Right, Middle, UM, S, max_speed) 00837 USE GlobalDeclarations 00838 REAL(KIND=qPrec), DIMENSION(3) :: Left, Right, Middle, Vacuum 00839 REAL(KIND=qPrec) :: S, CL, CR, PM, UM, SHL, SHR, UML, UMR, C,max_speed 00840 INTEGER :: solver 00841 Vacuum = (/REAL(0, KIND=qPrec), REAL(0,KIND=qPrec), REAL(0,KIND=qPrec) /) 00842 Right(1)=MAX(Right(1),0d0) 00843 Left(1)=MAX(Left(1),0d0) 00844 Right(2)=MAX(Right(2),MinDensity*MinTemp) 00845 Left(2)=MAX(Left(2),MinDensity*MinTemp) 00846 00847 IF (Left(1) > MinDensity .AND. Right(1) > MinDensity) THEN !check for vacuum generation 00848 CL = SQRT(gamma*Left(2)/Left(1)) 00849 CR = SQRT(gamma*Right(2)/Right(1)) 00850 IF(gamma2*(CL+CR) <= Right(3)-Left(3)) THEN 00851 solver = 1 00852 SHR = Right(3)+CR 00853 SHL = Left(3)-CL 00854 UMR = Right(3)-gamma2*CR 00855 UML = Left(3)+gamma2*CL 00856 max_speed = MAX(ABS(SHR), ABS(SHL)) 00857 IF (S<=SHL) THEN 00858 Middle = Left 00859 ELSE IF (S <= UML) THEN 00860 C=gamma9*(CL+gamma12*(Left(3)-S)) 00861 Middle=(/Left(1)*(C/CL)**gamma2, Left(2)*(C/CL)**gamma10, gamma9*(CL+gamma12*Left(3)+S)/) 00862 ELSE IF (S <= UMR) THEN 00863 Middle = Vacuum 00864 ELSE IF (S < SHR) THEN 00865 C=gamma9*(CR-gamma12*(Right(3)-S)) 00866 Middle = (/Right(1)*(C/CR)**gamma2, Right(2)*(C/CR)**gamma10,gamma9*(-CR+gamma12*Right(3)+S) /) 00867 ELSE 00868 Middle = Right 00869 END IF 00870 ELSE 00871 solver = 2 00872 CALL STARPU(Left,Right, CL, CR, PM, UM) 00873 CALL SAMPLE(PM, UM, S, Left, Right, Middle, CL, CR, max_speed) 00874 END IF 00875 00876 ELSE IF (Right(1) > MinDensity) THEN 00877 solver = 3 00878 CR = SQRT(gamma*Right(2)/Right(1)) 00879 SHR = Right(3)+CR 00880 UM = Right(3)-gamma2*CR 00881 max_speed = MAX(abs(SHR), abs(UM)) 00882 IF (S <= UM) THEN ! Sampled point in vacuum state 00883 Middle = Vacuum 00884 ELSE IF (S < SHR) THEN ! Sampled point in right fan 00885 C=gamma9*(CR-gamma12*(Right(3)-S)) 00886 Middle = (/Right(1)*(C/CR)**gamma2, Right(2)*(C/CR)**gamma10,gamma9*(-CR+gamma12*Right(3)+S) /) 00887 ELSE !Sampled point in right state 00888 Middle = Right 00889 END IF 00890 00891 ELSE IF (Left(1) > MinDensity) THEN 00892 solver = 4 00893 CL = SQRT(gamma*Left(2)/Left(1)) 00894 SHL = Left(3)-CL 00895 UM = Left(3)+gamma2*CL 00896 max_speed = MAX(abs(SHL), abs(UM)) 00897 IF (S >= UM) THEN ! Sampled point in vacuum state 00898 Middle = Vacuum 00899 ELSE IF (S > SHL) THEN ! Sampled point in left fan 00900 C=gamma9*(CL+gamma12*(Left(3)-S)) 00901 Middle=(/Left(1)*(C/CL)**gamma2, Left(2)*(C/CL)**gamma10,gamma9*(CL+gamma12*Left(3)+S)/) 00902 ELSE !Sampled point in left state 00903 Middle = Left 00904 END IF 00905 00906 ELSE !both vacuum 00907 solver = 5 00908 Middle=Vacuum 00909 max_speed=0d0 00910 END IF 00911 00912 IF (Middle(1)+Middle(2)+Middle(3) < 1e20) THEN 00913 ELSE 00914 write(*,*) Left, Right, Middle, solver 00915 STOP 00916 END IF 00917 END SUBROUTINE vacuum_solve 00918 00920 SUBROUTINE STARPU(Left, Right, CL, CR, PM, UM) 00921 REAL(KIND=qPrec), DIMENSION(3) :: Left, Right 00922 REAL(KIND=qPrec) :: CL, CR, CHANGE, FL, FLD, FR, FRD, POLD, PSTART, UDIFF, UM, PM 00923 INTEGER I, NRITER 00924 NRITER=100 00925 CALL guessp(Left, Right, CL, CR, PSTART) 00926 POLD = PSTART 00927 UDIFF = Right(3)-Left(3) 00928 DO I=1,NRITER 00929 CALL PREFUN(FL, FLD, POLD, Left, CL) 00930 CALL PREFUN(FR, FRD, POLD, Right, CR) 00931 PM=POLD-(FL+FR+UDIFF)/(FLD+FRD) 00932 IF (PM <= 0.0) PM=MinDensity*MinTemp 00933 CHANGE=2.0*ABS((PM-POLD)/(PM+POLD)) 00934 IF (CHANGE <= 1e-10) EXIT 00935 IF (I == NRITER ) THEN 00936 write(*,*) I 00937 write(*,*) "Left=", Left 00938 write(*,*) "Right=", Right 00939 write(*,*) "PSTART=", PSTART 00940 write(*,*) "PM=", PM 00941 write(*,*) "POLD=", POLD 00942 write(*,*) "Change=", CHANGE 00943 write(*,*) "CL=", CL 00944 write(*,*) "CR=", CR 00945 write(*,*) "STARPU reached max iterations in riemann_solvers.f90" 00946 STOP 00947 END IF 00948 POLD = PM 00949 END DO 00950 UM=half*(Left(3)+Right(3)+FR-FL) 00951 END SUBROUTINE STARPU 00952 00954 SUBROUTINE SAMPLE(PM, UM, S, Left, Right, Middle, CL, CR, max_speed) 00955 REAL(KIND=qPrec), DIMENSION(3) :: Left, Right, Middle, MiddleR, MiddleL 00956 REAL(KIND=qPrec) :: DL, UL, PL, CL, DR, UR, PR, CR, P, D, U, S, UM, PM 00957 REAL(KIND=qPrec) :: C, CML, CMR, PML, PMR, SHL, SHR, SL, SR, STL, STR,Smax 00958 REAL(KIND=qPREC) :: max_speed 00959 SHL = Left(3) - CL 00960 SHR = Right(3) + CR 00961 PML = PM/Left(2) 00962 PMR = PM/Right(2) 00963 SL = Left(3) - CL*SQRT(gamma8*PML+gamma3) 00964 SR = Right(3) + CR*SQRT(gamma8*PMR+gamma3) 00965 max_speed=0d0 00966 IF (PM<=Left(2)) THEN !Left rarefaction 00967 max_speed=MAX(max_speed,ABS(SHL)) 00968 ELSE !Left shock 00969 max_speed=MAX(max_speed,ABS(SL)) 00970 END IF 00971 IF (PM<=Right(2)) THEN !Right rarefaction 00972 max_speed=MAX(max_speed,ABS(SHR)) 00973 ELSE !Right shock 00974 max_speed=MAX(max_speed,ABS(SR)) 00975 END IF 00976 max_speed=MAX(max_speed,ABS(UM)) !Entropy wave 00977 00978 IF (S <= UM) THEN !Sample left of contact 00979 IF (PM <= Left(2)) THEN !Left rarefaction 00980 IF (S <= SHL) THEN !Sample in left state 00981 Middlel = Left 00982 ELSE 00983 CML=CL*(PML)**gamma3 00984 STL=UM-CML 00985 IF (S > STL) THEN !Sample in star left state 00986 Middlel = (/Left(1)*(PML)**(gamma6), PM, UM /) 00987 ELSE !Sample in left fan 00988 C=gamma9*(CL+gamma12*(Left(3)-S)) 00989 Middlel=(/Left(1)*(C/CL)**gamma2, Left(2)*(C/CL)**gamma10, gamma9*(CL+gamma12*Left(3)+S)/) 00990 END IF 00991 END IF 00992 ELSE !Left shock 00993 IF (S <= SL) THEN !Sample in left state 00994 Middlel = Left 00995 ELSE !Sample in star left state 00996 Middlel = (/Left(1)*(PML+gamma11)/(PML*gamma11+1d0), PM, UM/) 00997 END IF 00998 END IF 00999 middle=middlel 01000 END IF 01001 IF (S >= UM) THEN !Sample right of contact 01002 IF (PM > Right(2)) THEN !Right shock 01003 IF (S >= SR) THEN !Sample in right state 01004 MiddleR = Right 01005 ELSE !Sample in star right state 01006 MiddleR = (/Right(1)*(PMR+gamma11)/(PMR*gamma11+1d0), PM, UM/) 01007 END IF 01008 ELSE !Right rarefaction 01009 IF (S >= SHR) THEN !Sample in right state 01010 MiddleR = Right 01011 ELSE 01012 CMR = CR*(PMR)**gamma3 01013 STR = UM + CMR 01014 IF (S <= STR) THEN !Sample in star right state 01015 MiddleR = (/Right(1)*(PMR)**(gamma6), PM, UM/) 01016 ELSE !Sample in right fan 01017 C=gamma9*(CR-gamma12*(Right(3)-S)) 01018 MiddleR = (/Right(1)*(C/CR)**gamma2,Right(2)*(C/CR)**gamma10, gamma9*(-CR+gamma12*Right(3)+S) /) 01019 END IF 01020 END IF 01021 END IF 01022 middle=middleR 01023 END IF 01024 IF (S == UM) THEN 01025 Middle=.5d0*(middleL+middleR) 01026 END IF 01027 END SUBROUTINE SAMPLE 01028 01030 SUBROUTINE guessp(Left,Right,CL,CR,PM) 01031 REAL(KIND=qPrec) :: PM, PQ, UM, GEL, GER, CL, CR, UDIFF, PPV, PMIN, PMAX, QMAX, QUSER, PTL, PTR 01032 REAL(KIND=qPrec), DIMENSION(3) :: Left, Right 01033 QUSER=2.0 01034 UDIFF=Right(3)-Left(3) 01035 PPV=half*(Left(2)+Right(2))-0.125*UDIFF*(Left(1)+Right(1))*(CL+CR) 01036 PMIN=MIN(Left(2),Right(2)) 01037 PMAX=MAX(Left(2),Right(2)) 01038 QMAX=PMAX/PMIN 01039 IF (QMAX <= QUSER .AND. (PMIN <=PPV .AND. PPV <= PMAX)) THEN !PVRS 01040 PM=PPV 01041 ELSE 01042 IF (PPV < PMIN) THEN !TRRS 01043 PQ=(Left(2)/Right(2))**gamma3 01044 UM=(PQ*Left(3)/CL+Right(3)/CR+gamma2*(PQ-1d0))/(PQ/CL+1d0/CR) 01045 PTL=1d0+gamma12*(Left(3)-UM)/CL 01046 PTR=1d0+gamma12*(UM-Right(3))/CR 01047 PM=half*(Left(2)*PTL**gamma10+Right(2)*PTR**gamma10) 01048 ELSE !TSRS 01049 GEL=SQRT((gamma9/Left(1))/(gamma11*Left(2)+PPV)) 01050 GER=SQRT((gamma9/Right(1))/(gamma11*Right(2)+PPV)) 01051 PM=(GEL*Left(2)+GER*Right(2)-UDIFF)/(GEL+GER) 01052 END IF 01053 END IF 01054 PM = MAX(PM,MinDensity*MinTemp) 01055 END SUBROUTINE guessp 01056 01058 SUBROUTINE PREFUN(F, FD, Pstar, WK, CK) 01059 REAL(KIND=qPrec) :: Pstar, F, FD, CK, PRAT, AK, BK, QRT 01060 REAL(KIND=qPrec), DIMENSION(3) :: WK 01061 IF (Pstar <= WK(2)) THEN 01062 PRAT=Pstar/WK(2) 01063 F=gamma2*CK*(PRAT**gamma3-1.0) 01064 FD=(1.0/(WK(1)*CK))*PRAT**(-gamma8) 01065 ELSE 01066 AK=gamma9/WK(1) 01067 BK=gamma11*WK(2) 01068 QRT=SQRT(AK/(BK+Pstar)) 01069 F=(Pstar-WK(2))*QRT 01070 FD=(1d0-half*(Pstar-WK(2))/(BK+Pstar))*QRT 01071 END IF 01072 END SUBROUTINE PREFUN 01073 01074 END MODULE RiemannSolvers 01075