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
|
#include <cmath>
#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>
#include <loopcontrol.h>
#include "mask_excluded.hh"
namespace CarpetMask {
using namespace std;
/**
* Set the weight in the excluded regions to zero.
*/
void
CarpetExcludedSetup (CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
CCTK_REAL * const weight =
static_cast <CCTK_REAL *>
(CCTK_VarDataPtr (cctkGH, 0, "CarpetReduce::weight"));
if (not weight) {
CCTK_WARN (CCTK_WARN_ABORT,
"CarpetReduce is not active, or CarpetReduce::mask does not have storage");
}
for (int n = 0; n < 10; ++ n) {
CCTK_REAL const r0 = excluded_radius[n];
if (r0 >= 0.0) {
CCTK_REAL const x0 = excluded_centre_x[n];
CCTK_REAL const y0 = excluded_centre_y[n];
CCTK_REAL const z0 = excluded_centre_z[n];
CCTK_REAL const r2 = pow (r0, 2);
bool const exterior = exclude_exterior[n];
#pragma omp parallel
LC_LOOP3(CarpetExcludedSetup,
i,j,k,
0,0,0, cctk_lsh[0],cctk_lsh[1],cctk_lsh[2],
cctk_lsh[0],cctk_lsh[1],cctk_lsh[2])
{
int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
CCTK_REAL const dx2 = pow (x[ind] - x0, 2);
CCTK_REAL const dy2 = pow (y[ind] - y0, 2);
CCTK_REAL const dz2 = pow (z[ind] - z0, 2);
if (exterior ?
dx2 + dy2 + dz2 >= r2 :
dx2 + dy2 + dz2 <= r2)
{
weight[ind] = 0.0;
}
} LC_ENDLOOP3(CarpetExcludedSetup);
} // if r>=0
} // for n
}
} // namespace CarpetMask
|