aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-10-14 08:35:29 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2011-10-14 08:35:29 +0000
commit17788d73436659db4c7731458edfc5c73b1ee4f1 (patch)
treefafc9ca154112702f28b34764e3eaf04224e8665 /src
parent7a689eb93032bc54ddd3d7935188baef49b0f759 (diff)
correct formula used to compute eps in EoSChangeK
this only shows up if we use a general eos (not a polytype) during evolution git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@282 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/GRHydro_EoSChangeGamma.F902
-rw-r--r--src/GRHydro_Prim2Con.F906
2 files changed, 4 insertions, 4 deletions
diff --git a/src/GRHydro_EoSChangeGamma.F90 b/src/GRHydro_EoSChangeGamma.F90
index 432b469..2820c1b 100644
--- a/src/GRHydro_EoSChangeGamma.F90
+++ b/src/GRHydro_EoSChangeGamma.F90
@@ -195,7 +195,7 @@ subroutine GRHydro_EoSChangeK(CCTK_ARGUMENTS)
end if
press = local_k * rho**local_gamma
- eps = (local_gamma - 1.d0) * local_k * rho**local_gamma
+ eps = local_k / (local_gamma - 1.d0) * rho**(local_gamma-1)
call Primitive2ConservativePolyCells(CCTK_ARGUMENTS)
diff --git a/src/GRHydro_Prim2Con.F90 b/src/GRHydro_Prim2Con.F90
index b5e1a22..3915576 100644
--- a/src/GRHydro_Prim2Con.F90
+++ b/src/GRHydro_Prim2Con.F90
@@ -541,14 +541,14 @@ subroutine Prim2ConservativePolytype(CCTK_ARGUMENTS)
avg_detl = SPATIAL_DETERMINANT(g11l,g12l,g13l,g22l, g23l,g33l)
avg_detr = SPATIAL_DETERMINANT(g11r,g12r,g13r,g22r, g23r,g33r)
- call prim2conpolytype(GRHydro_eos_handle, g11l,g12l,g13l,&
+ call prim2conpolytype(GRHydro_polytrope_handle, g11l,g12l,g13l,&
g22l,g23l,g33l, &
avg_detl,densminus(i,j,k),sxminus(i,j,k),&
syminus(i,j,k),szminus(i,j,k),tauminus(i,j,k),rhominus(i,j,k), &
velxminus(i,j,k),velyminus(i,j,k),velzminus(i,j,k),&
epsminus(i,j,k),pressminus(i,j,k),w_lorentzminus(i, j, k))
- call prim2conpolytype(GRHydro_eos_handle, g11r,g12r,g13r,&
+ call prim2conpolytype(GRHydro_polytrope_handle, g11r,g12r,g13r,&
g22r,g23r,g33r, &
avg_detr, densplus(i,j,k),sxplus(i,j,k),&
syplus(i,j,k),szplus(i,j ,k),tauplus(i,j,k),&
@@ -705,7 +705,7 @@ subroutine Primitive2ConservativePolyCells(CCTK_ARGUMENTS)
det = SPATIAL_DETERMINANT(g11(i,j,k),g12(i,j,k),g13(i,j,k),\
g22(i,j,k),g23(i,j,k),g33(i,j,k))
- call prim2conpolytype(GRHydro_eos_handle,g11(i,j,k),g12(i,j,k),&
+ call prim2conpolytype(GRHydro_polytrope_handle,g11(i,j,k),g12(i,j,k),&
g13(i,j,k),g22(i,j,k),g23(i,j,k),g33(i,j,k),&
det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
tau(i,j,k),rho(i,j,k),velx(i,j,k),vely(i,j,k),velz(i,j,k),&