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
106
|
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
#include "util_ErrorCodes.h"
#include "util_Table.h"
#include "stencil.h"
void DiffGv2 ( const CCTK_POINTER_TO_CONST cctkGH_, const CCTK_INT dir,
const CCTK_REAL *var, CCTK_REAL *dvar2,
const CCTK_INT table_handle )
{
cGH const * restrict const cctkGH = cctkGH_;
DECLARE_CCTK_PARAMETERS
DECLARE_CCTK_ARGUMENTS
CCTK_INT ni, nj, nk, gsize, loc_order;
CCTK_REAL delta;
CCTK_INT bb[2];
int onesided[6];
int nelements;
void CCTK_FCALL CCTK_FNAME(deriv2_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 *dvar2);
void CCTK_FCALL CCTK_FNAME(deriv2_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 *dvar2);
ni = cctk_lsh[0]; nj = cctk_lsh[1]; nk = cctk_lsh[2];
if ( table_handle >=0 ) {
nelements = Util_TableGetInt ( table_handle, &loc_order, "order" );
if ( nelements == UTIL_ERROR_TABLE_NO_SUCH_KEY ) {
loc_order = order;
} else if ( nelements != 1) {
CCTK_WARN (0, "The options table has an entry \"order\", but it does not have the right properties");
}
} else {
loc_order = order;
}
SBP_determine_onesided_stencil (cctkGH, onesided);
switch(dir) {
case 0: {
delta = CCTK_DELTA_SPACE(0);
bb[0] = onesided[0];
bb[1] = onesided[1];
gsize = cctk_nghostzones[0];
break;
}
case 1: {
delta = CCTK_DELTA_SPACE(1);
bb[0] = onesided[2];
bb[1] = onesided[3];
gsize = cctk_nghostzones[1];
break;
}
case 2: {
delta = CCTK_DELTA_SPACE(2);
bb[0] = onesided[4];
bb[1] = onesided[5];
gsize = cctk_nghostzones[2];
break;
}
default:
CCTK_WARN (0, "Wrong direction specified");
}
if ( CCTK_Equals(norm_type,"Diagonal") ) {
switch(loc_order) {
case 2: {
CCTK_FNAME(deriv2_gf_2_1)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar2);
break;
}
case 4: {
CCTK_FNAME(deriv2_gf_4_2)(var,&ni,&nj,&nk,&dir,bb,&gsize,&delta,dvar2);
break;
}
default:
CCTK_WARN (0, "Unknown 2nd derivative stencil specified");
}
} else {
switch(loc_order) {
default:
CCTK_WARN (0, "Unknown stencil specified");
}
}
}
|