aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Con2PrimM_pt_EOSOmni.c
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2012-05-13 19:12:00 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2012-05-13 19:12:00 +0000
commit39ae36d214aba68151dc8ac45faf312c8bc988e1 (patch)
tree230daabf855d4bb19a28331c7f39e0fe2429e08a /src/GRHydro_Con2PrimM_pt_EOSOmni.c
parent73c3adae7b5e7e50f1e5f83caedc15ada6fb55b8 (diff)
Adding Josh's latest changes to the general EOS Con2PrimM routine.
patch by Philipp Moesta git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@324 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_Con2PrimM_pt_EOSOmni.c')
-rw-r--r--src/GRHydro_Con2PrimM_pt_EOSOmni.c10
1 files changed, 5 insertions, 5 deletions
diff --git a/src/GRHydro_Con2PrimM_pt_EOSOmni.c b/src/GRHydro_Con2PrimM_pt_EOSOmni.c
index 051db33..f0b6a39 100644
--- a/src/GRHydro_Con2PrimM_pt_EOSOmni.c
+++ b/src/GRHydro_Con2PrimM_pt_EOSOmni.c
@@ -89,7 +89,7 @@ static CCTK_INT twod_newton_raphson( CCTK_REAL x[], CCTK_REAL gammaeos, struct L
CCTK_REAL [], CCTK_REAL [][2],
CCTK_REAL *, CCTK_REAL *, CCTK_REAL, struct LocGlob *) );
-static CCTK_INT threed_newton_raphson( CCTK_REAL x[], struct eosomnivars *eosvars, struct LocGlob *lgp,
+static CCTK_INT threed_newton_raphson_omni( CCTK_REAL x[], struct eosomnivars *eosvars, struct LocGlob *lgp,
void (*funcd) (CCTK_REAL [], CCTK_REAL [],
CCTK_REAL [], CCTK_REAL [][3],
CCTK_REAL *, CCTK_REAL *,
@@ -386,7 +386,7 @@ void CCTK_FCALL CCTK_FNAME(GRHydro_Con2PrimM_pt) (
} else {
//USE 3d NR for non-polytropes!
x_3d[2] = u;
- *retval = 1.0*threed_newton_raphson( x_3d, &eosvars, &lg, func_vsq_eosomni ) ;
+ *retval = 1.0*threed_newton_raphson_omni( x_3d, &eosvars, &lg, func_vsq_eosomni ) ;
}
W = x_3d[0];
@@ -646,14 +646,14 @@ static CCTK_INT twod_newton_raphson( CCTK_REAL x[], CCTK_REAL gammaeos, struct L
/************************************************************
- threed_newton_raphson():
+ threed_newton_raphson_omni():
-- performs Newton-Rapshon method on an 2d system for polytropes.
-- inspired in part by Num. Rec.'s routine newt();
*****************************************************************/
-static CCTK_INT threed_newton_raphson( CCTK_REAL x[], struct eosomnivars *eosvars, struct LocGlob *lgp,
+static CCTK_INT threed_newton_raphson_omni( CCTK_REAL x[], struct eosomnivars *eosvars, struct LocGlob *lgp,
void (*funcd) (CCTK_REAL [], CCTK_REAL [], CCTK_REAL [],
CCTK_REAL [][3], CCTK_REAL *,
CCTK_REAL *, struct eosomnivars *, struct LocGlob *) )
@@ -900,7 +900,7 @@ static void func_vsq_eosomni(CCTK_REAL x[], CCTK_REAL dx[], CCTK_REAL resid[],
j01 = -1.0*B2plusW_sq;
- //Eq. 5: Qdotn + B^2(1+vsq)/2 - QdotBsq/2/W^2 + vsq W + rho_0(vsq) + u = 0
+ //Eq. 5: -Qdotn - B^2(1+vsq)/2 + QdotBsq/2/W^2 - vsq W - rho_0(vsq) - u = 0
//rho0 = D * sqrt(1-vsq)
resid[1] = -lgp->Qdotn - lgp->half_Bsq*(1.0+vsq) + 0.5*QB2Winv2 - vsq*W - rho0 - u;