aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_FluxSplit.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_FluxSplit.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_FluxSplit.F90')
-rw-r--r--src/GRHydro_FluxSplit.F90154
1 files changed, 40 insertions, 114 deletions
diff --git a/src/GRHydro_FluxSplit.F90 b/src/GRHydro_FluxSplit.F90
index c5fe6d7..d54e446 100644
--- a/src/GRHydro_FluxSplit.F90
+++ b/src/GRHydro_FluxSplit.F90
@@ -11,6 +11,7 @@
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
#include "cctk_Functions.h"
+#include "GRHydro_Macros.h"
/*@@
@routine GRHydro_FSAlpha
@@ -66,22 +67,11 @@ subroutine GRHydro_FSAlpha(CCTK_ARGUMENTS)
do j = 1, ny
do i = 1, nx
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- 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))
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
if (shift_state .eq. 0) then
beta = 0.d0
@@ -123,22 +113,11 @@ subroutine GRHydro_FSAlpha(CCTK_ARGUMENTS)
do j = 1, ny
do i = 1, nx
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- 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))
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
if (shift_state .eq. 0) then
beta = 0.d0
@@ -180,22 +159,11 @@ subroutine GRHydro_FSAlpha(CCTK_ARGUMENTS)
do j = 1, ny
do i = 1, nx
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
- 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))
- 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det,&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
if (shift_state .eq. 0) then
beta = 0.d0
@@ -347,28 +315,14 @@ subroutine GRHydro_SplitFlux(CCTK_ARGUMENTS)
else
dummy = betax(:,j,k)
end if
- if (conformal_state .eq. 0) then
- do i = 1, cctk_lsh(1)
- 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(i))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(i),&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(i,j,k))
- upper(i) = uxx
- end do
- else
- do i = 1, cctk_lsh(1)
- 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(i))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(i),&
- 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))
- upper(i) = uxx
- end do
- end if
+ do i = 1, cctk_lsh(1)
+ 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(i),&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
+ upper(i) = uxx
+ end do
call GRHydro_SplitFlux_1D(GRHydro_eos_handle, nx,&
fs_alpha1, fs_alpha2, fs_alpha3, fs_alpha4, fs_alpha5, &
@@ -398,28 +352,14 @@ subroutine GRHydro_SplitFlux(CCTK_ARGUMENTS)
else
dummy = betay(i,:,k)
end if
- if (conformal_state .eq. 0) then
- do j = 1, cctk_lsh(2)
- 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(j))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(j),&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(i,j,k))
- upper(j) = uyy
- end do
- else
- do j = 1, cctk_lsh(2)
- 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(j))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(j),&
- 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))
- upper(j) = uyy
- end do
- end if
+ do j = 1, cctk_lsh(2)
+ 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(j),&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
+ upper(j) = uyy
+ end do
call GRHydro_SplitFlux_1D(GRHydro_eos_handle, ny,&
fs_alpha1, fs_alpha2, fs_alpha3, fs_alpha4, fs_alpha5, &
@@ -450,28 +390,14 @@ subroutine GRHydro_SplitFlux(CCTK_ARGUMENTS)
dummy = betaz(i,j,:)
end if
- if (conformal_state .eq. 0) then
- do k = 1, cctk_lsh(3)
- 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(k))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(k),&
- gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
- gyz(i,j,k),gzz(i,j,k))
- upper(k) = uzz
- end do
- else
- do k = 1, cctk_lsh(3)
- 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(k))
- call UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(k),&
- 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))
- upper(k) = uzz
- end do
- end if
+ do k = 1, cctk_lsh(3)
+ 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 UpperMetric(uxx,uxy,uxz,uyy,uyz,uzz,det(k),&
+ gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),gyy(i,j,k),&
+ gyz(i,j,k),gzz(i,j,k))
+ upper(k) = uzz
+ end do
call GRHydro_SplitFlux_1D(GRHydro_eos_handle, nz,&
fs_alpha1, fs_alpha2, fs_alpha3, fs_alpha4, fs_alpha5, &