aboutsummaryrefslogtreecommitdiff
path: root/src/StepSize.c
blob: c32812caba72098bf99069804ac784c38a4f8abb (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
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
 /*@@
   @file      StepSize.c
   @date      Tue Sep 07 2004
   @author    Erik Schnetter
   @desc 
   Control the time step size.
   @enddesc 
   @version   $Header$
 @@*/

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

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

#include "ExternalVariables.h"

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

CCTK_FILEVERSION(CactusBase_MoL_StepSize_c);

/********************************************************************
 *********************     Local Data Types   ***********************
 ********************************************************************/

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

/********************************************************************
 ***************** Scheduled Routine Prototypes *********************
 ********************************************************************/

void MoL_StartLoop(CCTK_ARGUMENTS);

void MoL_InitAdaptiveError(CCTK_ARGUMENTS);
void MoL_FindAdaptiveError(CCTK_ARGUMENTS);
void MoL_ReduceAdaptiveError(CCTK_ARGUMENTS);

void MoL_SetEstimatedDt(CCTK_ARGUMENTS);

void MoL_FinishLoop(CCTK_ARGUMENTS);

/********************************************************************
 ********************* Other Routine Prototypes *********************
 ********************************************************************/

/********************************************************************
 *********************     Local Data   *****************************
 ********************************************************************/

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

 /*@@
   @routine    MoL_StartLoop
   @date       Tue Sep 07 2004
   @author     Erik Schnetter
   @desc 
   Start the step size control loop, so that at least one iteration is done.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void
MoL_StartLoop(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  *MoL_Stepsize_Bad = 1;

  if (adaptive_stepsize)
  {
    *EstimatedDt = cctkGH->cctk_delta_time;
  }
  
}

 /*@@
   @routine    MoL_InitAdaptiveError
   @date       Tue Sep 07 2004
   @author     Erik Schnetter
   @desc 
   Initialize error counters for adaptive stepsize control
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

static inline CCTK_REAL
square (CCTK_REAL const x)
{
  return x * x;
}

void MoL_InitAdaptiveError(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  /* Initialise global error */
  *Error = 0;
  *Count = 0;
}

 /*@@
   @routine    MoL_FindAdaptiveError
   @date       Thu Jan 27 10:22:26 2005
   @author     Erik Schnetter
   @desc 
   Compute the error local to this component/patch/...
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_FindAdaptiveError(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  /* TODO: exclude symmetry boundaries */
  assert (cctk_dim <= 3);
  int imin[3], imax[3];
  for (int d = 0; d < cctk_dim; d++)
  {
    imin[d] = cctk_bbox[2*d] ? 0 : cctk_nghostzones[d];
    imax[d] = cctk_lsh[d] - (cctk_bbox[2*d+1] ? 0 : cctk_nghostzones[d]);
  }

  /* Calculate absolute error */
  CCTK_REAL local_error = 0.0;
  for (int var = 0; var < MoLNumEvolvedVariables; var++)
  {

    CCTK_REAL const * restrict const
      UpdateVar = CCTK_VarDataPtrI(cctkGH, 0, EvolvedVariableIndex[var]);
    CCTK_REAL const * restrict const
      RHSVar = CCTK_VarDataPtrI(cctkGH, 0, RHSVariableIndex[var]);
    CCTK_REAL const * restrict const
      ErrorVar
      = CCTK_VarDataPtrI(cctkGH, 0, 
                         CCTK_FirstVarIndex("MOL::ERRORESTIMATE") + var);

    CCTK_REAL const rhs_relative_error
      = maximum_relative_error * RHS_error_weight * (*Original_Delta_Time);

    assert (cctk_dim == 3);
#pragma omp parallel for reduction(+: local_error)
    for (int k = imin[2]; k < imax[2]; k++)
    {
      for (int j = imin[1]; j < imax[1]; j++)
      {
        for (int i = imin[0]; i < imax[0]; i++)
        {
          int const index = CCTK_GFINDEX3D(cctkGH, i, j, k);
          CCTK_REAL const scale
            = (square(maximum_absolute_error)
               + square(maximum_relative_error * UpdateVar[index])
               + square(rhs_relative_error * RHSVar[index]));
          local_error += square(ErrorVar[index]) / scale;
        }
      }
    }

  } /* for var */

  *Error += local_error;
  *Count
    += (MoLNumEvolvedVariables
        * (imax[0] - imin[0]) * (imax[1] - imin[1]) * (imax[2] - imin[2]));
}

 /*@@
   @routine    MoL_ReduceAdaptiveError
   @date       Thu Jan 27 10:23:14 2005
   @author     Erik Schnetter
   @desc 
   Find the global error estimate. 
   Change the timestep based on the error estimate.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_ReduceAdaptiveError(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  int redop;

  CCTK_REAL red_local[2], red_global[2];
  CCTK_REAL p1, p2;
  int ierr;

  /* Get global result over all processors */
  redop = CCTK_ReductionHandle ("sum");
  assert (redop >= 0);

  red_local[0] = *Error;
  red_local[1] = *Count;
  ierr = CCTK_ReduceLocArrayToArray1D
    (cctkGH, -1, redop, red_local, red_global, 2, CCTK_VARIABLE_REAL);
  assert (ierr == 0);
  *Error = red_global[0];
  *Count = red_global[1];

  /* Calculate L2-norm */
  *Error = sqrt(*Error / *Count);
  if (! CCTK_EQUALS(verbose, "none"))
  {
    CCTK_VInfo (CCTK_THORNSTRING, "Integration accuracy quotient is %g", (double)*Error);
  }
  if (! isfinite(*Error))
  {
    CCTK_VWarn (CCTK_WARN_ALERT,__LINE__,__FILE__,CCTK_THORNSTRING,
                "Integration accuracy quotient is %g, which is not a finite number -- reducing the step size", (double)*Error);
  }

  if ( CCTK_EQUALS(ODE_Method,"RK45") )
  {
    p1 = 1.0/5.0;
    p2 = 1.0/4.0;
  }
  else if ( CCTK_EQUALS(ODE_Method,"RK65") )
  {
    p1 = 1.0/6.0;
    p2 = 1.0/5.0;
  }
  else if ( CCTK_EQUALS(ODE_Method,"RK87") )
  {
    p1 = 1.0/8.0;
    p2 = 1.0/7.0;
  }
  else
  {
    CCTK_WARN (CCTK_WARN_ABORT, "unsupported ODE_Method in stepsize control");
    /* Avoid compiler warnings */
    p1 = 0;
    p2 = 0;
  }

  /* Decide whether to accept this step */
  *MoL_Stepsize_Bad = ! isfinite(*Error) || *Error > 1;

  if (*MoL_Stepsize_Bad)
  {
    /* The error is too large; reject the time step and reduce the
       step size */
    cctkGH->cctk_time -= cctkGH->cctk_delta_time;
    if (! CCTK_EQUALS(verbose, "none"))
    {
      CCTK_VInfo (CCTK_THORNSTRING, "*** REJECTING TIME STEP ***");
    }

    if (isfinite(*Error))
    {
      cctkGH->cctk_delta_time
        = safety_factor * (*Original_Delta_Time) / pow(*Error, p1);
    }
    else
    {
      cctkGH->cctk_delta_time = (*Original_Delta_Time) / maximum_decrease;
    }
    /* if (! CCTK_EQUALS(verbose, "none")) */
    /* { */
    /*   CCTK_VInfo (CCTK_THORNSTRING, "Setting time step to %g", (double)cctkGH->cctk_delta_time); */
    /* } */

    if (cctkGH->cctk_delta_time < (*Original_Delta_Time) / maximum_decrease)
    {
      /* No more than a factor of 10 decrease */
      cctkGH->cctk_delta_time = (*Original_Delta_Time) / maximum_decrease;
      if (! CCTK_EQUALS(verbose, "none"))
      {
        CCTK_VInfo (CCTK_THORNSTRING, "   Time step reduction too large; clamping time step to %g", (double)cctkGH->cctk_delta_time);
      }
    }
    if (cctkGH->cctk_delta_time == (CCTK_REAL)0.0)
      /* yes, we want to compare to zero exactly, to catch underflows */
      {
        CCTK_WARN (CCTK_WARN_ABORT, "New step size would be zero -- aborting");
      }

    cctkGH->cctk_time += cctkGH->cctk_delta_time;
  }
  else
  {
    /* The error is acceptable; estimate the next step size */
    *EstimatedDt = 
      safety_factor * (*Original_Delta_Time) / pow(*Error, p2);

    if (*EstimatedDt > (*Original_Delta_Time) * maximum_increase)
    {
      /* No more than a factor of 5 increase */
      *EstimatedDt = (*Original_Delta_Time) * maximum_increase;
      if (! CCTK_EQUALS(verbose, "none"))
      {
        CCTK_VInfo (CCTK_THORNSTRING, "   Time step increase too large; clamping time step to %g", (double)(*EstimatedDt));
      }
    }
  }
}

 /*@@
   @routine    MoL_SetEstimatedDt
   @date       Thu Jan 27 14:05:08 2005
   @author     Ian Hawke
   @desc 
   Actually set the timestep in PostStep. 
   This avoids problems when the timestep is changed in the middle of the
   evolution loop.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_SetEstimatedDt(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  cctkGH->cctk_delta_time = *EstimatedDt;

  if (! CCTK_EQUALS(verbose, "none"))
  {
    CCTK_VInfo (CCTK_THORNSTRING, "Setting time step to %g", 
                (double)cctkGH->cctk_delta_time);
  }

}

 /*@@
   @routine    MoL_FinishLoop
   @date       Thu Jan 27 10:24:19 2005
   @author     Erik Schnetter
   @desc 
   Loop control if adaptive timestepping is not used.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_FinishLoop(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  /* Keep time step size unchanged */
  *MoL_Stepsize_Bad = 0;
  
  // Set these flags to ONE outside of MoL-loop!
  *MoL_SlowPostStep = 1;
  *MoL_SlowStep = 1;
}

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