diff options
Diffstat (limited to 'src/GRHydro_TVDReconstruct_drv.F90')
-rw-r--r-- | src/GRHydro_TVDReconstruct_drv.F90 | 45 |
1 files changed, 18 insertions, 27 deletions
diff --git a/src/GRHydro_TVDReconstruct_drv.F90 b/src/GRHydro_TVDReconstruct_drv.F90 index ae220ea..880e35e 100644 --- a/src/GRHydro_TVDReconstruct_drv.F90 +++ b/src/GRHydro_TVDReconstruct_drv.F90 @@ -52,6 +52,23 @@ subroutine GRHydro_TVDReconstruct_drv(CCTK_ARGUMENTS) DECLARE_CCTK_PARAMETERS DECLARE_CCTK_FUNCTIONS + ! save memory when MP is not used + CCTK_INT :: GRHydro_UseGeneralCoordinates + + if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then + call TVDReconstruct_drv(gaa, gab, gac, gbb, gbc, gcc, lvel, lBvec) + else + call TVDReconstruct_drv(gxx, gxy, gxz, gyy, gyz, gzz, vel, Bvec) + end if + +contains +subroutine TVDReconstruct_drv(g11, g12, g13, g22, g23, g33, vup, Bprim) + + implicit none + + CCTK_REAL, DIMENSION(cctk_ash1,cctk_ash2,cctk_ash3) :: g11, g12, g13, g22, g23, g33 + CCTK_REAL, DIMENSION(cctk_ash1,cctk_ash2,cctk_ash3,3) :: vup, Bprim + integer :: nx, ny, nz, i, j, k, itracer logical, dimension(:,:,:), allocatable :: trivial_rp @@ -67,33 +84,6 @@ subroutine GRHydro_TVDReconstruct_drv(CCTK_ARGUMENTS) CCTK_REAL, DIMENSION(:,:,:,:), ALLOCATABLE :: wvel ! vel*w_lorentz CCTK_REAL :: agxx,agxy,agxz,agyy,agyz,agzz ! metric components - ! save memory when MP is not used - CCTK_INT :: GRHydro_UseGeneralCoordinates - CCTK_REAL, DIMENSION(cctk_ash1,cctk_ash2,cctk_ash3,3) :: vup, Bprim - pointer (pvup,vup), (pBprim,Bprim) - CCTK_REAL, DIMENSION(cctk_ash1,cctk_ash2,cctk_ash3) :: g11, g12, g13, g22, g23, g33 - pointer (pg11,g11), (pg12,g12), (pg13,g13), (pg22,g22), (pg23,g23), (pg33,g33) - - if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then - pg11 = loc(gaa) - pg12 = loc(gab) - pg13 = loc(gac) - pg22 = loc(gbb) - pg23 = loc(gbc) - pg33 = loc(gcc) - pvup = loc(lvel) - pBprim = loc(lBvec) - else - pg11 = loc(gxx) - pg12 = loc(gxy) - pg13 = loc(gxz) - pg22 = loc(gyy) - pg23 = loc(gyz) - pg33 = loc(gzz) - pvup = loc(vel) - pBprim = loc(Bvec) - end if - allocate(trivial_rp(cctk_lsh(1),cctk_lsh(2),cctk_lsh(3)),STAT=ierr) if (ierr .ne. 0) then call CCTK_WARN(0, "Allocation problems with trivial_rp") @@ -317,6 +307,7 @@ subroutine GRHydro_TVDReconstruct_drv(CCTK_ARGUMENTS) deallocate(trivial_rp) ! Fortran 90 will automatically deallocate wvel for us +end subroutine end subroutine GRHydro_TVDReconstruct_drv |