aboutsummaryrefslogtreecommitdiff
path: root/src/surface_integral.c
blob: d84347de7eebe436d687542831474e0a94695084 (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
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
#include <math.h>

#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>

CCTK_INT find_closest(const cGH *cctkGH, CCTK_INT *cctk_lsh,
                      CCTK_REAL *x, CCTK_REAL x_min, CCTK_INT dir)
{
    CCTK_INT i, ijk, min_i = -1;
    CCTK_REAL min = 1.e100;
    for(i=3; i<cctk_lsh[dir]-3; i++)
    {
        ijk = CCTK_GFINDEX3D(cctkGH, (dir==0)?i:0, (dir==1)?i:0, (dir==2)?i:0);
        if (fabs(x[ijk] - x_min) < min)
        {
            min = fabs(x[ijk] - x_min);
            min_i = i;
        }
    }
    return min_i;
}

void ADMMass_Surface(CCTK_ARGUMENTS)
{
    DECLARE_CCTK_ARGUMENTS
    DECLARE_CCTK_PARAMETERS

    CCTK_INT i,j,k, ijk;
    CCTK_REAL x_min, x_max, y_min, y_max, z_min, z_max;
    CCTK_INT  x_min_i, x_max_i, y_min_j, y_max_j, z_min_k, z_max_k, min;
    CCTK_REAL detg, idetg, ds[3];
    CCTK_REAL u[3][3], dg[3][3][3];
    CCTK_REAL local_sum = 0.0;

    /* grid-function strides for ADMMacros */
    const CCTK_INT di = 1;
    const CCTK_INT dj = cctk_lsh[0];
    const CCTK_INT dk = cctk_lsh[0]*cctk_lsh[1];

    if (ADMMass_boundary_distance[*ADMMass_LoopCounter] > 0.0)
    {
        x_min = ADMMass_x_pos[*ADMMass_LoopCounter] -
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
        x_max = ADMMass_x_pos[*ADMMass_LoopCounter] +
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
        y_min = ADMMass_y_pos[*ADMMass_LoopCounter] -
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
        y_max = ADMMass_y_pos[*ADMMass_LoopCounter] +
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
        z_min = ADMMass_z_pos[*ADMMass_LoopCounter] -
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
        z_max = ADMMass_z_pos[*ADMMass_LoopCounter] +
                ADMMass_boundary_distance[*ADMMass_LoopCounter];
    }
    else
    {
        x_min = ADMMass_x_min[*ADMMass_LoopCounter];
        x_max = ADMMass_x_max[*ADMMass_LoopCounter];
        y_min = ADMMass_y_min[*ADMMass_LoopCounter];
        y_max = ADMMass_y_max[*ADMMass_LoopCounter];
        z_min = ADMMass_z_min[*ADMMass_LoopCounter];
        z_max = ADMMass_z_max[*ADMMass_LoopCounter];
    }
    x_min_i = find_closest(cctkGH, cctk_lsh, x, x_min, 0);
    x_max_i = find_closest(cctkGH, cctk_lsh, x, x_max, 0);
    y_min_j = find_closest(cctkGH, cctk_lsh, y, y_min, 1);
    y_max_j = find_closest(cctkGH, cctk_lsh, y, y_max, 1);
    z_min_k = find_closest(cctkGH, cctk_lsh, z, z_min, 2);
    z_max_k = find_closest(cctkGH, cctk_lsh, z, z_max, 2);

#include "CactusEinstein/ADMMacros/src/macro/UPPERMET_declare.h"
#include "CactusEinstein/ADMMacros/src/macro/DG_declare.h"
#include "CactusEinstein/ADMMacros/src/macro/ADM_Spacing_declare.h"

    for(i=3; i<cctk_lsh[0]-3; i++)
     for(j=3; j<cctk_lsh[1]-3; j++)
      for(k=3; k<cctk_lsh[2]-3; k++)
      {
          ijk = CCTK_GFINDEX3D(cctkGH, i, j, k);
          ADMMass_SurfaceMass_GF[ijk] = 0.0;
          if ((i >= x_min_i) &&
              (i <= x_max_i) &&
              (j >= y_min_j) &&
              (j <= y_max_j) &&
              (k >= z_min_k) &&
              (k <= z_max_k))
          {
              ds[0] = 0.0;
              ds[1] = 0.0;
              ds[2] = 0.0;
#include "CactusEinstein/ADMMacros/src/macro/ADM_Spacing.h"
              if ((i == x_min_i) && (j != y_max_j) && (k != z_max_k))
                  ds[0] = -dy*dz;
              if ((i == x_max_i) && (j != y_max_j) && (k != z_max_k))
                  ds[0] = dy*dz;
              if ((j == y_min_j) && (i != x_max_i) && (k != z_max_k))
                  ds[1] = -dx*dz;
              if ((j == y_max_j) && (i != x_max_i) && (k != z_max_k))
                  ds[1] = dx*dz;
              if ((k == z_min_k) && (j != y_max_j) && (i != x_max_i))
                  ds[2] = -dx*dy;
              if ((k == z_max_k) && (j != y_max_j) && (i != x_max_i))
                  ds[2] = dx*dy;
              if ((ds[0] != 0.0) || (ds[1] != 0.0) || (ds[2] != 0.0))
              {
#include "CactusEinstein/ADMMacros/src/macro/UPPERMET_guts.h"
#include "CactusEinstein/ADMMacros/src/macro/DG_guts.h"
                  u[0][0] = UPPERMET_UXX;
                  u[0][1] = UPPERMET_UXY;
                  u[0][2] = UPPERMET_UXZ;
                  u[1][0] = UPPERMET_UXY;
                  u[1][1] = UPPERMET_UYY;
                  u[1][2] = UPPERMET_UYZ;
                  u[2][0] = UPPERMET_UXZ;
                  u[2][1] = UPPERMET_UYZ;
                  u[2][2] = UPPERMET_UZZ;
                  dg[0][0][0] = DXDG_DXDGXX;
                  dg[0][0][1] = DYDG_DYDGXX;
                  dg[0][0][2] = DZDG_DZDGXX;
                  dg[0][1][0] = DXDG_DXDGXY;
                  dg[0][1][1] = DYDG_DYDGXY;
                  dg[0][1][2] = DZDG_DZDGXY;
                  dg[0][2][0] = DXDG_DXDGXZ;
                  dg[0][2][1] = DYDG_DYDGXZ;
                  dg[0][2][2] = DZDG_DZDGXZ;
                  dg[1][0][0] = DXDG_DXDGXY;
                  dg[1][0][1] = DYDG_DYDGXY;
                  dg[1][0][2] = DZDG_DZDGXY;
                  dg[1][1][0] = DXDG_DXDGYY;
                  dg[1][1][1] = DYDG_DYDGYY;
                  dg[1][1][2] = DZDG_DZDGYY;
                  dg[1][2][0] = DXDG_DXDGYZ;
                  dg[1][2][1] = DYDG_DYDGYZ;
                  dg[1][2][2] = DZDG_DZDGYZ;
                  dg[2][0][0] = DXDG_DXDGXZ;
                  dg[2][0][1] = DYDG_DYDGXZ;
                  dg[2][0][2] = DZDG_DZDGXZ;
                  dg[2][1][0] = DXDG_DXDGYZ;
                  dg[2][1][1] = DYDG_DYDGYZ;
                  dg[2][1][2] = DZDG_DZDGYZ;
                  dg[2][2][0] = DXDG_DXDGZZ;
                  dg[2][2][1] = DYDG_DYDGZZ;
                  dg[2][2][2] = DZDG_DZDGZZ;
                  for (int ti = 0; ti < 3; ti++)
                   for (int tj = 0; tj < 3; tj++)
                    for (int tk = 0; tk < 3; tk++)
                     for (int tl = 0; tl < 3; tl++)
                         ADMMass_SurfaceMass_GF[ijk] +=
                             u[ti][tj] * u[tk][tl] *
                             ( dg[ti][tk][tj] - dg[ti][tj][tk] ) * ds[tk];
                  ADMMass_SurfaceMass_GF[ijk] *= sqrt(DETG_DETG) / 16 / 3.14159265;
                  local_sum += ADMMass_SurfaceMass_GF[ijk];
              }
          }
      }
    #include "CactusEinstein/ADMMacros/src/macro/UPPERMET_undefine.h"
    #include "CactusEinstein/ADMMacros/src/macro/DG_undefine.h"
    #include "CactusEinstein/ADMMacros/src/macro/ADM_Spacing_undefine.h"
    /*CCTK_VInfo(CCTK_THORNSTRING,
     *           "ADM mass local to one processor (surface): %g\n", local_sum);
     */
}

void ADMMass_Surface_Global(CCTK_ARGUMENTS)
{
    DECLARE_CCTK_ARGUMENTS
    DECLARE_CCTK_PARAMETERS

    CCTK_INT reduction_handle;

    reduction_handle = CCTK_ReductionHandle("sum");
    if (reduction_handle < 0)
        CCTK_WARN(0, "Unable to ge reduction handle.");

    if (CCTK_Reduce(cctkGH, -1, reduction_handle, 1,
                    CCTK_VARIABLE_REAL,
                    &ADMMass_SurfaceMass[*ADMMass_LoopCounter], 1,
                    CCTK_VarIndex("ADMMass::ADMMass_SurfaceMass_GF")))
        CCTK_WARN(0, "Error while reducing ADMMass_SurfaceMass_GF");
    
    CCTK_VInfo(CCTK_THORNSTRING,
               "ADM mass, surface, detector %d: %g %g\n", *ADMMass_LoopCounter,
               ADMMass_boundary_distance[*ADMMass_LoopCounter],
               ADMMass_SurfaceMass[*ADMMass_LoopCounter]);
}

void ADMMass_Surface_Lapse(CCTK_ARGUMENTS)
{
    DECLARE_CCTK_ARGUMENTS
    DECLARE_CCTK_PARAMETERS

    CCTK_INT i,j,k, ijk;
    CCTK_REAL local_sum = 0.0;

    for(i=3; i<cctk_lsh[0]-3; i++)
     for(j=3; j<cctk_lsh[1]-3; j++)
      for(k=3; k<cctk_lsh[2]-3; k++)
      {
          ijk = CCTK_GFINDEX3D(cctkGH, i, j, k);
          ADMMass_SurfaceMass_GF[ijk] *= alp[ijk];
          local_sum += ADMMass_SurfaceMass_GF[ijk];
      }
}

void ADMMass_Surface_Lapse_Global(CCTK_ARGUMENTS)
{
    DECLARE_CCTK_ARGUMENTS
    DECLARE_CCTK_PARAMETERS

    CCTK_INT reduction_handle;

    reduction_handle = CCTK_ReductionHandle("sum");
    if (reduction_handle < 0)
        CCTK_WARN(0, "Unable to ge reduction handle.");

    if (CCTK_Reduce(cctkGH, -1, reduction_handle, 1,
                    CCTK_VARIABLE_REAL,
                    &ADMMass_SurfaceMass_Lapse[*ADMMass_LoopCounter], 1,
                    CCTK_VarIndex("ADMMass::ADMMass_SurfaceMass_GF")))
        CCTK_WARN(0, "Error while reducing ADMMass_SurfaceMass_GF");
    
    CCTK_VInfo(CCTK_THORNSTRING,
           "ADM mass, surface lapse, detector %d: %g %g\n",
           *ADMMass_LoopCounter,
           ADMMass_boundary_distance[*ADMMass_LoopCounter],
           ADMMass_SurfaceMass_Lapse[*ADMMass_LoopCounter]);
}