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

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

#include "carpet.hh"

static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Poison.cc,v 1.3 2001/11/02 10:58:58 schnetter Exp $";



namespace Carpet {
  
  using namespace std;
  
  
  
  void Poison (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*)cgh, group)) {
	PoisonGroup (cgh, group, where);
      } // if has storage
    } // for group
  }
  
  
  
  void PoisonGroup (cGH* cgh, const int group, const checktimes where)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (group>=0 && group<CCTK_NumGroups());
    
    if (! poison_new_timelevels) return;
    
    Checkpoint ("%*sPoisonGroup %s", 2*reflevel, "", CCTK_GroupName(group));
    
    if (! CCTK_QueryGroupStorageI((cGH*)cgh, group)) {
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
		  "Cannot poison group \"%s\" because it has no storage",
		  CCTK_GroupName(group));
      return;
    }
    
    for (int var=0; var<CCTK_NumVarsInGroupI(group); ++var) {
      
      const int n = CCTK_FirstVarIndexI(group) + var;
      const int sz = CCTK_VarTypeSize(CCTK_VarTypeI(n));
      assert (sz>0);
      
      const int num_tl = CCTK_NumTimeLevelsFromVarI(n);
      assert (num_tl>0);
      const int min_tl = mintl(where, num_tl);
      const int max_tl = maxtl(where, num_tl);
      
      for (int tl=min_tl; tl<=max_tl; ++tl) {
	
	BEGIN_COMPONENT_LOOP(cgh) {
	  if (hh->is_local(reflevel,component)) {
	    int np = 1;
	    const int gpdim = arrdata[group].info.dim;
	    for (int d=0; d<gpdim; ++d) {
	      np *= *CCTK_ArrayGroupSizeI((cGH*)cgh, d, group);
	    }
	    memset (cgh->data[n][tl], poison_value, np*sz);
	  }
	} END_COMPONENT_LOOP(cgh);
	
      } // for tl
      
    } // for var
  }
  
  
  
  void PoisonCheck (cGH* cgh, const checktimes where)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (! check_for_poison) return;
    
    Checkpoint ("%*sPoisonCheck", 2*reflevel, "");
    
    for (int group=0; group<CCTK_NumGroups(); ++group) {
      if (CCTK_QueryGroupStorageI((cGH*)cgh, group)) {
	for (int var=0; var<CCTK_NumVarsInGroupI(group); ++var) {
	  
	  const int n = CCTK_FirstVarIndexI(group) + var;
	  
	  const int num_tl = CCTK_NumTimeLevelsFromVarI(n);
	  assert (num_tl>0);
	  const int min_tl = mintl(where, num_tl);
	  const int max_tl = maxtl(where, num_tl);
	  
	  for (int tl=min_tl; tl<=max_tl; ++tl) {
	    
	    BEGIN_COMPONENT_LOOP(cgh) {
	      if (hh->is_local(reflevel,component)) {
		vect<int,dim> size(1);
		const int gpdim = arrdata[group].info.dim;
		for (int d=0; d<gpdim; ++d) {
		  size[d] = *CCTK_ArrayGroupSizeI((cGH*)cgh, d, group);
		}
		const int tp = CCTK_VarTypeI(n);
		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 = CCTK_GFINDEX3D(cgh,i,j,k);
		      bool poisoned=false;
		      switch (tp) {
#define TYPECASE(N,T)							\
		      case N: {						\
			T worm;						\
			memset (&worm, poison_value, sizeof(worm));	\
			poisoned = ((const T*)data)[idx] == worm;	\
			break;						\
		      }
#include "typecase"
#undef TYPECASE
		      default:
			UnsupportedVarType(n);
		      }
		      if (poisoned) {
			++numpoison;
			if (numpoison<=10) {
			  char* fullname = CCTK_FullName(n);
			  CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
				      "The variable \"%s\" contains poison at [%d,%d,%d] in timelevel %d",
				      fullname, i,j,k, tl);
			  free (fullname);
			}
		      } // if poisoned
		    } // for i
		  } // for j
		} // for k
		if (numpoison>10) {
		  char* fullname = CCTK_FullName(n);
		  CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
			      "The variable \"%s\" contains poison at %d locations in timelevel %d; not all locations were printed.",
			      fullname, numpoison, tl);
		  free (fullname);
		}
	      } // if is local
	    } END_COMPONENT_LOOP(cgh);
	    
	  } // for tl
	  
	} // for var
      } // if has storage
    } // for group
  }
  
} // namespace Carpet