aboutsummaryrefslogtreecommitdiff
path: root/src/call_derivs.c
blob: e8c051ab5d7c1e47106bc02587d4f5e6a6b6f18b (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
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

#include <assert.h>

void DiffGv ( const CCTK_POINTER cctkGH_, const CCTK_INT dir,
              const CCTK_REAL *var, CCTK_REAL *dvar )
{
  cGH const * restrict const cctkGH = cctkGH_;
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_ARGUMENTS

  CCTK_INT ni, nj, nk, gsize, ic;
  CCTK_REAL delta;
  CCTK_INT ierr;
  CCTK_INT lsh[3], bbox[6], bb[2], nghostzones[3];
  void CCTK_FCALL CCTK_FNAME(deriv_gf_2_1)(const 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)(const 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)(const 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_8_4)(const 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);

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

  switch(dir) {
  case 0: {
    delta = CCTK_DELTA_SPACE(0);
    bb[0] = cctk_bbox[0]; bb[1] = cctk_bbox[1];
    gsize = cctk_nghostzones[0];
    break;
  }
  case 1: {
    delta = CCTK_DELTA_SPACE(1);
    bb[0] = cctk_bbox[2]; bb[1] = cctk_bbox[3];
    gsize = cctk_nghostzones[1];
    break;
  }
  case 2: {
    delta = CCTK_DELTA_SPACE(2);
    bb[0] = cctk_bbox[4]; bb[1] = cctk_bbox[5];
    gsize = cctk_nghostzones[2];
    break;
  }
  default:
      assert(0);
  }

  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;
  }
  case 8: {
    CCTK_FNAME(deriv_gf_8_4)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar);
    break;
  }
  default:
    assert(0);
  }
}