aboutsummaryrefslogtreecommitdiff
path: root/Carpet/Carpet/src/Poison.cc
blob: 7f1a00f2ecc2be81b69e810c773c53fff369e307 (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
#include <assert.h>
#include <stdlib.h>
#include <string.h>

#include "cctk.h"
#include "cctk_Parameters.h"

#include "carpet.hh"

extern "C" {
  static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Poison.cc,v 1.17 2004/03/23 19:30:14 schnetter Exp $";
  CCTK_FILEVERSION(Carpet_Carpet_Poison_cc);
}



namespace Carpet {
  
  using namespace std;
  
  
  
  void Poison (const cGH* cgh, const checktimes where)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (! poison_new_timelevels) return;
    
    for (int group=0; group<CCTK_NumGroups(); ++group) {
      if (CCTK_QueryGroupStorageI(cgh, group)) {
	PoisonGroup (cgh, group, where);
      } // if has storage
    } // for group
  }
  
  
  
  void PoisonGroup (const cGH* cgh, const int group, const checktimes where)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (group>=0 && group<CCTK_NumGroups());
    
    if (! poison_new_timelevels) return;
    
    if (! CCTK_QueryGroupStorageI(cgh, group)) {
      char * const groupname = CCTK_GroupName(group);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
		  "Cannot poison group \"%s\" because it has no storage",
		  groupname);
      free (groupname);
      return;
    }
    
    const int n0 = CCTK_FirstVarIndexI(group);
    assert (n0>=0);
    const int nvar = CCTK_NumVarsInGroupI(group);
    const int sz = CCTK_VarTypeSize(CCTK_VarTypeI(n0));
    assert (sz>0);
    
    const int num_tl = CCTK_NumTimeLevelsFromVarI(n0);
    assert (num_tl>0);
    const int min_tl = mintl(where, num_tl);
    const int max_tl = maxtl(where, num_tl);
    
    if (min_tl <= max_tl) {
      
      {
        char * const groupname = CCTK_GroupName(group);
        Checkpoint ("PoisonGroup \"%s\"", groupname);
        free (groupname);
      }
      
      const int grouptype = CCTK_GroupTypeI(group);
      
      BEGIN_MAP_LOOP(cgh, grouptype) {
        BEGIN_LOCAL_COMPONENT_LOOP(cgh, grouptype) {
          
          ivect size(1);
          const int gpdim = groupdata.at(group).info.dim;
          for (int d=0; d<gpdim; ++d) {
            size[d] = groupdata.at(group).info.lsh[d];
          }
          const int np = prod(size);
          
          for (int var=0; var<nvar; ++var) {
            const int n = n0 + var;
            for (int tl=min_tl; tl<=max_tl; ++tl) {
              memset (cgh->data[n][tl], poison_value, np*sz);
            } // for tl
          } // for var
          
        } END_LOCAL_COMPONENT_LOOP;
      } END_MAP_LOOP;
      
    } // if tl
  }
  
  
  
  void PoisonCheck (const cGH* cgh, const checktimes where)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (! check_for_poison) return;
    
    Checkpoint ("PoisonCheck");
    
    for (int group=0; group<CCTK_NumGroups(); ++group) {
      if (CCTK_QueryGroupStorageI(cgh, group)) {
        
        const int grouptype = CCTK_GroupTypeI(group);
        const int n0 = CCTK_FirstVarIndexI(group);
        assert (n0>=0);
        const int nvar = CCTK_NumVarsInGroupI(group);
        const int tp = CCTK_VarTypeI(n0);
        const int gpdim = groupdata.at(group).info.dim;
        
        const int num_tl = CCTK_NumTimeLevelsFromVarI(n0);
        assert (num_tl>0);
        const int min_tl = mintl(where, num_tl);
        const int max_tl = maxtl(where, num_tl);
        
        BEGIN_MAP_LOOP(cgh, grouptype) {
          BEGIN_LOCAL_COMPONENT_LOOP(cgh, grouptype) {
            
            ivect size(1);
            for (int d=0; d<gpdim; ++d) {
              size[d] = groupdata.at(group).info.lsh[d];
            }
            const int np = prod(size);
            
            for (int var=0; var<nvar; ++var) {
              const int n = n0 + var;
              
              for (int tl=min_tl; tl<=max_tl; ++tl) {
                
                const void* const data = cgh->data[n][tl];
                int numpoison=0;
                for (int k=0; k<size[2]; ++k) {
                  for (int j=0; j<size[1]; ++j) {
                    for (int i=0; i<size[0]; ++i) {
                      const int idx = i + size[0] * (j + size[1] * k);
                      bool poisoned=false;
                      switch (tp) {
#define TYPECASE(N,T)                                                   \
                      case N: {                                         \
                        T worm;                                         \
                        memset (&worm, poison_value, sizeof worm);      \
                        const T & val = ((const T*)data)[idx];          \
                        poisoned = memcmp (&worm, &val, sizeof worm) == 0; \
                        break;                                          \
                      }
#include "typecase"
#undef TYPECASE
                      default:
                        UnsupportedVarType(n);
                      }
                      if (poisoned) {
                        ++numpoison;
                        if (max_poison_locations==-1
                            || numpoison<=max_poison_locations) {
                          char* fullname = CCTK_FullName(n);
                          CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                                      "Timelevel %d, component %d, map %d, refinement level %d of the variable \"%s\" contains poison at [%d,%d,%d]",
                                      tl, component, map, reflevel, fullname, i,j,k);
                          free (fullname);
                        }
                      } // if poisoned
                    } // for i
                  } // for j
                } // for k
                if (max_poison_locations!=-1 && numpoison>max_poison_locations) {
                  char* fullname = CCTK_FullName(n);
                  CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                              "Timelevel %d, component %d, map %d, refinement level %d of the variable \"%s\" contains poison at %d of %d locations; not all locations were printed",
                              tl, component, map, reflevel, fullname, numpoison, np);
                  free (fullname);
                } else if (numpoison>0) {
                  CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                              "Found poison at %d of %d locations",
                              numpoison, np);
                }
                
              } // for tl
            } // for var
          } END_LOCAL_COMPONENT_LOOP;
        } END_MAP_LOOP;
        
      } // if has storage
    } // for group
  }
  
} // namespace Carpet