aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_UpdateMask.F90
diff options
context:
space:
mode:
authorcott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-08-27 20:30:51 +0000
committercott <cott@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-08-27 20:30:51 +0000
commit39d5f5e568fecb8cf28f7a13418c472b14a85966 (patch)
tree702d5448cecb79d96d7fbc0ef758a493de569d98 /src/GRHydro_UpdateMask.F90
parentb9d5cef4e0c1a57d0b83961350d68b556f72e2c1 (diff)
* remove dependence on StaticConformal
* change calculation of the determinant of the 3-metric from a subroutine call to a macro. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@152 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_UpdateMask.F90')
-rw-r--r--src/GRHydro_UpdateMask.F9030
1 files changed, 9 insertions, 21 deletions
diff --git a/src/GRHydro_UpdateMask.F90 b/src/GRHydro_UpdateMask.F90
index 3753797..d7cdff0 100644
--- a/src/GRHydro_UpdateMask.F90
+++ b/src/GRHydro_UpdateMask.F90
@@ -11,6 +11,7 @@
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
+#include "GRHydro_Macros.h"
#include "SpaceMask.h"
#define velx(i,j,k) vel(i,j,k,1)
@@ -209,27 +210,14 @@ subroutine GRHydro_AtmosphereReset(CCTK_ARGUMENTS)
velx(i,j,k) = 0.0d0
vely(i,j,k) = 0.0d0
velz(i,j,k) = 0.0d0
- if (conformal_state .eq. 0) then
- call SpatialDeterminant(gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
- gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det)
- call prim2conpolytype(GRHydro_polytrope_handle, &
- gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
- gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det, &
- dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
- tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
- velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
- else
- psi4pt = psi(i,j,k)**4
- call SpatialDeterminant(psi4pt*gxx(i,j,k), psi4pt*gxy(i,j,k), &
- psi4pt*gxz(i,j,k), psi4pt*gyy(i,j,k), psi4pt*gyz(i,j,k), &
- psi4pt*gzz(i,j,k), det)
- call prim2conpolytype(GRHydro_polytrope_handle, &
- psi4pt*gxx(i,j,k), psi4pt*gxy(i,j,k), psi4pt*gxz(i,j,k), &
- psi4pt*gyy(i,j,k), psi4pt*gyz(i,j,k), psi4pt*gzz(i,j,k), &
- det, dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
- tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
- velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
- end if
+ det = SPATIAL_DETERMINANT(gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), \
+ gyy(i,j,k), gyz(i,j,k), gzz(i,j,k))
+ call prim2conpolytype(GRHydro_polytrope_handle, &
+ gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
+ gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det, &
+ dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
+ tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
+ velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
if (wk_atmosphere .eq. 0) then
atmosphere_mask(i, j, k) = 0
SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits,\