aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authoreschnett <eschnett@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-12-23 02:47:41 +0000
committereschnett <eschnett@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-12-23 02:47:41 +0000
commitb0598b903ba8f3bb5ed89c29e0807cf4cb108110 (patch)
treec4a97fe09279c3c1b573e4299231b2b2671cf72b /src
parent115eaef5522b79393771474f761edb8044154a75 (diff)
Loop in most prim2con routines over all grid points, not only the
evolved ones. In this way, the conserved variables are defined on all grid points. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@198 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/GRHydro_Prim2Con.F9036
1 files changed, 18 insertions, 18 deletions
diff --git a/src/GRHydro_Prim2Con.F90 b/src/GRHydro_Prim2Con.F90
index e142798..048e6ae 100644
--- a/src/GRHydro_Prim2Con.F90
+++ b/src/GRHydro_Prim2Con.F90
@@ -318,9 +318,9 @@ subroutine Primitive2ConservativeCells(CCTK_ARGUMENTS)
if(evolve_temper.ne.1) then
!$OMP PARALLEL DO PRIVATE(i, j, det)
- do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
- do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
- do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ 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))
@@ -337,9 +337,9 @@ subroutine Primitive2ConservativeCells(CCTK_ARGUMENTS)
!$OMP END PARALLEL DO
else
!$OMP PARALLEL DO PRIVATE(i, j, det)
- do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
- do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
- do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ 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))
@@ -360,9 +360,9 @@ subroutine Primitive2ConservativeCells(CCTK_ARGUMENTS)
if(evolve_Y_e.ne.0) then
!$OMP PARALLEL DO PRIVATE(i, j)
- do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
- do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
- do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ do i = 1,cctk_lsh(1)
Y_e_con(i,j,k) = Y_e(i,j,k) * dens(i,j,k)
enddo
enddo
@@ -551,9 +551,9 @@ subroutine Primitive2ConservativePolyCells(CCTK_ARGUMENTS)
CCTK_REAL :: det
!$OMP PARALLEL DO PRIVATE(i, j, det)
- do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
- do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
- do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ 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))
@@ -571,9 +571,9 @@ subroutine Primitive2ConservativePolyCells(CCTK_ARGUMENTS)
if(evolve_Y_e.ne.0) then
!$OMP PARALLEL DO PRIVATE(i, j)
- do k = GRHydro_stencil,cctk_lsh(3)-GRHydro_stencil+1
- do j = GRHydro_stencil,cctk_lsh(2)-GRHydro_stencil+1
- do i = GRHydro_stencil,cctk_lsh(1)-GRHydro_stencil+1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ do i = 1,cctk_lsh(1)
Y_e_con(i,j,k) = Y_e(i,j,k) * dens(i,j,k)
enddo
enddo
@@ -764,9 +764,9 @@ subroutine Primitive2ConservativeCellsGeneral(CCTK_ARGUMENTS)
ierr = EOS_SetGFs(cctkGH, EOS_Prim2ConCellsCall)
- do k = GRHydro_stencil, cctk_lsh(3) - GRHydro_stencil + 1
- do j = GRHydro_stencil, cctk_lsh(2) - GRHydro_stencil + 1
- do i = GRHydro_stencil, cctk_lsh(1) - GRHydro_stencil + 1
+ do k = 1,cctk_lsh(3)
+ do j = 1,cctk_lsh(2)
+ do i = 1,cctk_lsh(1)
gxxpt = gxx(i,j,k)
gxypt = gxy(i,j,k)