From 2bd3c3a7d785398654dc3fae9790a0ca33957504 Mon Sep 17 00:00:00 2001 From: diener Date: Wed, 1 Dec 2004 13:44:39 +0000 Subject: Routine to check if grid sizes and ghost_zones are enough for the order of finite differencing requires. Has not been tested thoroughly. git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@16 f69c4107-0314-4c4f-9ad4-17e986b73f4a --- src/CheckGridSizes.F90 | 54 ++++++++++++++++++++++++++++++++++++++++++++++++++ src/make.code.defn | 2 +- 2 files changed, 55 insertions(+), 1 deletion(-) create mode 100644 src/CheckGridSizes.F90 (limited to 'src') diff --git a/src/CheckGridSizes.F90 b/src/CheckGridSizes.F90 new file mode 100644 index 0000000..3240f6a --- /dev/null +++ b/src/CheckGridSizes.F90 @@ -0,0 +1,54 @@ +! Check if there are enough grid cells and ghost zones to compute the finite +! differences. +! $Header$ + +#include "cctk.h" +#include "cctk_Parameters.h" +#include "cctk_Arguments.h" +#include "cctk_Functions.h" + +subroutine SBP_CheckGridSizes(CCTK_ARGUMENTS) + + implicit none + + DECLARE_CCTK_PARAMETERS + DECLARE_CCTK_ARGUMENTS + DECLARE_CCTK_FUNCTIONS + + integer :: i + integer, dimension(3) :: bc + integer, dimension(3,3) :: np + character(len=128) :: info_message + character(len=1), dimension(3) :: dir = (/ 'x', 'y', 'z' /) + + bc(1) = cctk_bbox(1)+cctk_bbox(2)+1 ! 0 processor boundaries -> bc = 3 + bc(2) = cctk_bbox(3)+cctk_bbox(4)+1 ! 1 processor boundary -> bc = 2 + bc(3) = cctk_bbox(5)+cctk_bbox(6)+1 ! 2 processor boundaries -> bc = 1 + + if ( any (cctk_nghostzones < order / 2) .and. CCTK_nProcs(cctkGH) > 1 ) then + write(info_message,'(a22,i1,a13,i1,a25)') 'You need ghostsize >= ', & + order/2, ' to run with ', order, & + ' order finite differences' + call CCTK_WARN(0,info_message) + end if + + if ( CCTK_EQUALS(norm_type,"Diagonal") ) then + np(1,:) = 2*cctk_nghostzones + np(2,:) = 3*order/2 + np(3,:) = 2*order + else + np(1,:) = 2*cctk_nghostzones + np(2,:) = 3*order/2 + 1 + np(3,:) = 2*(order+1) + end if + + do i = 1, 3 + if ( cctk_lsh(i) < np(bc(i),i) ) then + write(info_message,'(a27,a1,a35,i1,a25)') & + 'Your local gridsize in the ', dir(i), & + '-direction is to small to run with ', order, & + ' order finite differences' + call CCTK_WARN(0,info_message) + end if + end do +end subroutine SBP_CheckGridSizes diff --git a/src/make.code.defn b/src/make.code.defn index 893ee3f..5810dc0 100644 --- a/src/make.code.defn +++ b/src/make.code.defn @@ -2,7 +2,7 @@ # $Header$ # Source files in this directory -SRCS = Derivatives_2_1.F90 Derivatives_4_2.F90 Derivatives_6_3.F90 Derivatives_8_4.F90 call_derivs.c Derivatives_4_3.F90 call_derivs_name.c Get_Coeff.F90 call_norm_mask.c set_norm_mask.F90 +SRCS = Derivatives_2_1.F90 Derivatives_4_2.F90 Derivatives_6_3.F90 Derivatives_8_4.F90 call_derivs.c Derivatives_4_3.F90 call_derivs_name.c Get_Coeff.F90 call_norm_mask.c set_norm_mask.F90 CheckGridSizes.F90 # Subdirectories containing source files SUBDIRS = -- cgit v1.2.3