aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Con2PrimM.F90
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-05-11 06:29:34 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-05-11 06:29:34 +0000
commit007c18aa5e6aed27c7edbf71b1f03267c0ffef9d (patch)
tree6d9c4331bcf20e8fc4f7fc2ab4fbb6802b2293a0 /src/GRHydro_Con2PrimM.F90
parentd73ce95f867e1a041cc51b82a954b6578e24610c (diff)
GRHydro: Fixing Con2PrimM to call the proper point-wise routines
This re-introduces routines that work for hybrid/hot EOS and corresponding changes in pointwise routines for hot EOS error checking and temperature treatment by adding old EOSOmni pointwise routine. From: Philipp Moesta <pmoesta@tapir.caltech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@511 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_Con2PrimM.F90')
-rw-r--r--src/GRHydro_Con2PrimM.F9026
1 files changed, 18 insertions, 8 deletions
diff --git a/src/GRHydro_Con2PrimM.F90 b/src/GRHydro_Con2PrimM.F90
index 0960b0d..720e826 100644
--- a/src/GRHydro_Con2PrimM.F90
+++ b/src/GRHydro_Con2PrimM.F90
@@ -339,8 +339,19 @@ subroutine Conservative2PrimitiveM(CCTK_ARGUMENTS)
!of the comoving B-field, b^{\mu} b_{\mu}. It is overwritten
!in this routine, but we may need to find a better notation
!avoid future confusions.
- call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, keytemp, &
- GRHydro_eos_rf_prec, local_gam(1), dens(i,j,k), &
+ if(GRHydro_eos_handle .eq. 1 .or. GRHydro_eos_handle .eq. 2) then
+
+ call GRHydro_Con2PrimM_ptold(GRHydro_eos_handle, local_gam(1), dens(i,j,k), &
+ scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3), tau(i,j,k), &
+ Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),rho_tmp, &
+ velx_tmp,vely_tmp,velz_tmp,eps_tmp,press_tmp, &
+ Bvecx_tmp,Bvecy_tmp,Bvecz_tmp,b2, w_lorentz_tmp,&
+ g11(i,j,k),g12(i,j,k),g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k), &
+ uxx,uxy,uxz,uyy,uyz,uzz,det, &
+ epsnegative,GRHydro_C2P_failed(i,j,k))
+
+ else
+ call GRHydro_Con2PrimM_pt2(GRHydro_eos_handle, keytemp, GRHydro_eos_rf_prec, GRHydro_perc_ptol, local_gam(1), dens(i,j,k), &
scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3), tau(i,j,k), &
Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3),xye(1), &
xtemp(1),rho_tmp,velx_tmp,vely_tmp,velz_tmp,&
@@ -349,6 +360,7 @@ subroutine Conservative2PrimitiveM(CCTK_ARGUMENTS)
g22(i,j,k),g23(i,j,k),g33(i,j,k), &
uxx,uxy,uxz,uyy,uyz,uzz,det, &
epsnegative,GRHydro_C2P_failed(i,j,k))
+ endif
if(evolve_entropy.ne.0) then
if(GRHydro_C2P_failed(i,j,k).ne.0) then
@@ -492,8 +504,7 @@ subroutine Conservative2PrimitiveM(CCTK_ARGUMENTS)
keytemp = 0
- call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, keytemp, &
- GRHydro_eos_rf_prec, local_gam(1), dens(i,j,k), &
+ call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, GRHydro_reflevel, i ,j ,k, x(i,j,k), y(i,j,k), z(i,j,k), keytemp, GRHydro_eos_rf_prec, GRHydro_perc_ptol, local_gam(1), dens(i,j,k), &
scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3), tau(i,j,k), &
Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3), Y_e(i,j,k), &
temperature(i,j,k),rho_tmp,velx_tmp,vely_tmp,velz_tmp,&
@@ -519,8 +530,7 @@ subroutine Conservative2PrimitiveM(CCTK_ARGUMENTS)
Bvecx_tmp = Bprim(i,j,k,1)
Bvecy_tmp = Bprim(i,j,k,2)
Bvecz_tmp = Bprim(i,j,k,3)
- call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, keytemp, &
- local_perc_ptol, local_gam(1), dens(i,j,k), &
+ call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, GRHydro_reflevel, i ,j ,k, x(i,j,k), y(i,j,k), z(i,j,k), keytemp, GRHydro_eos_rf_prec, GRHydro_perc_ptol, local_gam(1), dens(i,j,k), &
scon(i,j,k,1),scon(i,j,k,2),scon(i,j,k,3), tau(i,j,k), &
Bcons(i,j,k,1),Bcons(i,j,k,2),Bcons(i,j,k,3), Y_e(i,j,k), &
temperature(i,j,k),rho_tmp,velx_tmp,vely_tmp,velz_tmp,&
@@ -864,7 +874,7 @@ subroutine Conservative2PrimitiveBoundsM(CCTK_ARGUMENTS)
Y_e(i,j,k) = Y_e_con(i,j,k) / dens(i,j,k)
endif
- call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, keytemp, GRHydro_eos_rf_prec, local_gam(1),densminus(i,j,k), &
+ call GRHydro_Con2PrimM_pt2(GRHydro_eos_handle, keytemp, GRHydro_eos_rf_prec, GRHydro_perc_ptol, local_gam(1), densminus(i,j,k), &
sxminus(i,j,k),syminus(i,j,k),szminus(i,j,k), tauminus(i,j,k), &
Bconsxminus(i,j,k), Bconsyminus(i,j,k), Bconszminus(i,j,k), xye(1), xtemp(1), rhominus(i,j,k),&
velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k),epsminus(i,j,k),pressminus(i,j,k),&
@@ -922,7 +932,7 @@ subroutine Conservative2PrimitiveBoundsM(CCTK_ARGUMENTS)
epsnegative = 0
- call GRHydro_Con2PrimM_pt(GRHydro_eos_handle, keytemp, GRHydro_eos_rf_prec, local_gam(1),densplus(i,j,k), &
+ call GRHydro_Con2PrimM_pt2(GRHydro_eos_handle, keytemp, GRHydro_eos_rf_prec, GRHydro_perc_ptol, local_gam(1), densplus(i,j,k), &
sxplus(i,j,k),syplus(i,j,k),szplus(i,j,k), tauplus(i,j,k), &
Bconsxplus(i,j,k), Bconsyplus(i,j,k), Bconszplus(i,j,k), xye(1), xtemp(1), rhoplus(i,j,k),&
velxplus(i,j,k),velyplus(i,j,k),velzplus(i,j,k),epsplus(i,j,k),pressplus(i,j,k),&