aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetReduce/src/mask_test.c
blob: cbf57b3b36a4782940602aec4d9bfbe568d17f92 (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
#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>

#include <assert.h>
#include <math.h>



void
MaskBase_TestMask (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  if (verbose) {
    CCTK_INFO ("Testing weight");
  }
  
  
  
  int const sum = CCTK_ReductionHandle ("sum");
  assert (sum >= 0);
  
  int const proc = 0;
  
  int const weight_var = CCTK_VarIndex ("CarpetReduce::weight");
  int const one_var    = CCTK_VarIndex ("CarpetReduce::one");
  assert (weight_var >= 0);
  
  CCTK_REAL sum_weight;
  
  {
    int const ierr = CCTK_Reduce (cctkGH,
                                  proc,
                                  sum,
                                  1, CCTK_VARIABLE_REAL, &sum_weight,
                                  1, one_var);
    assert (ierr >= 0);
  }
  
  
  
  if (CCTK_MyProc(cctkGH) == proc) {
    
    CCTK_REAL physical_min[cctk_dim];
    CCTK_REAL physical_max[cctk_dim];
    CCTK_REAL interior_min[cctk_dim];
    CCTK_REAL interior_max[cctk_dim];
    CCTK_REAL exterior_min[cctk_dim];
    CCTK_REAL exterior_max[cctk_dim];
    CCTK_REAL spacing     [cctk_dim];
    int const ierr = GetDomainSpecification (cctk_dim,
                                             physical_min,
                                             physical_max,
                                             interior_min,
                                             interior_max,
                                             exterior_min,
                                             exterior_max,
                                             spacing);
    assert (!ierr);
    
    CCTK_REAL domain_volume = 1.0;
    for (int d=0; d<cctk_dim; ++d) {
      domain_volume *= (physical_max[d] - physical_min[d]) / spacing[d];
    }
    
    
    
    int const there_is_a_problem =
      fabs(sum_weight - domain_volume) > 1.0e-12 * (sum_weight + domain_volume);
    
    
    
    if (verbose || there_is_a_problem) {
      CCTK_VInfo (CCTK_THORNSTRING,
                  "Simulation domain volume: %.15g", (double)domain_volume);
      CCTK_VInfo (CCTK_THORNSTRING,
                  "Reduction weight sum:     %.15g", (double)sum_weight);
    }
    
    if (there_is_a_problem) {
      CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Simulation domain volume and reduction weight sum differ");
    }
    
  }
  
  
  
  int const iret1 = CCTK_DisableGroupStorage (cctkGH, "CarpetReduce::iweight");
  int const iret2 = CCTK_DisableGroupStorage (cctkGH, "CarpetReduce::one");
  assert (iret1==1);
  assert (iret2==1);
}