aboutsummaryrefslogtreecommitdiff
path: root/src/surface_integral.c
diff options
context:
space:
mode:
Diffstat (limited to 'src/surface_integral.c')
-rw-r--r--src/surface_integral.c191
1 files changed, 134 insertions, 57 deletions
diff --git a/src/surface_integral.c b/src/surface_integral.c
index d84347d..7927a06 100644
--- a/src/surface_integral.c
+++ b/src/surface_integral.c
@@ -1,20 +1,27 @@
#include <math.h>
+#include <stdio.h>
+#include <stdlib.h>
+#include <assert.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)
+#define PI 3.141592653589793238462643383279502884197169399375105820974944592307816406286208998628034825342117068
+
+CCTK_INT find_closest(const cGH *cctkGH, CCTK_INT *cctk_lsh, CCTK_REAL *cctk_delta_space,
+ CCTK_INT ghost, CCTK_REAL *coord, CCTK_REAL coord_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++)
+
+ for(i=ghost; i<cctk_lsh[dir]-ghost; i++)
{
ijk = CCTK_GFINDEX3D(cctkGH, (dir==0)?i:0, (dir==1)?i:0, (dir==2)?i:0);
- if (fabs(x[ijk] - x_min) < min)
+
+ if (fabs(coord[ijk] - coord_min) < min)
{
- min = fabs(x[ijk] - x_min);
+ min = fabs(coord[ijk] - coord_min);
min_i = i;
}
}
@@ -26,60 +33,122 @@ 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_INT i,j,k, ijk, ierr, ghost;
+ CCTK_INT x_min_i, x_max_i, y_min_j, y_max_j, z_min_k, z_max_k;
CCTK_REAL detg, idetg, ds[3];
CCTK_REAL u[3][3], dg[3][3][3];
- CCTK_REAL local_sum = 0.0;
+
+ CCTK_REAL physical_min[3];
+ CCTK_REAL physical_max[3];
+ CCTK_REAL interior_min[3];
+ CCTK_REAL interior_max[3];
+ CCTK_REAL exterior_min[3];
+ CCTK_REAL exterior_max[3];
+ CCTK_REAL spacing[3];
+
+ char coordinate_parameter_type[9];
/* 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];
+
+ x_min_i = -10000000;
+ x_max_i = 10000000;
+ y_min_j = -10000000;
+ y_max_j = 10000000;
+ z_min_k = -10000000;
+ z_max_k = 10000000;
- 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];
- }
+ if (ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter] > 0.0)
+ {
+ if (strcmp(CCTK_ParameterValString("type", "cartgrid3d"), "coordbase"))
+ {
+ CCTK_WARN(0,"this thorns used with the ADMMass_distance_from_grid_boundar parameter requires to set coordinates through coordbase.");
+ }
+
+ /* Find the physical coordinates of the boundaries */
+ ierr = GetDomainSpecification
+ (3, physical_min, physical_max, interior_min, interior_max,
+ exterior_min, exterior_max, spacing);
+ if (ierr)
+ CCTK_VWarn ( 0, __LINE__, __FILE__, "CartGrid3D", "error returned from function GetDomainSpecification");
+ *ADMMass_box_x_min = physical_min[0] + ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ *ADMMass_box_x_max = physical_max[0] - ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ *ADMMass_box_y_min = physical_min[1] + ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ *ADMMass_box_y_max = physical_max[1] - ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ *ADMMass_box_z_min = physical_min[2] + ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ *ADMMass_box_z_max = physical_max[2] - ADMMass_distance_from_grid_boundary[*ADMMass_LoopCounter];
+ }
+ else if (ADMMass_surface_distance[*ADMMass_LoopCounter] > 0.0)
+ {
+ *ADMMass_box_x_min = ADMMass_x_pos[*ADMMass_LoopCounter] -
+ ADMMass_surface_distance[*ADMMass_LoopCounter];
+ *ADMMass_box_x_max = ADMMass_x_pos[*ADMMass_LoopCounter] +
+ ADMMass_surface_distance[*ADMMass_LoopCounter];
+ *ADMMass_box_y_min = ADMMass_y_pos[*ADMMass_LoopCounter] -
+ ADMMass_surface_distance[*ADMMass_LoopCounter];
+ *ADMMass_box_y_max = ADMMass_y_pos[*ADMMass_LoopCounter] +
+ ADMMass_surface_distance[*ADMMass_LoopCounter];
+ *ADMMass_box_z_min = ADMMass_z_pos[*ADMMass_LoopCounter] -
+ ADMMass_surface_distance[*ADMMass_LoopCounter];
+ *ADMMass_box_z_max = ADMMass_z_pos[*ADMMass_LoopCounter] +
+ ADMMass_surface_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);
+ {
+ *ADMMass_box_x_min = ADMMass_x_min[*ADMMass_LoopCounter];
+ *ADMMass_box_x_max = ADMMass_x_max[*ADMMass_LoopCounter];
+ *ADMMass_box_y_min = ADMMass_y_min[*ADMMass_LoopCounter];
+ *ADMMass_box_y_max = ADMMass_y_max[*ADMMass_LoopCounter];
+ *ADMMass_box_z_min = ADMMass_z_min[*ADMMass_LoopCounter];
+ *ADMMass_box_z_max = ADMMass_z_max[*ADMMass_LoopCounter];
+ }
+
+ if (CCTK_IsThornActive("PUGH"))
+ {
+ const void *ghost_ptr = CCTK_ParameterGet("ghost_size","PUGH",NULL);
+ assert( ghost_ptr != NULL );
+ ghost = *(const int *)ghost_ptr;
+ }
+ else /* carpet */
+ {
+ const void *ghost_ptr = CCTK_ParameterGet("ghost_size","Carpet",NULL);
+ assert( ghost_ptr != NULL );
+ ghost = *(const int *)ghost_ptr;
+ }
+
+ ijk = CCTK_GFINDEX3D(cctkGH, cctk_ubnd[0]-cctk_lbnd[0], cctk_ubnd[1]-cctk_lbnd[1],cctk_ubnd[2]-cctk_lbnd[2] );
+
+ if ( (x[0] < *ADMMass_box_x_min) && (*ADMMass_box_x_min < x[ijk]) )
+ x_min_i = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, x, *ADMMass_box_x_min, 0);
+
+ if ( (x[0] < *ADMMass_box_x_max) && (*ADMMass_box_x_max < x[ijk]) )
+ x_max_i = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, x, *ADMMass_box_x_max, 0);
+
+ if ( (y[0] < *ADMMass_box_y_min) && (*ADMMass_box_y_min < y[ijk]) )
+ y_min_j = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, y, *ADMMass_box_y_min, 1);
+
+ if ( (y[0] < *ADMMass_box_y_max) && (*ADMMass_box_y_max < y[ijk]) )
+ y_max_j = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, y, *ADMMass_box_y_max, 1);
+
+ if ( (z[0] < *ADMMass_box_z_min) && (*ADMMass_box_z_min < z[ijk]) )
+ z_min_k = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, z, *ADMMass_box_z_min, 2);
+
+ if ( (z[0] < *ADMMass_box_z_max) && (*ADMMass_box_z_max < z[ijk]) )
+ z_max_k = find_closest(cctkGH, cctk_lsh, cctk_delta_space, ghost, z, *ADMMass_box_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++)
+ for(i=ghost; i<cctk_lsh[0]-ghost; i++)
+ for(j=ghost; j<cctk_lsh[1]-ghost; j++)
+ for(k=ghost; k<cctk_lsh[2]-ghost; k++)
{
ijk = CCTK_GFINDEX3D(cctkGH, i, j, k);
ADMMass_SurfaceMass_GF[ijk] = 0.0;
- if ((i >= x_min_i) &&
+ if ((i >= x_min_i) &&
(i <= x_max_i) &&
(j >= y_min_j) &&
(j <= y_max_j) &&
@@ -149,11 +218,12 @@ void ADMMass_Surface(CCTK_ARGUMENTS)
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];
+ ADMMass_SurfaceMass_GF[ijk] *= sqrt(DETG_DETG);
+ /* 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"
@@ -171,18 +241,23 @@ void ADMMass_Surface_Global(CCTK_ARGUMENTS)
reduction_handle = CCTK_ReductionHandle("sum");
if (reduction_handle < 0)
- CCTK_WARN(0, "Unable to ge reduction handle.");
+ CCTK_WARN(0, "Unable to get 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");
-
+
+ ADMMass_SurfaceMass[*ADMMass_LoopCounter] /= 16.0*PI;
+
CCTK_VInfo(CCTK_THORNSTRING,
- "ADM mass, surface, detector %d: %g %g\n", *ADMMass_LoopCounter,
- ADMMass_boundary_distance[*ADMMass_LoopCounter],
- ADMMass_SurfaceMass[*ADMMass_LoopCounter]);
+ "detector %d: ADM mass %g, surface (cube): xmin %g, xmax %g, ymin %g, ymax %g, zmin %g, zmax %g\n",
+ *ADMMass_LoopCounter,
+ ADMMass_SurfaceMass[*ADMMass_LoopCounter],
+ *ADMMass_box_x_min,*ADMMass_box_x_max,
+ *ADMMass_box_y_min,*ADMMass_box_y_max,
+ *ADMMass_box_z_min,*ADMMass_box_z_max);
}
void ADMMass_Surface_Lapse(CCTK_ARGUMENTS)
@@ -191,15 +266,13 @@ void ADMMass_Surface_Lapse(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];
}
}
@@ -219,11 +292,15 @@ void ADMMass_Surface_Lapse_Global(CCTK_ARGUMENTS)
&ADMMass_SurfaceMass_Lapse[*ADMMass_LoopCounter], 1,
CCTK_VarIndex("ADMMass::ADMMass_SurfaceMass_GF")))
CCTK_WARN(0, "Error while reducing ADMMass_SurfaceMass_GF");
-
+
+ ADMMass_SurfaceMass_Lapse[*ADMMass_LoopCounter] /= 16.0*PI;
+
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]);
+ "detector %d: ADM mass with lapse %g, surface (cube): xmin %g, xmax %g, ymin %g, ymax %g, zmin %g, zmax %g\n",
+ *ADMMass_LoopCounter,
+ ADMMass_SurfaceMass_Lapse[*ADMMass_LoopCounter],
+ *ADMMass_box_x_min,*ADMMass_box_x_max,
+ *ADMMass_box_y_min,*ADMMass_box_y_max,
+ *ADMMass_box_z_min,*ADMMass_box_z_max);
}