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
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
|
// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetSlab/src/slab.cc,v 1.5 2003/05/12 16:25:11 schnetter Exp $
#include <assert.h>
#include <stdlib.h>
#include <string.h>
#include <vector>
#include "cctk.h"
#include "Carpet/CarpetLib/src/bbox.hh"
#include "Carpet/CarpetLib/src/bboxset.hh"
#include "Carpet/CarpetLib/src/dh.hh"
#include "Carpet/CarpetLib/src/gdata.hh"
#include "Carpet/CarpetLib/src/gh.hh"
#include "Carpet/CarpetLib/src/ggf.hh"
#include "Carpet/CarpetLib/src/vect.hh"
#include "Carpet/Carpet/src/carpet.hh"
#include "slab.hh"
extern "C" {
static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetSlab/src/slab.cc,v 1.5 2003/05/12 16:25:11 schnetter Exp $";
CCTK_FILEVERSION(Carpet_CarpetSlab_slab_cc);
}
namespace CarpetSlab {
using namespace Carpet;
void* GetSlab (cGH* const cgh,
const int dest_proc,
const int n,
const int ti,
const int hdim,
const int origin[/*vdim*/],
const int dirs[/*hdim*/],
const int stride[/*hdim*/],
const int length[/*hdim*/])
{
if (reflevel == -1) {
CCTK_WARN (0, "It is not possible to use hyperslabbing in global mode");
}
// Check global state
assert (reflevel>=0);
assert (mglevel>=0);
// Save global state
int saved_component = component;
component = -1;
// Check Cactus grid hierarchy
assert (cgh);
// Check destination processor
assert (dest_proc>=-1 && dest_proc<CCTK_nProcs(cgh));
// Check variable index
assert (n>=0 && n<CCTK_NumVars());
// Get info about variable
const int group = CCTK_GroupIndexFromVarI(n);
assert (group>=0);
const int n0 = CCTK_FirstVarIndexI(group);
assert (n0>=0);
const int var = n - n0;
assert (var>=0);
// Get info about group
cGroup gp;
CCTK_GroupData (group, &gp);
assert (gp.dim<=dim);
assert (CCTK_QueryGroupStorageI(cgh, group));
const int typesize = CCTK_VarTypeSize(gp.vartype);
assert (typesize>0);
// Check dimension
assert (hdim>=0 && hdim<=gp.dim);
// Check timelevel
const int num_tl = gp.numtimelevels;
assert (ti>=0 && ti<num_tl);
const int tl = -ti;
// Check origin
// for (int d=0; d<dim; ++d) {
// assert (origin[d]>=0 && origin[d]<=sizes[d]);
// }
// Check directions
for (int dd=0; dd<hdim; ++dd) {
assert (dirs[dd]>=1 && dirs[dd]<=dim);
}
// Check stride
for (int dd=0; dd<hdim; ++dd) {
assert (stride[dd]>0);
}
// Check length
for (int dd=0; dd<hdim; ++dd) {
assert (length[dd]>=0);
}
// Check extent
// for (int dd=0; dd<hdim; ++dd) {
// assert (origin[dirs[dd]-1] + length[dd] <= sizes[dirs[dd]]);
// }
// Get insider information about variable
const gh<dim>* myhh;
const dh<dim>* mydd;
const ggf<dim>* myff;
assert (group < (int)arrdata.size());
myhh = arrdata[group].hh;
assert (myhh);
mydd = arrdata[group].dd;
assert (mydd);
assert (var < (int)arrdata[group].data.size());
myff = arrdata[group].data[var];
assert (myff);
// Detemine collecting processor
const int collect_proc = dest_proc<0 ? 0 : dest_proc;
// Determine own rank
const int rank = CCTK_MyProc(cgh);
// Calculate global size
int totalsize = 1;
for (int dd=0; dd<hdim; ++dd) {
totalsize *= length[dd];
}
// Allocate memory
void* hdata = 0;
if (dest_proc==-1 || rank==dest_proc) {
hdata = malloc(totalsize * typesize);
assert (hdata);
memset (hdata, 0, totalsize * typesize);
}
if (hh->components(reflevel) > 0) {
// Only temporarily
component = 0;
// Get sample data
const gdata<dim>* mydata;
mydata = (*myff)(tl, reflevel, component, mglevel);
// Stride of data in memory
const vect<int,dim> str = mydata->extent().stride();
// Stride of collected data
vect<int,dim> hstr = str;
for (int dd=0; dd<hdim; ++dd) {
hstr[dirs[dd]-1] *= stride[dd];
}
// Lower bound of collected data
vect<int,dim> hlb(0);
for (int d=0; d<gp.dim; ++d) {
hlb[d] = origin[d] * str[d];
}
// Upper bound of collected data
vect<int,dim> hub = hlb;
for (int dd=0; dd<hdim; ++dd) {
hub[dirs[dd]-1] += (length[dd]-1) * hstr[dirs[dd]-1];
}
// Calculate extent to collect
const bbox<int,dim> hextent (hlb, hub, hstr);
assert (hextent.num_points() == totalsize);
// Create collector data object
void* myhdata = rank==collect_proc ? hdata : 0;
gdata<dim>* const alldata = mydata->make_typed();
alldata->allocate (hextent, collect_proc, myhdata);
// Done with the temporary stuff
mydata = 0;
component = -1;
// Loop over all components, copying data from them
assert (component == -1);
for (component=0; component<hh->components(reflevel); ++component) {
// Get data object
mydata = (*myff)(tl, reflevel, component, mglevel);
// Calculate overlapping extents
const bboxset<int,dim> myextents
= ((mydd->boxes[reflevel][component][mglevel].sync_not
| mydd->boxes[reflevel][component][mglevel].interior)
& hextent);
// Loop over overlapping extents
for (bboxset<int,dim>::const_iterator ext_iter = myextents.begin();
ext_iter != myextents.end();
++ext_iter) {
// Copy data
alldata->copy_from (mydata, *ext_iter);
}
} // Loop over components
component = -1;
// Copy result to all processors
if (dest_proc == -1) {
for (int proc=0; proc<CCTK_nProcs(cgh); ++proc) {
if (proc != collect_proc) {
void* myhdata = rank==proc ? hdata : 0;
gdata<dim>* const tmpdata = mydata->make_typed();
tmpdata->allocate (alldata->extent(), proc, myhdata);
tmpdata->copy_from (alldata, alldata->extent());
delete tmpdata;
}
}
} // Copy result
delete alldata;
} // if components>0
// Restore global state
component = saved_component;
// Success
return hdata;
}
int Hyperslab_GetHyperslab (cGH* const GH,
const int target_proc,
const int vindex,
const int vtimelvl,
const int hdim,
const int global_startpoint [/*vdim*/],
const int directions [/*vdim*/],
const int lengths [/*hdim*/],
const int downsample [/*hdim*/],
void** const hdata,
int hsize [/*hdim*/])
{
const int gpdim = CCTK_GroupDimFromVarI(vindex);
assert (gpdim>=1 && gpdim<=dim);
// Check some arguments
assert (hdim>=0 && hdim<=dim);
// Check output arguments
assert (hdata);
assert (hsize);
// Calculate more convenient representation of the direction
vector<int> dirs(hdim);
// The following if statement is written according to the
// definition of "dir".
if (hdim==1) {
// 1-dimensional hyperslab
int mydir = 0;
for (int d=0; d<dim; ++d) {
if (directions[d]!=0) {
mydir = d+1;
break;
}
}
assert (mydir>0);
for (int d=0; d<dim; ++d) {
if (d == mydir-1) {
assert (directions[d]!=0);
} else {
assert (directions[d]==0);
}
}
dirs[0] = mydir;
} else if (hdim==dim) {
// dim-dimensional hyperslab
for (int dd=0; dd<hdim; ++dd) {
dirs[dd] = dd+1;
}
} else if (hdim==2) {
// 2-dimensional hyperslab with dim==3
assert (dim==3);
int mydir = 0;
for (int d=0; d<dim; ++d) {
if (directions[d]==0) {
mydir = d+1;
break;
}
}
assert (mydir>0);
for (int d=0; d<dim; ++d) {
if (d == mydir-1) {
assert (directions[d]==0);
} else {
assert (directions[d]!=0);
}
}
int dd=0;
for (int d=0; d<dim; ++d) {
if (d != mydir-1) {
dirs[dd] = d+1;
++dd;
}
}
assert (dd==hdim);
} else {
assert (0);
}
// Calculate lengths
for (int dd=0; dd<hdim; ++dd) {
if (lengths[dd]<0) {
int gsh[dim];
int ierr = CCTK_GroupgshVI(GH, dim, gsh, vindex);
assert (!ierr);
const int totlen = gsh[dirs[dd]-1];
assert (totlen>=0);
// Partial argument check
assert (global_startpoint[dirs[dd]-1]>=0);
assert (global_startpoint[dirs[dd]-1]<=totlen);
assert (downsample[dd]>0);
hsize[dd] = (totlen - global_startpoint[dirs[dd]-1]) / downsample[dd];
} else {
hsize[dd] = lengths[dd];
}
assert (hsize[dd]>=0);
}
// Get the slab
*hdata = GetSlab (GH,
target_proc,
vindex,
vtimelvl,
hdim,
global_startpoint,
&dirs[0],
downsample,
hsize);
// Return with success
return 1;
}
} // namespace CarpetSlab
|