aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_PPMReconstruct_drv_opt.F90
diff options
context:
space:
mode:
authorrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:31 +0000
committerrhaas <rhaas@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2013-07-06 18:11:31 +0000
commit2d1ec54a5d8de927693b8a8f92d2bada3214a84e (patch)
treece6599c3ced1079a480209f1cf0e6f5d37b6e785 /src/GRHydro_PPMReconstruct_drv_opt.F90
parent4112b118acdcffdad4590b4aba56674fa37bc6ed (diff)
GRHydro: reimplement PPM in C
* add prototype of a re-implementation of PPM in C. So far we can do old PPM without and with temperature/Y_e reconstruction. * also add a version of the reconstruction driver that is optimized for the use of the new PPM routine. * no change to standard code version; one needs to manually replace source files to get the new stuff. * new version is factor 3 faster than old version (based on preliminary measurements). From: Christian Ott <cott@tapir.caltech.edu> git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@560 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src/GRHydro_PPMReconstruct_drv_opt.F90')
-rw-r--r--src/GRHydro_PPMReconstruct_drv_opt.F90212
1 files changed, 212 insertions, 0 deletions
diff --git a/src/GRHydro_PPMReconstruct_drv_opt.F90 b/src/GRHydro_PPMReconstruct_drv_opt.F90
new file mode 100644
index 0000000..4222965
--- /dev/null
+++ b/src/GRHydro_PPMReconstruct_drv_opt.F90
@@ -0,0 +1,212 @@
+ /*@@
+ @file GRHydro_PPMReconstruct_drv.F90
+ @date Tue Jul 19 13:22:03 EDT 2011
+ @author Bruno C. Mundim, Joshua Faber, Christian D. Ott
+ @desc
+ Driver routine to perform the PPM reconstructions.
+ @enddesc
+ @@*/
+
+#include "cctk.h"
+#include "cctk_Arguments.h"
+#include "cctk_Parameters.h"
+#include "cctk_Functions.h"
+
+#include "SpaceMask.h"
+
+#define velx(i,j,k) vup(i,j,k,1)
+#define vely(i,j,k) vup(i,j,k,2)
+#define velz(i,j,k) vup(i,j,k,3)
+#define sx(i,j,k) scon(i,j,k,1)
+#define sy(i,j,k) scon(i,j,k,2)
+#define sz(i,j,k) scon(i,j,k,3)
+
+
+ /*@@
+ @routine GRHydro_PPMReconstruct_drv
+ @date Tue Jul 19 13:24:34 EDT 2011
+ @author Luca Baiotti, Ian Hawke, Bruno C. Mundim, Joshua Faber, Christian D. Ott
+ @desc
+ A driver routine to do PPM reconstructions.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+subroutine GRHydro_PPMReconstruct_drv(CCTK_ARGUMENTS)
+
+ implicit none
+
+ ! save memory when MP is not used
+ ! TARGET as to be before DECLARE_CCTK_ARGUMENTS for gcc 4.1
+ TARGET gaa, gab, gac, gbb, gbc, gcc
+ TARGET gxx, gxy, gxz, gyy, gyz, gzz
+ TARGET betaa, betab, betac
+ TARGET betax, betay, betaz
+ TARGET lvel, vel
+
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+ DECLARE_CCTK_FUNCTIONS
+
+ integer :: nx, ny, nz, i, j, k
+
+ logical, dimension(:,:,:), allocatable :: trivial_rp
+
+ CCTK_INT :: type_bitsx, trivialx, not_trivialx, &
+ &type_bitsy, trivialy, not_trivialy, &
+ &type_bitsz, trivialz, not_trivialz
+
+ CCTK_INT :: ierr
+
+ ! save memory when MP is not used
+ CCTK_INT :: GRHydro_UseGeneralCoordinates
+ CCTK_REAL, DIMENSION(:,:,:), POINTER :: g11, g12, g13, g22, g23, g33
+ CCTK_REAL, DIMENSION(:,:,:), POINTER :: beta1, beta2, beta3
+ CCTK_REAL, DIMENSION(:,:,:,:), POINTER :: vup
+
+ logical :: apply_enhanced_ppm
+
+ if (GRHydro_UseGeneralCoordinates(cctkGH).ne.0) then
+ g11 => gaa
+ g12 => gab
+ g13 => gac
+ g22 => gbb
+ g23 => gbc
+ g33 => gcc
+ beta1 => betaa
+ beta2 => betab
+ beta3 => betac
+ vup => lvel
+ else
+ g11 => gxx
+ g12 => gxy
+ g13 => gxz
+ g22 => gyy
+ g23 => gyz
+ g33 => gzz
+ beta1 => betax
+ beta2 => betay
+ beta3 => betaz
+ vup => vel
+ 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")
+ end if
+
+ call SpaceMask_GetTypeBits(type_bitsx, "Hydro_RiemannProblemX")
+ call SpaceMask_GetStateBits(trivialx, "Hydro_RiemannProblemX", &
+ &"trivial")
+ call SpaceMask_GetStateBits(not_trivialx, "Hydro_RiemannProblemX", &
+ &"not_trivial")
+ call SpaceMask_GetTypeBits(type_bitsy, "Hydro_RiemannProblemY")
+ call SpaceMask_GetStateBits(trivialy, "Hydro_RiemannProblemY", &
+ &"trivial")
+ call SpaceMask_GetStateBits(not_trivialy, "Hydro_RiemannProblemY", &
+ &"not_trivial")
+ call SpaceMask_GetTypeBits(type_bitsz, "Hydro_RiemannProblemZ")
+ call SpaceMask_GetStateBits(trivialz, "Hydro_RiemannProblemZ", &
+ &"trivial")
+ call SpaceMask_GetStateBits(not_trivialz, "Hydro_RiemannProblemZ", &
+ &"not_trivial")
+
+ nx = cctk_lsh(1)
+ ny = cctk_lsh(2)
+ nz = cctk_lsh(3)
+
+ ! if use_enhanced_ppm, allow old PPM on one level
+ if (GRHydro_oppm_reflevel .eq. (-1) .or. &
+ GRHydro_reflevel .ne. GRHydro_oppm_reflevel) then
+ apply_enhanced_ppm = use_enhanced_ppm .ne. 0
+ else
+ apply_enhanced_ppm = .false.
+ end if
+
+
+
+!!$ PPM starts:
+ if (flux_direction == 1) then
+ !$OMP PARALLEL DO PRIVATE(i, j, k)
+ do k = GRHydro_stencil, nz - GRHydro_stencil + 1
+ do j = GRHydro_stencil, ny - GRHydro_stencil + 1
+
+ trivial_rp(:,j,k) = .false.
+ call PPMtyc(nx,CCTK_DELTA_SPACE(1),rho(:,j,k),velx(:,j,k),&
+ vely(:,j,k),velz(:,j,k),eps(:,j,k),temperature(:,j,k),Y_e(:,j,k),&
+ press(:,j,k),rhominus(:,j,k),&
+ velxminus(:,j,k),velyminus(:,j,k),velzminus(:,j,k),&
+ epsminus(:,j,k),tempminus(:,j,k), y_e_minus(:,j,k),&
+ rhoplus(:,j,k),&
+ velxplus(:,j,k),velyplus(:,j,k),velzplus(:,j,k),&
+ epsplus(:,j,k),tempplus(:,j,k),y_e_plus(:,j,k))
+ do i = 1, nx
+ if (trivial_rp(i,j,k)) then
+ SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, trivialx)
+ else
+ SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bitsx, not_trivialx)
+ end if
+ end do
+ end do
+ end do
+ !$OMP END PARALLEL DO
+ else if (flux_direction == 2) then
+ !$OMP PARALLEL DO PRIVATE(i, j, k)
+ do k = GRHydro_stencil, nz - GRHydro_stencil + 1
+ do j = GRHydro_stencil, nx - GRHydro_stencil + 1
+ trivial_rp(j,:,k) = .false.
+ call PPMtyc(ny,CCTK_DELTA_SPACE(2),rho(j,:,k),&
+ vely(j,:,k),velz(j,:,k),velx(j,:,k),&
+ eps(j,:,k),temperature(j,:,k),Y_e(j,:,k), &
+ press(j,:,k),rhominus(j,:,k),&
+ velyminus(j,:,k),velzminus(j,:,k),velxminus(j,:,k),&
+ epsminus(j,:,k),tempminus(j,:,k), y_e_minus(j,:,k),rhoplus(j,:,k),&
+ velyplus(j,:,k),velzplus(j,:,k),velxplus(j,:,k),&
+ epsplus(j,:,k),tempplus(j,:,k),y_e_plus(j,:,k))
+ do i = 1, ny
+ if (trivial_rp(j,i,k)) then
+ SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, trivialy)
+ else
+ SpaceMask_SetStateBitsF90(space_mask, j, i, k, type_bitsy, not_trivialy)
+ end if
+ end do
+ end do
+ end do
+ !$OMP END PARALLEL DO
+
+ else if (flux_direction == 3) then
+ !$OMP PARALLEL DO PRIVATE(i, j, k)
+ do k = GRHydro_stencil, ny - GRHydro_stencil + 1
+ do j = GRHydro_stencil, nx - GRHydro_stencil + 1
+ trivial_rp(j,k,:) = .false.
+ call PPMtyc(nz,CCTK_DELTA_SPACE(3),rho(j,k,:),&
+ velz(j,k,:),velx(j,k,:),vely(j,k,:),&
+ eps(j,k,:),temperature(j,k,:),Y_e(j,k,:), press(j,k,:),rhominus(j,k,:),&
+ velzminus(j,k,:),velxminus(j,k,:),velyminus(j,k,:),&
+ epsminus(j,k,:),tempminus(j,k,:),y_e_minus(j,k,:),rhoplus(j,k,:),&
+ velzplus(j,k,:),velxplus(j,k,:),velyplus(j,k,:),&
+ epsplus(j,k,:),tempplus(j,k,:),y_e_plus(j,k,:))
+ do i = 1, nz
+ if (trivial_rp(j,k,i)) then
+ SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, trivialz)
+ else
+ SpaceMask_SetStateBitsF90(space_mask, j, k, i, type_bitsz, not_trivialz)
+ end if
+ end do
+ end do
+ end do
+ !$OMP END PARALLEL DO
+ else
+ call CCTK_WARN(0, "Flux direction not x,y,z")
+ end if
+!!$ PPM ends.
+
+ deallocate(trivial_rp)
+
+end subroutine GRHydro_PPMReconstruct_drv
+