aboutsummaryrefslogtreecommitdiff
path: root/src/RK4.c
blob: 4e72389e986f380ea36a6f75ba305e390a4e02d0 (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
 /*@@
   @file      RK4.c
   @date      Fri July 14, 2006
   @author    Yosef Zlochower
   @desc 
   A routine to perform RK4 evolution. Mostly copied from
   genericRK.c
   @enddesc 
   @version   $Header$
 @@*/

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

#include "ExternalVariables.h"
#include "Operators.h"

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

CCTK_FILEVERSION(CactusBase_MoL_RK4_c);

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

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

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

void MoL_RK4Add(CCTK_ARGUMENTS);

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

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

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

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

@@*/

void MoL_RK4Add(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
    
  cGroupDynamicData arraydata;
  CCTK_INT groupindex, ierr;
  CCTK_INT arraytotalsize, 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;
  CCTK_INT index, var;
  CCTK_INT totalsize;
  CCTK_REAL alpha, beta;
  CCTK_REAL * restrict UpdateVar;
  CCTK_REAL * restrict OldVar;
  CCTK_REAL const * restrict RHSVar;
  CCTK_REAL * restrict ScratchVar;

  CCTK_INT arrayscratchlocation;

  /* Keep a running total of alpha as we perform the substeps, so that
     we know the "real" alpha (including round-off errors) when we
     calculate the final result. */
  CCTK_REAL const time_rhs = 1.0;
  CCTK_REAL const old_time = 0.0;
  static CCTK_REAL time, scratch_time;
  
  totalsize = 1;
  for (arraydim = 0; arraydim < cctk_dim; arraydim++)
  {
    totalsize *= cctk_ash[arraydim];
  }

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

  switch (MoL_Intermediate_Steps - (*MoL_Intermediate_Step))
  {
    case 0:
      alpha = 1.0 / 3.0;
      beta  = 1.0 / 2.0;
      break;
    case 1:
      alpha = 2.0 / 3.0;
      beta  = 1.0 / 2.0;
      break;
    case 2:
      alpha = 1.0 / 3.0;
      beta  = 1.0;
      break;
    case 3:
      alpha = 1.0;
      beta  = 1.0 / 6.0;
      break;
  }

  if (MoL_Intermediate_Steps == (*MoL_Intermediate_Step))
  {
    time = 0.0;
  }
  
  /* FIXME */

#ifdef MOLDOESCOMPLEX
  
  Complex_Delta_Time = CCTK_Cmplx((*Original_Delta_Time) / cctkGH->cctk_timefac, 0);
  Complex_beta = CCTK_Cmplx(beta, 0);

#endif

  /* Real GFs */

  for (var = 0; var < MoLNumEvolvedVariables; var++)
  {
    {
      int const nsrcs = 2;
      CCTK_INT const srcs[] =
        {EvolvedVariableIndex[var], RHSVariableIndex[var]};
      CCTK_INT const tls[] = {1, 0};
      CCTK_REAL const facts[] =
        {1.0, (*Original_Delta_Time) / cctkGH->cctk_timefac * beta};
      MoL_LinearCombination(cctkGH,
                            EvolvedVariableIndex[var], 0.0,
                            srcs, tls, facts, nsrcs);
      time = facts[0] * old_time + facts[1] * time_rhs;
    }

    /* scratch storage */
    if ((*MoL_Intermediate_Step) == MoL_Intermediate_Steps)
    {
      MoL_LinearCombination(cctkGH,
                            scratchspace_firstindex + var, 0.0,
                            NULL, NULL, NULL, 0);
      scratch_time = 0.0;
    }

    if ((*MoL_Intermediate_Step)>1)
    {
      int const nsrcs = 1;
      CCTK_INT const srcs[] = {EvolvedVariableIndex[var]};
      CCTK_INT const tls[] = {0};
      CCTK_REAL const facts[] = {alpha};
      MoL_LinearCombination(cctkGH,
                            scratchspace_firstindex + var, 1.0,
                            srcs, tls, facts, nsrcs);
      scratch_time += facts[0] * time;
    }
    else
    {
      int const nsrcs = 2;
      CCTK_INT const srcs[] =
        {scratchspace_firstindex + var, EvolvedVariableIndex[var]};
      CCTK_INT const tls[] = {0, 1};
      CCTK_REAL const facts[] = {1.0, -4.0/3.0};
      MoL_LinearCombination(cctkGH,
                            EvolvedVariableIndex[var], 1.0,
                            srcs, tls, facts, nsrcs);
      time += facts[0] * scratch_time + facts[1] * old_time;
    }
  }

  /* Complex GFs */

  /* FIXME */

#ifdef MOLDOESCOMPLEX
CCTK_WARN(0, "not done");
#endif

  /* Real arrays */

  arrayscratchlocation = 0;

  for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
  {
    
    UpdateVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                              EvolvedArrayVariableIndex[var]);
    OldVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 1, 
                                              EvolvedArrayVariableIndex[var]);
    RHSVar = (CCTK_REAL const *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                 RHSArrayVariableIndex[var]);
    
    groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]);
    ierr = CCTK_GroupDynamicData(cctkGH, groupindex,
                                 &arraydata);
    if (ierr)
    {
      CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING, 
                 "The driver does not return group information "
                 "for group '%s'.", 
                 CCTK_GroupName(groupindex));
    }
    arraytotalsize = 1;
    for (arraydim = 0; arraydim < arraydata.dim; arraydim++)
    {
      arraytotalsize *= arraydata.ash[arraydim];
    }

    ScratchVar = &ArrayScratchSpace[arrayscratchlocation];

#pragma omp parallel for
    for (index = 0; index < arraytotalsize; index++)
    {
      UpdateVar[index] = OldVar[index] +
        (*Original_Delta_Time) / cctkGH->cctk_timefac * beta * RHSVar[index];
    }

    if ((*MoL_Intermediate_Step) == MoL_Intermediate_Steps)
    {
#pragma omp parallel for
      for (index = 0; index < arraytotalsize; index++)
      {
        ScratchVar[index] = 0;
      }
    }
    
    if ((*MoL_Intermediate_Step)>1)
    {
#pragma omp parallel for
      for (index = 0; index < arraytotalsize; index++)
      {
        ScratchVar[index] += alpha * UpdateVar[index];
      }
    }
    else
    {
#pragma omp parallel for
      for (index = 0; index < arraytotalsize; index++)
      {
        UpdateVar[index] += ScratchVar[index] - 4.0 / 3.0 * OldVar[index];
      }
    }
    arrayscratchlocation += arraytotalsize;
  }

  return;
}