aboutsummaryrefslogtreecommitdiff
path: root/src/stencil.c
blob: e57f9ff6e75a3ad97382c95c64c0b407ffedb598 (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
66
#include "cctk.h"

#include "util_ErrorCodes.h"
#include "util_Table.h"

#include "stencil.h"



/* Determine whether a boundary with the symmetry handle symbnd is a
   "regular" symmetry boundary (where the SBP stencils should not be
   modified), or an outer boundary (or a multi-block boundary), where
   the SBP stencils need to be modified.  */

void SBP_determine_onesided_stencil (const cGH * cctkGH, int * onesided)
{
  int symtable;
  int pen_sym_handle;
  CCTK_INT symbnd[6];
  int ierr;
  int d;
  
  symtable = SymmetryTableHandleForGrid (cctkGH);
  if (symtable<0) {
    CCTK_WARN(0,"Cannot get symmetry table handle -- maybe thorn SymBase is not active?");
  }
  
  ierr = Util_TableGetIntArray (symtable, 6, symbnd, "symmetry_handle");
  if (ierr!=6) {
    CCTK_WARN(0,"Cannot get symmetry handles");
  }
  
  pen_sym_handle = SymmetryHandleOfName ( "multipatch" );
  
  for (d=0; d<6; ++d) {
    if (! cctkGH->cctk_bbox[d]) {
      /* This is an inter-processor boundary */
      onesided[d] = 0;
    } else {
      /* On an outer boundary (which is not a symmetry boundary),
         symbnd < 0.  */
      if (symbnd[d] < 0) {
        onesided[d] = 1;
      } else {
        if (pen_sym_handle >= 0) {
          /* If the symmetry boundary is a multi-block boundary, then
             symbnd = pen_sym_handle.  However, we can only check that
             thorn MultiPatch is active, i.e., whether pen_sym_handle
             >= 0.  */
          if (symbnd[d] == pen_sym_handle) {
            onesided[d] = 1;
          } else {
            onesided[d] = 0;
          }
        } else {
          onesided[d] = 0;
        }
      }
    }
  }
}

CCTK_FCALL void CCTK_FNAME (SBP_determine_onesided_stencil) (CCTK_POINTER_TO_CONST * cctkGH, int * onesided)
{
  SBP_determine_onesided_stencil (* cctkGH, onesided);
}