aboutsummaryrefslogtreecommitdiff
path: root/src/call_derivs_name.c
blob: c696e4ebd40000660d73ca231ed9d2a4dc0f1b32 (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
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
#include "cctk.h"
#include "cctk_Parameters.h"

#include <assert.h>

void DiffGf ( const CCTK_POINTER cctkGH, const CCTK_INT dir,
              const char *var_name, const char *dvar_name )
{
  CCTK_REAL *var, *dvar;
  CCTK_INT ni, nj, nk, gsize, ic;
  CCTK_REAL delta;
  CCTK_REAL phys_min[3], phys_max[3];
  CCTK_REAL int_min[3], int_max[3];
  CCTK_REAL ext_min[3], ext_max[3];
  CCTK_REAL spacing[3];
  CCTK_INT ierr;
  CCTK_INT lsh[3], bbox[6], bb[2], nghostzones[3];
  void CCTK_FCALL CCTK_FNAME(deriv_gf_2_1)(CCTK_REAL *var, const CCTK_INT *ni,
                                       const CCTK_INT *nj, const CCTK_INT *nk,
                                       const CCTK_INT *dir,
                                       const CCTK_INT *bb,
                                       const CCTK_INT *gsize,
                                       const CCTK_REAL *delta,
                                       CCTK_REAL *dvar);
  void CCTK_FCALL CCTK_FNAME(deriv_gf_4_2)(CCTK_REAL *var, const CCTK_INT *ni,
                                       const CCTK_INT *nj, const CCTK_INT *nk,
                                       const CCTK_INT *dir,
                                       const CCTK_INT *bb,
                                       const CCTK_INT *gsize,
                                       const CCTK_REAL *delta,
                                       CCTK_REAL *dvar);
  void CCTK_FCALL CCTK_FNAME(deriv_gf_6_3)(CCTK_REAL *var, const CCTK_INT *ni,
                                       const CCTK_INT *nj, const CCTK_INT *nk,
                                       const CCTK_INT *dir,
                                       const CCTK_INT *bb,
                                       const CCTK_INT *gsize,
                                       const CCTK_REAL *delta,
                                       CCTK_REAL *dvar);
  DECLARE_CCTK_PARAMETERS

  
  if (CCTK_GrouplshVN(cctkGH, 3, lsh, var_name) < 0)
  {
    CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Error getting local size of grid for variable '%s'",var_name);
  }
  if (CCTK_GroupbboxVN(cctkGH, 6, bbox, var_name) < 0)
  {
    CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Error getting bounding box for variable '%s'",var_name);
  }
  if (CCTK_GroupnghostzonesVN(cctkGH, 3, nghostzones, var_name) < 0)
  {
    CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Error getting number of ghostzones for variable '%s'",var_name);
  }

  ni = lsh[0]; nj = lsh[1]; nk = lsh[2];

  ierr = GetDomainSpecification ( 3, phys_min, phys_max, int_min, int_max,
                                                  ext_min, ext_max, spacing );
  
  switch(dir) {
  case 0: {
    delta = spacing[0];
    bb[0] = bbox[0]; bb[1] = bbox[1];
    gsize = nghostzones[0];
    break;
  }
  case 1: {
    delta = spacing[1];
    bb[0] = bbox[2]; bb[1] = bbox[3];
    gsize = nghostzones[1];
    break;
  }
  case 2: {
    delta = spacing[2];
    bb[0] = bbox[4]; bb[1] = bbox[5];
    gsize = nghostzones[2];
    break;
  }
  default:
      assert(0);
  }

  var = (CCTK_REAL *)(CCTK_VarDataPtr(cctkGH,0,var_name));
  dvar = (CCTK_REAL *)(CCTK_VarDataPtr(cctkGH,0,dvar_name));

  switch(order) {
  case 2: {
    CCTK_FNAME(deriv_gf_2_1)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar);
    break;
  }
  case 4: {
    CCTK_FNAME(deriv_gf_4_2)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar);
    break;
  }
  case 6: {
    CCTK_FNAME(deriv_gf_6_3)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar);
    break;
  }
  default:
    assert(0);
  }
}