summaryrefslogtreecommitdiff
path: root/src/main/CactusDefaultEvolve.c
blob: 9bb64cf837a51cc42ef2799eb7195dcb6cb89d56 (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
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
 /*@@
   @file      CactusDefaultEvolve.c
   @date      Thu Oct  8 17:28:46 1998
   @author    Tom Goodale
   @desc
              Default cactus evolution stuff.
   @enddesc
   @version   $Id$
 @@*/

/* #define DEBUG_CCTK 1 */

#include <stdio.h>
#include <stdlib.h>

#include "definethisthorn.h"

#include "cctk_Flesh.h"
#include "cctk_Parameters.h"
#include "cctk_Groups.h"
#include "cctk_WarnLevel.h"
#include "cctk_Termination.h"
#include "cctk_Main.h"
#include "cctk_Misc.h"
#include "cctk_IO.h"

#ifdef HAVE_TIME_GETTIMEOFDAY
#ifdef HAVE_SYS_TIME_H
#include <sys/time.h>
#endif
#include <unistd.h>
#endif

static const char *rcsid="$Header$";

CCTK_FILEVERSION(main_CactusDefaultEvolve_c);

/* Define some macros for convenience. */

#define ForallConvLevels(iteration, conv_level)                               \
        {                                                                     \
          int factor = 1;                                                     \
          for (conv_level = 0; conv_level < config->nGHs; conv_level++)       \
          {                                                                   \
            if (iteration % factor == 0)                                      \
            {

#define EndForallConvLevels                                                   \
            }                                                                 \
            factor *= 2;                                                      \
          }                                                                   \
        }

/********************************************************************
 ********************* Local Routine Prototypes *********************
 ********************************************************************/

static int DoneMainLoop (const cGH *GH, CCTK_REAL cctk_time, int iteration);
static void StepGH (cGH *GH);

/********************************************************************
 *********************     External Routines   **********************
 ********************************************************************/

int CactusDefaultEvolve (tFleshConfig *config);

/********************************************************************
 *********************     Local Routines   *************************
 ********************************************************************/

 /*@@
   @routine    CactusDefaultEvolve
   @date       Thu Oct  8 17:30:15 1998
   @author     Tom Goodale
   @desc
               The default cactus evolution routine.
   @enddesc
   @calls      CCTK_Traverse
               CCTK_OutputGH
               DoneMainLoop
               StepGH
   @history
   @hdate      Fri May 12 2000
   @hauthor    Thomas Radke
   @hdesc      Moved evolution loop termination check into DoneMainLoop()
   @endhistory

   @var        config
   @vdesc      pointer to flesh configuration structure
   @vtype      tFleshConfig *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               always returns 0
   @endreturndesc
@@*/
int CactusDefaultEvolve (tFleshConfig *config)
{
  int var;
  unsigned int convergence_level, iteration;

  /*** Call OUTPUT for this GH (this routine    ***/
  /*** checks if output is necessary) and makes ***/
  /*** a Traverse with CCTK_ANALYSIS      ***/
  iteration = CCTK_MainLoopIndex ();
  ForallConvLevels (iteration, convergence_level)
  {
    CCTK_Traverse (config->GH[convergence_level], "CCTK_ANALYSIS");
    CCTK_OutputGH (config->GH[convergence_level]);
  }
  EndForallConvLevels;

  while (! DoneMainLoop (config->GH[0], config->GH[0]->cctk_time, iteration))
  {
    if (iteration == 0)
    {
      /* Can only use CactusDefaultEvolve with one timelevel */
      for (var = CCTK_NumVars () - 1; var >= 0; var--)
      {
        if (CCTK_MaxTimeLevelsVI (var) > 1)
        {
          CCTK_VWarn (0,__LINE__,__FILE__,"Cactus",
                      "Variable '%s' has multiple timelevels, default Cactus "
                      "evolve routine cannot rotate",
                      CCTK_VarName (var));
        }
      }
    }
    /* HERE ROTATE TIMELEVELS FOR ALL CONVERGENCE LEVELS */

    iteration = CCTK_SetMainLoopIndex (++iteration);

    /* Step each convergence level */
    ForallConvLevels (CCTK_MainLoopIndex (), convergence_level)
    {
      StepGH (config->GH[convergence_level]);
      /*
      CCTK_InfoOutput (config->GH[convergence_level], convergence_level);
      */
    }
    EndForallConvLevels;

    /* Dump out checkpoint data on all levels */
    ForallConvLevels (iteration, convergence_level)
    {
      CCTK_Traverse (config->GH[convergence_level], "CCTK_CHECKPOINT");
    }
    EndForallConvLevels;

    /*** Call OUTPUT for this GH (this routine    ***/
    /*** checks if output is necessary) and makes ***/
    /*** an Traverse with CCTK_ANALYSIS      ***/
    ForallConvLevels (iteration, convergence_level)
    {
        CCTK_Traverse (config->GH[convergence_level], "CCTK_ANALYSIS");
        CCTK_OutputGH (config->GH[convergence_level]);
    }
    EndForallConvLevels;

  } /*** END OF MAIN ITERATION LOOP ***/

  return (0);
}


 /*@@
   @routine    DoneMainLoop
   @date       Fri May 12 2000
   @author     Thomas Radke
   @desc
               Check the termination conditions for the evolution loop
   @enddesc
   @calls      CCTK_TerminationReached
   @history
   @hdate      Thu 6 Nov 2002
   @hauthor    Thomas Radke
   @hdesc      Added max_runtime condition test
   @endhistory

   @var        GH
   @vdesc      pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        simulation_time
   @vdesc      current physical simulation time
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        iteration
   @vdesc      current iteration
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               true (1) or false (0) for done/not done with main loop
   @endreturndesc
@@*/
static int DoneMainLoop (const cGH *GH, CCTK_REAL simulation_time,int iteration)
{
  int max_iteration_reached, max_simulation_time_reached, max_runtime_reached;
  int retval;
#ifdef HAVE_TIME_GETTIMEOFDAY
  struct timeval runtime;
  static struct timeval starttime = {0, 0};
#endif
  DECLARE_CCTK_PARAMETERS


#ifdef HAVE_TIME_GETTIMEOFDAY
  /* on the first time through, get the start time */
  if (starttime.tv_sec == 0 && starttime.tv_usec == 0)
  {
    gettimeofday (&starttime, NULL);
  }
#endif

  retval = terminate_next || CCTK_TerminationReached (GH);
  if (! retval && ! CCTK_Equals (terminate, "never"))
  {
    max_iteration_reached = iteration >= cctk_itlast;

    if (cctk_initial_time < cctk_final_time)
    {
      max_simulation_time_reached = simulation_time >= cctk_final_time;
    }
    else
    {
      max_simulation_time_reached = simulation_time <= cctk_final_time;
    }

    max_runtime_reached = 0;
#ifdef HAVE_TIME_GETTIMEOFDAY
    if (max_runtime > 0)
    {
      /* get the elapsed runtime in minutes and compare with max_runtime */
      gettimeofday (&runtime, NULL);
      runtime.tv_sec -= starttime.tv_sec;
      max_runtime_reached = ((CCTK_REAL) runtime.tv_sec / 60.0) >= max_runtime;
    }
#endif

    if (CCTK_Equals (terminate, "iteration"))
    {
      retval = max_iteration_reached;
    }
    else if (CCTK_Equals (terminate, "time"))
    {
      retval = max_simulation_time_reached;
    }
    else if (CCTK_Equals (terminate, "runtime"))
    {
      retval = max_runtime_reached;
    }
    else if (CCTK_Equals (terminate, "any"))
    {
      retval = max_iteration_reached || max_simulation_time_reached ||
               max_runtime_reached;
    }
    else if (CCTK_Equals (terminate, "all"))
    {
      retval = max_iteration_reached && max_simulation_time_reached &&
               max_runtime_reached;
    }
    /* the following two conditions are deprecated in BETA14 */
    else if (CCTK_Equals (terminate, "either"))
    {
      retval = max_iteration_reached || max_simulation_time_reached;
    }
    else /* if (CCTK_Equals (terminate, "both")) */
    {
      retval = max_iteration_reached && max_simulation_time_reached;
    }
  }

  return (retval);
}


 /*@@
   @routine    StepGH
   @date       Fri Aug 14 12:39:49 1998
   @author     Gerd Lanfermann
   @desc
               The full set of routines used to execute all schedule points
               in the main iteration loop. Makes calls to the individual
               routines for each schedule point.
   @enddesc
   @calls      CCTK_Traverse

   @var        GH
   @vdesc      pointer to CCTK grid hierachy
   @vtype      cGH *
   @vio        in
   @endvar
 @@*/
static void StepGH (cGH *GH)
{
  GH->cctk_time = GH->cctk_time + GH->cctk_delta_time;
  GH->cctk_iteration++;

  CCTK_Traverse (GH, "CCTK_PRESTEP");
  CCTK_Traverse (GH, "CCTK_EVOL");
  CCTK_Traverse (GH, "CCTK_POSTSTEP");
}