aboutsummaryrefslogtreecommitdiff
path: root/src/SetTime.c
blob: 94380b130627c79d2522801b8472889561abeb47 (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
 /*@@
   @file      SetTime.c
   @date      Mon May 20 09:45:45 2002
   @author    Ian Hawke
   @desc 
   Sets the time and dt depending on the ODE method and 
   position in the loop.
   @enddesc 
   @version   $Header$
 @@*/

#include <stdlib.h>

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

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

CCTK_FILEVERSION(CactusBase_MoL_SetTime_c);

/* #define MOLDEBUG 1 */

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

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

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

int MoL_SetTime(CCTK_ARGUMENTS);

int MoL_ResetTime(CCTK_ARGUMENTS);

int MoL_ResetDeltaTime(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_SetTime
   @date       Mon May 20 09:48:55 2002
   @author     Ian Hawke
   @desc 
   Sets the time and timestep before the MoL evolution loop starts.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int MoL_SetTime(CCTK_ARGUMENTS)
{
  
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_REAL beta;

  *Original_Time = cctkGH->cctk_time;
  *Original_Delta_Time = cctkGH->cctk_delta_time;
  cctkGH->cctk_time -= cctkGH->cctk_delta_time / cctkGH->cctk_timefac;

  if (CCTK_EQUALS(ODE_Method,"ICN"))
  {
    cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"ICN-avg"))
  {
    cctkGH->cctk_delta_time = *Original_Delta_Time;
  }
  else if (CCTK_EQUALS(ODE_Method,"Generic"))
  {
    beta = RKBetaCoefficients[0];

    cctkGH->cctk_delta_time = beta*(*Original_Delta_Time);
  }
  
  return 0;
}

 /*@@
   @routine    MoL_ResetTime
   @date       Mon May 20 09:49:41 2002
   @author     Ian Hawke
   @desc 
   Sets the time during the MoL evolution loop.
   At the last time all methods should end up with the original
   values for time and timestep.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int MoL_ResetTime(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  CCTK_REAL *previous_times;
  CCTK_INT alphaindex, i, j;
    
  previous_times = (CCTK_REAL*)malloc(MoL_Intermediate_Steps*sizeof(CCTK_REAL));
  if (!previous_times)
  {
    CCTK_WARN(0,"Failed to allocate memory for very small array!");
  }

  if (*MoL_Intermediate_Step == 0) 
  {
    cctkGH->cctk_time = (*Original_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"ICN"))
  {
    cctkGH->cctk_time = (*Original_Time) - 
      0.5*(*Original_Delta_Time)/cctkGH->cctk_timefac;
  }
  else if (CCTK_EQUALS(ODE_Method,"ICN-avg"))
  {
    cctkGH->cctk_time = (*Original_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"Generic"))
  {
    previous_times[0] = (*Original_Time) -
      (*Original_Delta_Time)/cctkGH->cctk_timefac;
    for (i = MoL_Intermediate_Steps - 1; i > *MoL_Intermediate_Step - 1; i--)
    {
      previous_times[MoL_Intermediate_Steps - i] = 
        RKBetaCoefficients[MoL_Intermediate_Steps - i - 1] * 
        (*Original_Delta_Time)/cctkGH->cctk_timefac;
      for (j = MoL_Intermediate_Steps; j > i; j--)
      {
        alphaindex = (MoL_Intermediate_Steps - i - 1) * 
          MoL_Intermediate_Steps + MoL_Intermediate_Steps - j;
        previous_times[MoL_Intermediate_Steps - i] += 
          RKAlphaCoefficients[alphaindex] *
          previous_times[MoL_Intermediate_Steps - j];
#ifdef MOLDEBUG
        printf("i %d j %d is %d index %d t %g dt %g alpha %g beta %g\n",
               i, j, MoL_Intermediate_Steps, alphaindex,
               previous_times[MoL_Intermediate_Steps - i], 
               (*Original_Delta_Time)/cctkGH->cctk_timefac,
               RKAlphaCoefficients[alphaindex], 
               RKBetaCoefficients[MoL_Intermediate_Steps - i - 1]);
#endif
      }
    }
#ifdef MOLDEBUG
    printf("MoL says the previous times are ");
    for (i = 0; i < MoL_Intermediate_Steps - *MoL_Intermediate_Step
           + 1; i++)
    {
      printf("%g ", previous_times[i]);
    }
    printf("\n");
#endif
    cctkGH->cctk_time = previous_times[MoL_Intermediate_Steps - 
                                      *MoL_Intermediate_Step];
  }
  else if (CCTK_EQUALS(ODE_Method,"RK2"))
  {
    if (*MoL_Intermediate_Step == 1)
    {
      cctkGH->cctk_time = (*Original_Time);
    }
  }
  else if (CCTK_EQUALS(ODE_Method,"RK3"))
  {
    if (*MoL_Intermediate_Step == 2)
    {
      cctkGH->cctk_time = (*Original_Time);
    }
    else if (*MoL_Intermediate_Step == 1)
    {
      cctkGH->cctk_time = (*Original_Time) - 
        0.5*(*Original_Delta_Time)/cctkGH->cctk_timefac;
    }
  }
#ifdef MOLDEBUG
  printf("MoL has once more reset t (%d): %f.\n", 
         *MoL_Intermediate_Step, cctkGH->cctk_time);
  fflush(stdout);
#endif

  free(previous_times);
  previous_times = NULL;

  return 0;
}

 /*@@
   @routine    MoL_ResetDeltaTime
   @date       Mon May 20 09:49:41 2002
   @author     Ian Hawke
   @desc 
   Sets the timestep during the MoL evolution loop.
   At the last time all methods should end up with the original
   values for time and timestep.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;

  if (*MoL_Intermediate_Step == 0) 
  {
    cctkGH->cctk_delta_time = (*Original_Delta_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"ICN"))
  {
    if (*MoL_Intermediate_Step == 1)
    {
      cctkGH->cctk_delta_time = (*Original_Delta_Time);
    }
    else
    {
      cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
    }
  }
  else if (CCTK_EQUALS(ODE_Method,"ICN-avg"))
  {
    cctkGH->cctk_delta_time = (*Original_Delta_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"Generic"))
  {
    cctkGH->cctk_delta_time = RKBetaCoefficients[MoL_Intermediate_Steps - 
                                                (*MoL_Intermediate_Step)] *
      (*Original_Delta_Time);
  }
  else if (CCTK_EQUALS(ODE_Method,"RK2"))
  {
    if (*MoL_Intermediate_Step == 1)
    {
      cctkGH->cctk_delta_time = 0.5*(*Original_Delta_Time);
    }
  }
  else if (CCTK_EQUALS(ODE_Method,"RK3"))
  {
    if (*MoL_Intermediate_Step == 2)
    {
      cctkGH->cctk_delta_time = 0.25*(*Original_Delta_Time);
    }
    else if (*MoL_Intermediate_Step == 1)
    {
      cctkGH->cctk_delta_time = 2.0/3.0*(*Original_Delta_Time);
    }
  }
#ifdef MOLDEBUG
  printf("MoL has once more reset dt (%d): %f.\n", 
         *MoL_Intermediate_Step, 
         cctkGH->cctk_delta_time/cctkGH->cctk_timefac);
  fflush(stdout);
#endif

  return 0;
}

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