aboutsummaryrefslogtreecommitdiff
path: root/src/RK2-MR-2_1.c
blob: 9179cdd0d543f75ac0ea304a53e4770e28094e68 (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
 /*@@
   @file      RK2-MR-2_1.c
   @date      2012-03-25
   @author    Christian Reisswig
   @desc 
   A routine to perform RK2 2:1 MR evolution. Mostly copied from
   genericRK.c
   @enddesc 
   @version   $Header$
 @@*/

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

#include <stdio.h>
#include "ExternalVariables.h"

//#define MOLDEBUG

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

CCTK_FILEVERSION(CactusBase_MoL_RK2_MR_2_1_c);

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

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

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

void MoL_RK2_MR_2_1_Add(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_RK2_MR_2_1_Add
   @date       
   @author     
   @desc 
   Performs a single step of a RK2_MR_2_1 type time
   integration.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_RK2_MR_2_1_Add(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
    
  CCTK_INT arraydim;

  /* FIXME */

#ifdef MOLDOESCOMPLEX

CCTK_WARN(0, "not implemented");
  CCTK_COMPLEX Complex_alpha, Complex_beta, Complex_Delta_Time;
  CCTK_COMPLEX * restrict OldComplexVar;
  CCTK_COMPLEX * restrict UpdateComplexVar;
  CCTK_COMPLEX const * restrict RHSComplexVar;
  CCTK_COMPLEX * restrict ScratchComplexVar;

#endif

  static CCTK_INT scratchspace_firstindex = -99;
  static CCTK_INT scratchspace_firstindex_slow = -99;
  CCTK_INT index, var, scratchstep;
  CCTK_INT totalsize;
  CCTK_REAL alpha[4], beta[5];
  CCTK_REAL alpha_slow[4], beta_slow[5];
  CCTK_REAL * restrict UpdateVar;
  CCTK_REAL * restrict OldVar;
  CCTK_REAL const * restrict RHSVar;
  CCTK_REAL * restrict ScratchVar;

  totalsize = 1;
  for (arraydim = 0; arraydim < cctk_dim; arraydim++)
  {
     totalsize *= cctk_ash[arraydim];
  }

  if (scratchspace_firstindex == -99)
  {
    scratchspace_firstindex = CCTK_FirstVarIndex("MOL::SCRATCHSPACE");
  }
  
  if (scratchspace_firstindex_slow == -99)
  {
    scratchspace_firstindex_slow = CCTK_FirstVarIndex("MOL::SCRATCHSPACESLOW");
  }

  switch (MoL_Intermediate_Steps - (*MoL_Intermediate_Step))
  {
    case 0:
      alpha[0] = 1.0/2.0; 
      alpha[1] = 0;
      alpha[2] = 0;
      alpha[3] = 0;
      alpha_slow[0] = 1.0/2.0; 
      alpha_slow[1] = 0;
      alpha_slow[2] = 0;
      alpha_slow[3] = 0;
      break;
    case 1:
      alpha[0] = 1.0/4.0;
      alpha[1] = 1.0/4.0;
      alpha[2] = 0;
      alpha[3] = 0;
      alpha_slow[0] = 1.0/2.0; 
      alpha_slow[1] = 0;
      alpha_slow[2] = 0;
      alpha_slow[3] = 0;
      break;
    case 2:
      alpha[0] = 1.0/4.0; 
      alpha[1] = 1.0/4.0;
      alpha[2] = 1.0/2.0;
      alpha[3] = 0.0;
      alpha_slow[0] = 1.0; 
      alpha_slow[1] = 0;
      alpha_slow[2] = 0;
      alpha_slow[3] = 0;
      break;
    case 3:
      alpha[0] = 1.0/4.0;
      alpha[1] = 1.0/4.0; 
      alpha[2] = 1.0/4.0; 
      alpha[3] = 1.0/4.0;
      alpha_slow[0] = 1.0; 
      alpha_slow[1] = 0;
      alpha_slow[2] = 0;
      alpha_slow[3] = 0;
  }

  beta[0] = 1.0/4.0; 
  beta[1] = 1.0/4.0; 
  beta[2] = 1.0/4.0; 
  beta[3] = 1.0/4.0;
  beta[4] = 0.0;

  beta_slow[0] = 1.0/2.0; 
  beta_slow[1] = 0.0; 
  beta_slow[2] = 0.0; 
  beta_slow[3] = 0.0;
  beta_slow[4] = 1.0/2.0;


  /* FIXME */


  /* Real GFs, the "fast" part */

  for (var = 0; var < MoLNumEvolvedVariables; var++)
  {
    
    UpdateVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                              EvolvedVariableIndex[var]);
    OldVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 1, 
                                              EvolvedVariableIndex[var]);
    RHSVar = (CCTK_REAL const *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                 RHSVariableIndex[var]);
/* #define MOLDEBUG 1 */
#ifdef MOLDEBUG
    printf("In multirate RK. Variable %d (%s). RHS %d (%s). beta %g.\n",
           EvolvedVariableIndex[var],
           CCTK_VarName(EvolvedVariableIndex[var]),
           RHSVariableIndex[var],
           CCTK_VarName(RHSVariableIndex[var]),
           beta[MoL_Intermediate_Steps - (*MoL_Intermediate_Step)]);
#endif

    ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                      scratchspace_firstindex
                                      + var
                                      + MoLNumEvolvedVariables * (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));

#pragma omp parallel for
    for (index = 0; index < totalsize; index++)
    {
      ScratchVar[index] = (*Original_Delta_Time) / cctkGH->cctk_timefac * RHSVar[index];
    
#ifdef MOLDEBUG
      if (CCTK_EQUALS(verbose,"extreme"))
      {
        printf("Variable: %d. Index: %d. dt: %f. beta %f. RHS: %f. q: %f.\n",
               var, index, (*Original_Delta_Time) / cctkGH->cctk_timefac, beta, RHSVar[index], 
               UpdateVar[index]);
      }
#endif
    }
    
    
#pragma omp parallel for
    for (index = 0; index < totalsize; index++)
    {
      UpdateVar[index] = OldVar[index];
    }

    if ((*MoL_Intermediate_Step)>1)
    {
      //printf("Step %d \n", MoL_Intermediate_Steps - (*MoL_Intermediate_Step));
      for (scratchstep = 0; scratchstep <= MoL_Intermediate_Steps - (*MoL_Intermediate_Step); scratchstep++)
      {
         
         //printf("Scratch Step %d, alpha %g \n", scratchstep, alpha[scratchstep]);
         
         ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                       scratchspace_firstindex
                                       + var
                                       + MoLNumEvolvedVariables * scratchstep);
      
#pragma omp parallel for
         for (index = 0; index < totalsize; index++)
         {
            UpdateVar[index] += alpha[scratchstep] * ScratchVar[index];
         }
      }
    }
    else
    {
      //printf("Final Step!\n");
    
      for (scratchstep = 0; scratchstep < MoL_Intermediate_Steps; scratchstep++)
      {
         
         //printf("Scratch Step %d, beta %g \n", scratchstep, beta[scratchstep]);
         
         ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                       scratchspace_firstindex
                                       + var
                                       + MoLNumEvolvedVariables * scratchstep);
      
#pragma omp parallel for
         for (index = 0; index < totalsize; index++)
         {
            UpdateVar[index] += beta[scratchstep] * ScratchVar[index];
         }
      }
    }

  }


  for (var = 0; var < MoLNumEvolvedVariablesSlow; var++)
  {
    
    UpdateVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                              EvolvedVariableIndexSlow[var]);
    OldVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 1, 
                                              EvolvedVariableIndexSlow[var]);
    RHSVar = (CCTK_REAL const *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                 RHSVariableIndexSlow[var]);
/* #define MOLDEBUG 1 */
#ifdef MOLDEBUG
    printf("In multirate RK. SLOW Variable %d (%s). RHS %d (%s). beta %g.\n",
           EvolvedVariableIndexSlow[var],
           CCTK_VarName(EvolvedVariableIndexSlow[var]),
           RHSVariableIndexSlow[var],
           CCTK_VarName(RHSVariableIndexSlow[var]),
           beta[MoL_Intermediate_Steps - (*MoL_Intermediate_Step)]);
#endif

    ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                      scratchspace_firstindex_slow
                                      + var
                                      + MoLNumEvolvedVariablesSlow * (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));

#pragma omp parallel for
    for (index = 0; index < totalsize; index++)
    {
      ScratchVar[index] = (*Original_Delta_Time) / cctkGH->cctk_timefac * RHSVar[index];
    
#ifdef MOLDEBUG
      if (CCTK_EQUALS(verbose,"extreme"))
      {
        printf("SLOW Variable: %d. Index: %d. dt: %f. beta %f. RHS: %f. q: %f.\n",
               var, index, (*Original_Delta_Time) / cctkGH->cctk_timefac, beta, RHSVar[index], 
               UpdateVar[index]);
      }
#endif
    }
    
    
#pragma omp parallel for
    for (index = 0; index < totalsize; index++)
    {
      UpdateVar[index] = OldVar[index];
    }

    if ((*MoL_Intermediate_Step)>1)
    {
      //printf("Step %d \n", MoL_Intermediate_Steps - (*MoL_Intermediate_Step));
      for (scratchstep = 0; scratchstep <= /*MoL_Intermediate_Steps - (*MoL_Intermediate_Step)*/ 0; scratchstep++)
      {
         
         //printf("Scratch Step %d, alpha %g \n", scratchstep, alpha_slow[scratchstep]);
         
         ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                       scratchspace_firstindex_slow
                                       + var
                                       + MoLNumEvolvedVariablesSlow * scratchstep);
      
#pragma omp parallel for
         for (index = 0; index < totalsize; index++)
         {
            UpdateVar[index] += alpha_slow[scratchstep] * ScratchVar[index];
         }
      }
    }
    else
    {
      //printf("Final Step!\n");
    
      for (scratchstep = 0; scratchstep < MoL_Intermediate_Steps; scratchstep+=4)
      {
         
         //printf("Scratch Step %d, beta %g \n", scratchstep, beta_slow[scratchstep]);
         
         ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                       scratchspace_firstindex_slow
                                       + var
                                       + MoLNumEvolvedVariablesSlow * scratchstep);
      
#pragma omp parallel for
         for (index = 0; index < totalsize; index++)
         {
            UpdateVar[index] += beta_slow[scratchstep] * ScratchVar[index];
         }
      }
    }

  }

  return;
}