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
|
#include <assert.h>
#include <stdlib.h>
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Termination.h"
#include "Carpet/CarpetLib/src/th.hh"
#include "carpet.hh"
static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Evolve.cc,v 1.7 2002/01/09 17:45:39 schnetter Exp $";
namespace Carpet {
using namespace std;
static bool do_terminate (const cGH *cgh, CCTK_REAL time, int iteration)
{
DECLARE_CCTK_PARAMETERS;
// Early shortcut
if (terminate_next || CCTK_TerminationReached(cgh)) return false;
bool term_iter = iteration >= cctk_itlast;
bool term_time = (cctk_initial_time < cctk_final_time
? time >= cctk_final_time
: time <= cctk_final_time);
if (CCTK_Equals(terminate, "never")) {
return true;
} else if (CCTK_Equals(terminate, "iteration")) {
return term_iter;
} else if (CCTK_Equals(terminate, "time")) {
return term_time;
} else if (CCTK_Equals(terminate, "either")) {
return term_iter || term_time;
} else if (CCTK_Equals(terminate, "both")) {
return term_iter && term_time;
} else {
abort();
}
}
int Evolve (tFleshConfig* fc)
{
DECLARE_CCTK_PARAMETERS;
Waypoint ("starting Evolve...");
const int convlev = 0;
cGH* cgh = fc->GH[convlev];
// Main loop
while (! do_terminate(cgh, cgh->cctk_time, cgh->cctk_iteration)) {
// Advance time
++cgh->cctk_iteration;
cgh->cctk_time += base_delta_time / maxreflevelfact;
Waypoint ("Evolving iteration %d...", cgh->cctk_iteration);
BEGIN_MGLEVEL_LOOP(cgh) {
if ((cgh->cctk_iteration-1) % mglevelfact == 0) {
Waypoint ("Evolving multigrid level %d...", mglevel);
BEGIN_REFLEVEL_LOOP(cgh) {
const int do_every = mglevelfact * maxreflevelfact/reflevelfact;
if ((cgh->cctk_iteration-1) % do_every == 0) {
// Cycle time levels
CycleTimeLevels (cgh);
// Advance level times
tt->advance_time (reflevel, mglevel);
tt0->advance_time (reflevel, mglevel);
for (int group=0; group<CCTK_NumGroups(); ++group) {
switch (CCTK_GroupTypeI(group)) {
case CCTK_SCALAR:
break;
case CCTK_ARRAY:
arrdata[group].tt->advance_time (reflevel, mglevel);
break;
case CCTK_GF:
break;
default:
abort();
}
}
// Checking
CalculateChecksums (cgh, allbutcurrenttime);
Poison (cgh, currenttimebutnotifonly);
// Evolve
Waypoint ("%*sScheduling PRESTEP", 2*reflevel, "");
CCTK_ScheduleTraverse ("CCTK_PRESTEP", cgh, CallFunction);
Waypoint ("%*sScheduling EVOL", 2*reflevel, "");
CCTK_ScheduleTraverse ("CCTK_EVOL", cgh, CallFunction);
Waypoint ("%*sScheduling POSTSTEP", 2*reflevel, "");
CCTK_ScheduleTraverse ("CCTK_POSTSTEP", cgh, CallFunction);
// Checking
PoisonCheck (cgh, currenttimebutnotifonly);
// Recompose grid hierarchy
Recompose (cgh);
}
} END_REFLEVEL_LOOP(cgh);
}
} END_MGLEVEL_LOOP(cgh);
BEGIN_MGLEVEL_LOOP(cgh) {
if (cgh->cctk_iteration % mglevelfact == 0) {
BEGIN_REVERSE_REFLEVEL_LOOP(cgh) {
const int do_every = mglevelfact * maxreflevelfact/reflevelfact;
if (cgh->cctk_iteration % do_every == 0) {
// Restrict
Restrict (cgh);
// Checking
CalculateChecksums (cgh, currenttime);
// Waypoint
Waypoint ("%*sScheduling CHECKPOINT", 2*reflevel, "");
CCTK_ScheduleTraverse ("CCTK_CHECKPOINT", cgh, CallFunction);
// Analysis
Waypoint ("%*sScheduling ANALYSIS", 2*reflevel, "");
CCTK_ScheduleTraverse ("CCTK_ANALYSIS", cgh, CallFunction);
// Output
Waypoint ("%*sOutputGH", 2*reflevel, "");
CCTK_OutputGH (cgh);
// Checking
CheckChecksums (cgh, alltimes);
}
} END_REVERSE_REFLEVEL_LOOP(cgh);
}
} END_MGLEVEL_LOOP(cgh);
} // main loop
Waypoint ("done with Evolve.");
return 0;
}
} // namespace Carpet
|