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
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
|
#include <assert.h>
#include <stdlib.h>
#include <iostream>
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctki_GHExtensions.h"
#include "cctki_ScheduleBindings.h"
#include "cctki_WarnLevel.h"
#include "carpet.hh"
extern "C" {
static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Initialise.cc,v 1.35 2004/01/25 14:57:27 schnetter Exp $";
CCTK_FILEVERSION(Carpet_Carpet_Initialise_cc);
}
namespace Carpet {
using namespace std;
int Initialise (tFleshConfig* fc)
{
DECLARE_CCTK_PARAMETERS;
// Initialise stuff
const int convlev = 0;
cGH* const cgh = CCTK_SetupGH (fc, convlev);
CCTKi_AddGH (fc, convlev, cgh);
// Delay checkpoint until MPI has been initialised
Waypoint ("Starting initialisation");
// Initialise stuff
cgh->cctk_iteration = 0;
global_time = cctk_initial_time;
delta_time = 1.0;
cgh->cctk_time = global_time;
cgh->cctk_delta_time = delta_time;
do_global_mode = true;
do_meta_mode = true;
// Enable storage and communtication
CCTKi_ScheduleGHInit (cgh);
// Initialise stuff
CCTKi_InitGHExtensions (cgh);
BEGIN_MGLEVEL_LOOP(cgh) {
do_global_mode = true;
do_meta_mode = mglevel==mglevels-1;
// Register coordinates
Checkpoint ("Scheduling CCTK_WRAGH");
CCTK_ScheduleTraverse ("CCTK_WRAGH", cgh, CallFunction);
// Check parameters
Checkpoint ("Scheduling PARAMCHECK");
CCTK_ScheduleTraverse ("CCTK_PARAMCHECK", cgh, CallFunction);
CCTKi_FinaliseParamWarn();
} END_MGLEVEL_LOOP;
if (fc->recovered) {
// if recovering
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
cgh->cctk_time = global_time;
// Set up the grids
Checkpoint ("Scheduling BASEGRID");
CCTK_ScheduleTraverse ("CCTK_BASEGRID", cgh, CallFunction);
// Recover
Checkpoint ("Scheduling RECOVER_VARIABLES");
CCTK_ScheduleTraverse ("CCTK_RECOVER_VARIABLES", cgh, CallFunction);
Checkpoint ("Scheduling POST_RECOVER_VARIABLES");
CCTK_ScheduleTraverse
("CCTK_POST_RECOVER_VARIABLES", cgh, CallFunction);
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
} // for rl
} else {
// if not reconvering
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
cgh->cctk_time = global_time;
Waypoint ("Initialisation at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Checking
Poison (cgh, alltimes);
// Set up the grids
Checkpoint ("Scheduling BASEGRID");
CCTK_ScheduleTraverse ("CCTK_BASEGRID", cgh, CallFunction);
const int num_tl = init_each_timelevel ? 3 : 1;
// Rewind
for (int m=0; m<maps; ++m) {
vtt[m]->set_delta
(reflevel, mglevel, - vtt[m]->get_delta (reflevel, mglevel));
FlipTimeLevels (cgh);
for (int tl=0; tl<num_tl; ++tl) {
vtt[m]->advance_time (reflevel, mglevel);
CycleTimeLevels (cgh);
}
vtt[m]->set_delta
(reflevel, mglevel, - vtt[m]->get_delta (reflevel, mglevel));
FlipTimeLevels (cgh);
}
const bool outer_do_global_mode = do_global_mode;
for (int tl=num_tl-1; tl>=0; --tl) {
do_global_mode = outer_do_global_mode && tl==0;
// Advance times
for (int m=0; m<maps; ++m) {
vtt[m]->advance_time (reflevel, mglevel);
}
cgh->cctk_time
= global_time - tl * delta_time * mglevelfact / reflevelfact;
CycleTimeLevels (cgh);
// Set up the initial data
Checkpoint ("Scheduling INITIAL");
CCTK_ScheduleTraverse ("CCTK_INITIAL", cgh, CallFunction);
} // for tl
do_global_mode = outer_do_global_mode;
// Checking
PoisonCheck (cgh, currenttime);
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
// Regrid
Checkpoint ("Regrid");
Regrid (cgh, rl, rl+1, prolongate_initial_data);
} // for rl
for (int rl=reflevels-1; rl>=0; --rl) {
BEGIN_REVERSE_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
// Restrict
Restrict (cgh);
Checkpoint ("Scheduling POSTRESTRICTINITIAL");
CCTK_ScheduleTraverse
("CCTK_POSTRESTRICTINITIAL", cgh, CallFunction);
// Poststep
Checkpoint ("Scheduling POSTINITIAL");
CCTK_ScheduleTraverse ("CCTK_POSTINITIAL", cgh, CallFunction);
Checkpoint ("Scheduling POSTSTEP");
CCTK_ScheduleTraverse ("CCTK_POSTSTEP", cgh, CallFunction);
// Checking
PoisonCheck (cgh, alltimes);
CalculateChecksums (cgh, allbutcurrenttime);
leave_level_mode (cgh);
} END_REVERSE_MGLEVEL_LOOP;
} // for rl
if (init_3_timelevels) {
// Use Scott Hawley's algorithm for getting two extra
// timelevels of data
Waypoint ("Initialising three timelevels");
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
// Advance times
for (int m=0; m<maps; ++m) {
vtt[m]->advance_time (reflevel, mglevel);
}
cgh->cctk_time
= global_time + delta_time * mglevelfact / reflevelfact;
CycleTimeLevels (cgh);
Waypoint ("Initialisation 3TL evolution I (a) (forwards) at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Checking
CalculateChecksums (cgh, allbutcurrenttime);
Poison (cgh, currenttimebutnotifonly);
// Evolve forward
Checkpoint ("Scheduling PRESTEP");
CCTK_ScheduleTraverse ("CCTK_PRESTEP", cgh, CallFunction);
Checkpoint ("Scheduling EVOL");
CCTK_ScheduleTraverse ("CCTK_EVOL", cgh, CallFunction);
// Checking
PoisonCheck (cgh, currenttime);
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
} // for rl
delta_time *= -1;
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
// Flip time levels
Waypoint ("Flipping timelevels");
FlipTimeLevels (cgh);
cgh->cctk_time
= global_time + delta_time * mglevelfact / reflevelfact;
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
} // for rl
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
Waypoint ("Initialisation 3TL evolution I (b) (backwards) at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Checking
CalculateChecksums (cgh, allbutcurrenttime);
Poison (cgh, currenttimebutnotifonly);
// Evolve backward
Checkpoint ("Scheduling PRESTEP");
CCTK_ScheduleTraverse ("CCTK_PRESTEP", cgh, CallFunction);
Checkpoint ("Scheduling EVOL");
CCTK_ScheduleTraverse ("CCTK_EVOL", cgh, CallFunction);
// Checking
PoisonCheck (cgh, alltimes);
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
} // for rl
Waypoint ("Hourglass structure in place");
// Evolve each level "backwards" one more timestep
// Starting with the finest level and proceeding to the coarsest
for (int rl=reflevels-1; rl>=0; --rl) {
BEGIN_REVERSE_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
Waypoint ("Initialisation 3TL evolution II (b) (backwards) at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Restrict
Restrict (cgh);
Checkpoint ("Scheduling POSTRESTRICT");
CCTK_ScheduleTraverse ("CCTK_POSTRESTRICT", cgh, CallFunction);
Checkpoint ("Scheduling POSTSTEP");
CCTK_ScheduleTraverse ("CCTK_POSTSTEP", cgh, CallFunction);
// Checking
PoisonCheck (cgh, alltimes);
// Advance times
for (int m=0; m<maps; ++m) {
vtt[m]->advance_time (reflevel, mglevel);
}
cgh->cctk_time
= global_time + 2 * delta_time * mglevelfact / reflevelfact;
CycleTimeLevels (cgh);
Waypoint ("Initialisation 3TL evolution I (c) (backwards) at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Checking
CalculateChecksums (cgh, allbutcurrenttime);
Poison (cgh, currenttimebutnotifonly);
// Evolve backward
Checkpoint ("Scheduling PRESTEP");
CCTK_ScheduleTraverse ("CCTK_PRESTEP", cgh, CallFunction);
Checkpoint ("Scheduling EVOL");
CCTK_ScheduleTraverse ("CCTK_EVOL", cgh, CallFunction);
Checkpoint ("Scheduling POSTSTEP");
CCTK_ScheduleTraverse ("CCTK_POSTSTEP", cgh, CallFunction);
// Checking
PoisonCheck (cgh, alltimes);
leave_level_mode (cgh);
} END_REVERSE_MGLEVEL_LOOP;
} // for rl
delta_time *= -1;
for (int rl=0; rl<reflevels; ++rl) {
BEGIN_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
// Flip time levels back
Waypoint ("Flipping timelevels back");
FlipTimeLevels (cgh);
// Invert level times back
for (int m=0; m<maps; ++m) {
vtt[m]->set_delta
(reflevel, mglevel, - vtt[m]->get_delta (reflevel, mglevel));
vtt[m]->advance_time (reflevel, mglevel);
vtt[m]->advance_time (reflevel, mglevel);
vtt[m]->set_delta
(reflevel, mglevel, - vtt[m]->get_delta (reflevel, mglevel));
}
cgh->cctk_time = global_time;
leave_level_mode (cgh);
} END_MGLEVEL_LOOP;
} // for rl
Waypoint ("Finished initialising three timelevels");
} // if init_3_timelevels
} // if not recovering
for (int rl=reflevels-1; rl>=0; --rl) {
BEGIN_REVERSE_MGLEVEL_LOOP(cgh) {
enter_level_mode (cgh, rl);
do_global_mode = reflevel==0;
do_meta_mode = do_global_mode && mglevel==mglevels-1;
Waypoint ("Initialisation II at iteration %d time %g%s%s",
cgh->cctk_iteration, (double)cgh->cctk_time,
(do_global_mode ? " (global)" : ""),
(do_meta_mode ? " (meta)" : ""));
// Checkpoint
Checkpoint ("Scheduling CPINITIAL");
CCTK_ScheduleTraverse ("CCTK_CPINITIAL", cgh, CallFunction);
// Analysis
Checkpoint ("Scheduling ANALYSIS");
CCTK_ScheduleTraverse ("CCTK_ANALYSIS", cgh, CallFunction);
// Output
Checkpoint ("OutputGH");
CCTK_OutputGH (cgh);
// Checking
PoisonCheck (cgh, alltimes);
CheckChecksums (cgh, allbutcurrenttime);
leave_level_mode (cgh);
} END_REVERSE_MGLEVEL_LOOP;
} // for rl
Waypoint ("Done with initialisation");
return 0;
}
} // namespace Carpet
|