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

/*#define DEBUG_CCTK*/

#include <stdio.h>

#include "cctk_Flesh.h"
#include "cctk.h"
#include "cctk_parameters.h"
#include "rfrConstants.h"
#include "rfrInterface.h"

#include "CactusIOFunctions.h"

static char *rcsid="$Header$";

/* 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;                  \
                                        };                             \
                                     }

 /* Quick stuff for testing purposes. */
#define EVOLUTION 1
#define OUTPUT    2
int cactus_terminate;
static int cactus_terminate_global = 0;
#define TERMINATION_RAISED_BRDCAST 4

/* Local function prototypes. */

int CCTK_StepGH(cGH *GH);
 

/* the iteration counter used in the evolution loop */
static int iteration = 0;


 /*@@
   @routine    CCTK_SetMainLoopIndex
   @date       Sep 22 1999
   @author     Thomas Radke
   @desc 
               Sets the iteration counter variable of the evolution loop.
               This is used for recovery.
   @enddesc 
   @calls     
   @calledby   

@@*/
int CCTK_SetMainLoopIndex (int main_loop_index)
{
  iteration = main_loop_index;
  return 0;
}


 /*@@
   @routine    CCTK_MainLoopIndex
   @date       Sep 22 1999
   @author     Thomas Radke
   @desc 
               Returns the iteration counter variable of the evolution loop.
               This is used for checkpointing.
   @enddesc 
   @calls     
   @calledby   

@@*/
int CCTK_MainLoopIndex (void)
{
  return (iteration);
}


 /*@@
   @routine    CactusDefaultEvolve
   @date       Thu Oct  8 17:30:15 1998
   @author     Tom Goodale
   @desc 
   The default cactus evolution routine.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
int CactusDefaultEvolve(tFleshConfig *config)
{

  DECLARE_CCTK_PARAMETERS

  int convergence_level;

#ifdef DEBUG_CCTK
  CCTK_PRINTSEPARATOR
  printf("In CactusDefaultEvolve\n----------------------\n");
  printf("  Initializing iteration = %d\n",iteration);
  CCTK_PRINTSEPARATOR
#endif

#if 0
  CactusStartTimer(config->timer[OUTPUT]);
#endif

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

#if 0      
  CactusStopTimer(config->timer[OUTPUT]);

  CactusStartTimer(config->timer[EVOLUTION]);
#endif
  /*
  CCTK_InfoHeader(config);
  */


  while (iteration<cctk_itlast || (cctk_final_time>cctk_initial_time?config->GH[0]->cctk_time<cctk_final_time:0)) 
  {

#ifdef DEBUG_CCTK
    CCTK_PRINTSEPARATOR
    printf("In CactusDefaultEvolve\n----------------------\n");
    printf("  Advancing iteration %d = %d + 1\n",iteration+1,
         iteration); 
    CCTK_PRINTSEPARATOR
#endif

    iteration++;

    /* Step each convergence level */


    ForallConvLevels(iteration, convergence_level)
    {

      CCTK_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_rfrTraverse(config->GH[convergence_level],CCTK_CHECKPOINT);
    }
    EndForallConvLevels;

    /* Output perhaps */
#if 0
    CactusStartTimer(config->timer[OUTPUT]);
#endif
    /*** Call OUTPUT for this GH (this routine    ***/
    /*** checks if output is necessary) and makes ***/
    /*** an rfrTraverse with CCTK_ANALYSIS      ***/
    ForallConvLevels(iteration, convergence_level)
    {
        CCTK_rfrTraverse(config->GH[convergence_level],CCTK_ANALYSIS);
        CCTK_OutputGH(config->GH[convergence_level]);
    }
    EndForallConvLevels;
      
#if 0
    CactusStopTimer(config->timer[OUTPUT]);
#endif

#if 0
    ConvergenceReport(config->GH, iteration);

    TerminationStepper(config->GH[0]);

#endif
    /* Termination has been raised and broadcasted, exit loop*/
    if (cactus_terminate==TERMINATION_RAISED_BRDCAST) break;

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

#if 0
  CactusStopTimer(config->timer[EVOLUTION]);
#endif

  return 0;
}

/************************************************************************/

/* The following routines have been nicked from 3.0 for the moment. */





 /*@@
   @routine    CStepper
   @date       Fri Aug 14 12:39:49 1998
   @author     Gerd Lanfermann
   @desc 
     The full set of routines used to execute all rfr steps 
     int the main iteration loop. Makes calls to the individual 
     routines for each rfr step.
   @enddesc 
   @calls  PreStepper, EvolStepper, PostStepper
   @calledby main   
 @@*/

int CCTK_StepGH(cGH *GH) 
{

  void PreStepper(cGH *GH);
  void EvolStepper(cGH *GH);
  void PostStepper(cGH *GH);

  /* Advance GH->iteration BEFORE evolving */
#ifdef DEBUG_CCTK
  CCTK_PRINTSEPARATOR
  printf("In CCTK_StepGH\n--------------\n");
  printf("  Advancing GH->iteration to %lu = %lu + 1\n",(GH->cctk_iteration+1),
         GH->cctk_iteration);
  CCTK_PRINTSEPARATOR
#endif

  GH->cctk_iteration++;

  PreStepper(GH);
  EvolStepper(GH);

  /* Advance GH->time AFTER evolving */
#ifdef DEBUG_CCTK
  CCTK_PRINTSEPARATOR
  printf("In CCTK_StepGH\n--------------\n");
  printf("  Advancing GH->cctk_time %f = %f + %f\n",GH->cctk_time+GH->cctk_delta_time,
         GH->cctk_time,GH->cctk_delta_time);
  CCTK_PRINTSEPARATOR
#endif

  GH->cctk_time = GH->cctk_time + GH->cctk_delta_time;

  PostStepper(GH);


  return 0;
}

 /*@@
   @routine    PreStepper
   @date       Fri Aug 14 12:43:20 1998
   @author     Gerd Lanfermann
   @desc 
     calls RFR-PRESTEP
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void PreStepper(cGH *GH) {
  int Rstep;  

  /* Call the rfr with CCTK_PRESTEP */
  CCTK_rfrTraverse(GH, CCTK_PRESTEP);
}
 /*@@
   @routine    EvolStepper
   @date       Fri Aug 14 12:44:00 1998
   @author     Gerd Lanfermann
   @desc 
     calls RFR-EVOLUTION, checks for nans, increases physical time
   @enddesc 
   @calls     
   @calledby   
   @history

   @endhistory
@@*/

void EvolStepper(cGH *GH) 
{

  /* Call the rfr with Evolution */
  CCTK_rfrTraverse(GH, CCTK_EVOL);
  /* after Evolution check for NANs */

}

 /*@@
   @routine    PostStepper
   @date       Fri Aug 14 12:45:39 1998
   @author     Gerd Lanfermann
   @desc 
     calls the routines rgeistered as CCTK_POSSTEPS
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/
void PostStepper(cGH *GH) {
  int Rstep;  
   /* Call the rfr with post step */
    CCTK_rfrTraverse(GH, CCTK_POSTSTEP); 
}

 /*@@
   @routine    TerminationStepper
   @date       Fri Aug 14 13:07:11 1998
   @author     Gerd Lanfermann
   @desc 
     catctus_terminate is a global variable with these values: 
     TERMINATION_NOT_RAISED    : not signaled yet (cactus_initial.c)
     TERMINATION_RAISED_LOCAL  : signaled on one PE, not reduced (MPI_LOR) 
                                 to all PEs yet (main.c)
     TERMINATION_RAISED_BRDCAST: reduced -> can now be used to terminate 
                                 (chkpnt_terminate.c) by RFR 
     the raised termiantion signal is caught on 1 PE only and has to be recduced
     on all PEs before a termination sequenced can be launched (I like that)
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void TerminationStepper(cGH *GH) {
  int cactus_terminate_global; 
  
  cactus_terminate_global=cactus_terminate;

#if 0

#ifdef MPI
  MPI_Allreduce(&cactus_terminate,&cactus_terminate_global,1,
                MPI_INT,MPI_LOR,GH->PUGH_COMM_WORLD);
#endif
#endif
  if (cactus_terminate_global) { 
    cactus_terminate=TERMINATION_RAISED_BRDCAST;
    printf("RECEIVED GLOBAL TERMINATION SIGNAL \n");
  }
  CCTK_rfrTraverse(GH,CCTK_TERMINATE);
}