aboutsummaryrefslogtreecommitdiff
path: root/src/CheckGridSizes.F90
blob: b4366d9f27f802c6291a728084bc16d6a26e9c3f (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
! 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

  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 ( ( use_dissipation > 0 ) .and. &
       ( CCTK_EQUALS(dissipation_type,"Kreiss-Oliger" ) ) .and. &
       ( any (cctk_nghostzones < order / 2 + 1) ) .and. &
       ( CCTK_nProcs(cctkGH) > 1 ) ) then
    write(info_message,'(a22,i1,a13,i2,a32)') 'You need ghostsize >= ', &
                                    order/2+1, ' to run with ', order+2, &
                                    ' order Kreiss-Oliger dissipation'
    call CCTK_WARN(0,info_message) 
  end if

  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)') &

         '-direction is to small to run with ', order, &
         ' order finite differences'
      ! This is a warning only since this thorn may be active, but
      ! never called to calculate finite differences
      call CCTK_WARN(1,info_message)
    end if
  end do
end subroutine SBP_CheckGridSizes