aboutsummaryrefslogtreecommitdiff
path: root/src/CheckGridSizes.F90
diff options
context:
space:
mode:
authordiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2004-12-01 13:44:39 +0000
committerdiener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a>2004-12-01 13:44:39 +0000
commit2bd3c3a7d785398654dc3fae9790a0ca33957504 (patch)
treedc8c82a8a822ce295e0ca6712da4ea67c59c7631 /src/CheckGridSizes.F90
parente424c9749e261bbadab712bcaf1cad190143291a (diff)
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
Diffstat (limited to 'src/CheckGridSizes.F90')
-rw-r--r--src/CheckGridSizes.F9054
1 files changed, 54 insertions, 0 deletions
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