aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Con2PrimM_pt.c
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:12 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:12 +0000
commit669d5b37c19d16380ba6c63a7195cbcea8d08f35 (patch)
tree43c35e9b941be21d80ba34711d3bbfe7825b3b6d /src/GRHydro_Con2PrimM_pt.c
parentfcbaf698714e9594ef35fbba18530a4a1ab164ba (diff)
GRHydro: add grid function for sqrt(detg)
* add new 1 tl grid function sdetg that stores the sqrt of the determinent of the 3-metric. * replace lots of re-computation of det by use of this grid function From: Christian Ott <cott@bethe.tapir.caltech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@555 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_Con2PrimM_pt.c')
-rw-r--r--src/GRHydro_Con2PrimM_pt.c9
1 files changed, 4 insertions, 5 deletions
diff --git a/src/GRHydro_Con2PrimM_pt.c b/src/GRHydro_Con2PrimM_pt.c
index bf947d3..06edca8 100644
--- a/src/GRHydro_Con2PrimM_pt.c
+++ b/src/GRHydro_Con2PrimM_pt.c
@@ -119,7 +119,7 @@ void CCTK_FCALL CCTK_FNAME(GRHydro_Con2PrimM_ptold) (
CCTK_REAL *gyy, CCTK_REAL *gyz, CCTK_REAL *gzz,
CCTK_REAL *uxx, CCTK_REAL *uxy, CCTK_REAL *uxz,
CCTK_REAL *uyy, CCTK_REAL *uyz, CCTK_REAL *uzz,
- CCTK_REAL *det,
+ CCTK_REAL *sdet,
CCTK_INT *epsnegative,
CCTK_REAL *retval);
@@ -201,7 +201,7 @@ void CCTK_FCALL CCTK_FNAME(GRHydro_Con2PrimM_ptold) (
CCTK_REAL *gyy, CCTK_REAL *gyz, CCTK_REAL *gzz,
CCTK_REAL *uxx, CCTK_REAL *uxy, CCTK_REAL *uxz,
CCTK_REAL *uyy, CCTK_REAL *uyz, CCTK_REAL *uzz,
- CCTK_REAL *det,
+ CCTK_REAL *sdet,
CCTK_INT *epsnegative,
CCTK_REAL *retval)
@@ -213,8 +213,7 @@ void CCTK_FCALL CCTK_FNAME(GRHydro_Con2PrimM_ptold) (
CCTK_REAL QdotB;
CCTK_REAL rho0,u,p,w,gammasq,gamma,gtmp,W_last,W,vsq;
CCTK_REAL g_o_WBsq, QdB_o_W;
- CCTK_REAL detg = (*det);
- CCTK_REAL sqrt_detg = sqrt(detg);
+ CCTK_REAL sqrt_detg = *sdet;
CCTK_REAL inv_sqrt_detg = 1./sqrt_detg;
CCTK_INT i_increase;
@@ -260,7 +259,7 @@ void CCTK_FCALL CCTK_FNAME(GRHydro_Con2PrimM_ptold) (
fprintf(stdout," *uyy = %26.16e \n", *uyy );
fprintf(stdout," *uyz = %26.16e \n", *uyz );
fprintf(stdout," *uzz = %26.16e \n", *uzz );
- fprintf(stdout," *det = %26.16e \n", *det );
+ fprintf(stdout," *sdet = %26.16e \n", *sdet );
fprintf(stdout," *epsnegative = %10d \n", *epsnegative );
fprintf(stdout," *retval = %26.16e \n", *retval );
fflush(stdout);