aboutsummaryrefslogtreecommitdiff
path: root/src/RK4-RK2.c
blob: 6e5b85e8065bb6866b476a9d06c5c85ce6e52174 (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
 /*@@
   @file      RK4-RK2.c
   @date      2012-03-25
   @author    Christian Reisswig
   @desc 
   A routine to perform homegrown RK4RK2 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 */

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


void MoL_RK4_RK2_Add(CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
    
#ifdef MOLDOESCOMPLEX
  CCTK_WARN(CCTK_WARN_ABORT, "not implemented");
#endif
  
  static int scratchspace_firstindex      = -1;
  static int scratchspace_firstindex_slow = -1;
  if (scratchspace_firstindex < 0) {
    scratchspace_firstindex      = CCTK_FirstVarIndex("MOL::SCRATCHSPACE");
    scratchspace_firstindex_slow = CCTK_FirstVarIndex("MOL::SCRATCHSPACESLOW");
  }
  
  int const step = MoL_Intermediate_Steps - *MoL_Intermediate_Step;
  
  int totalsize = 1;
  for (int d=0; d<cctk_dim; ++d) totalsize *= cctk_lsh[d];
  
  CCTK_REAL const dt = *Original_Delta_Time / cctkGH->cctk_timefac;
  
  
  
  int const allvar1 = MoLNumEvolvedVariables + MoLNumEvolvedVariablesSlow;
#pragma omp parallel for schedule(dynamic)
  for (int var1=0; var1<allvar1; ++var1) {
    
    if (var1 < MoLNumEvolvedVariables) {
      /* a fast variable */
      int const var = var1;
      
      CCTK_REAL *restrict const UpdateVar =
        CCTK_VarDataPtrI(cctkGH, 0, EvolvedVariableIndex[var]);
      CCTK_REAL const *restrict const OldVar =
        CCTK_VarDataPtrI(cctkGH, 1, EvolvedVariableIndex[var]);
      CCTK_REAL const *restrict const RHSVar =
        CCTK_VarDataPtrI(cctkGH, 0, RHSVariableIndex[var]);
      
#define SCRATCHINDEX(step)                                              \
      (scratchspace_firstindex + var + MoLNumEvolvedVariables * (step))
      CCTK_REAL *restrict const ScratchVar =
        CCTK_VarDataPtrI(cctkGH, 0, SCRATCHINDEX(0));
      
      switch (step) {
        
      case 0:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          ScratchVar[i] = OldVar[i] + 1.0/3.0 * scaled_rhs;
          UpdateVar[i] = OldVar[i] + 0.5 * scaled_rhs;
        }
        break;
        
      case 1:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          ScratchVar[i] += 1.0/6.0 * scaled_rhs;
          UpdateVar[i] = OldVar[i] + 0.5 * scaled_rhs;
        }
        break;
        
      case 2:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          ScratchVar[i] += 1.0/6.0 * scaled_rhs;
          UpdateVar[i] = OldVar[i] + scaled_rhs;
        }
        break;
        
      case 3:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          /* ScratchVar contains OldVar */
          UpdateVar[i] = ScratchVar[i] + 1.0/3.0 * scaled_rhs;
        }
        break;
        
      default:
        assert(0);
      }
#undef SCRATCHINDEX
      
    } else {
      /* a slow variable */
      int const var = var1 - MoLNumEvolvedVariables;
      
      CCTK_REAL *restrict const UpdateVar =
        CCTK_VarDataPtrI(cctkGH, 0, EvolvedVariableIndexSlow[var]);
      CCTK_REAL const *restrict const OldVar =
        CCTK_VarDataPtrI(cctkGH, 1, EvolvedVariableIndexSlow[var]);
      CCTK_REAL const *restrict const RHSVar =
        CCTK_VarDataPtrI(cctkGH, 0, RHSVariableIndexSlow[var]);
      
#define SCRATCHINDEX(step)                                              \
      (scratchspace_firstindex_slow + var + MoLNumEvolvedVariablesSlow * (step))
      CCTK_REAL *restrict const ScratchVar =
        CCTK_VarDataPtrI(cctkGH, 0, SCRATCHINDEX(0));
      
      switch (step) {
        
      case 0:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          CCTK_REAL const scratchval = OldVar[i] + scaled_rhs;
          ScratchVar[i] = scratchval;
          UpdateVar[i] = scratchval;
        }
        break;
        
      case 1:
      case 2:
        for (int i=0; i<totalsize; ++i) {
          /* This is the same value as for the previous MoL step.
             However, MoL_PostStep may have modified it (e.g. enforced
             a constraint), so we need to recreate the original value
             here for consistency. */
          /* ScratchVar contains OldVar */
          UpdateVar[i] = ScratchVar[i];
        }
        break;
        
      case 3:
        for (int i=0; i<totalsize; ++i) {
          CCTK_REAL const scaled_rhs = dt * RHSVar[i];
          /* ScratchVar contains OldVar */
          UpdateVar[i] =
            0.5 * OldVar[i] + 0.5 * ScratchVar[i] + 0.5 * scaled_rhs;
        }
        break;
        
      default:
        assert(0);
      }
#undef SCRATCHINDEX
      
    } /* if fast or slow */
  }   /* for var */
  
}