aboutsummaryrefslogtreecommitdiff
path: root/src/GenericRK.c
blob: 52dc1f48c30231eae3e7a5034d3cecd96f6b6dcc (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
 /*@@
   @file      GenericRK.c
   @date      Sun May 26 03:47:15 2002
   @author    Ian Hawke
   @desc 
   This routine performs a generic Runge-Kutta type integration
   given the set of coefficients defined in the RKAlphaCoefficients
   and RKBetaCoefficients arrays. See the article by Shu referenced
   in the documentation for more details.
   @enddesc 
   @version   $Header$
 @@*/

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

#include "ExternalVariables.h"

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

CCTK_FILEVERSION(CactusMoL2_MoL2_GenericRK_c);

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

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

CCTK_INT AlphaIndex(CCTK_INT Step_Number, CCTK_INT Scratch_Level);

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

void MoL_GenericRKAdd(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_GenericRKAdd
   @date       Sun May 26 03:50:44 2002
   @author     Ian Hawke
   @desc 
   Performs a single step of a generic Runge-Kutta type time
   integration.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_GenericRKAdd(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INT index, var, scratchstep, alphaindex, scratchindex;
  CCTK_INT totalsize;
  CCTK_REAL alpha, beta;
  CCTK_REAL *UpdateVar;
  CCTK_REAL *RHSVar;
  CCTK_REAL *ScratchVar;

  totalsize = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
  
  beta = RKBetaCoefficients[MoL_Intermediate_Steps - 
                           (*MoL_Intermediate_Step)];

  for (var = 0; var < MoLNumEvolvedVariables; var++)
  {
    
    UpdateVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                              EvolvedVariableIndex[var]);
    RHSVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                           RHSVariableIndex[var]);

    /*
    printf("In generic RK. Variable %d (%s). RHS %d (%s).\n",
           EvolvedVariableIndex[var],
           CCTK_VarName(EvolvedVariableIndex[var]),
           RHSVariableIndex[var],
           CCTK_VarName(RHSVariableIndex[var]));
    */

    for (index = 0; index < totalsize; index++)
    {
      UpdateVar[index] = (*Original_Delta_Time) * beta * RHSVar[index];
      /*      UpdateVar[index] = CCTK_DELTA_TIME * beta * RHSVar[index];*/
      /*      UpdateVar[index] = CCTK_DELTA_TIME * RHSVar[index];*/
      /*
      printf("Variable: %d. Index: %d. dt: %f. beta %f. RHS: %f. q: %f.\n",
            var, index, (*Original_Delta_Time), beta, RHSVar[index], 
             UpdateVar[index]);
      */
      /*
      if (EvolvedVariableIndex[var]==CCTK_VarIndex("admbase::alp"))
      {
        
        printf("RK: Index: %d. dt: %f. beta %f. RHS: %f (%s). q: %f.\n",
            index, (*Original_Delta_Time), beta, RHSVar[index], 
               CCTK_VarName(RHSVariableIndex[var]),UpdateVar[index]);
      }
      */
    }
    
    for (scratchstep = 0; 
         scratchstep < MoL_Intermediate_Steps - (*MoL_Intermediate_Step) + 1;
         scratchstep++)
    {

      alphaindex = AlphaIndex(*MoL_Intermediate_Step, scratchstep);
      scratchindex = scratchstep - 1;

      alpha = RKAlphaCoefficients[alphaindex];

      if (scratchstep) 
      {
        ScratchVar = &ScratchSpace[(var * MoL_Num_Scratch_Levels + 
                                   scratchindex) * totalsize];
        /*
        printf("Reading from scratch space, initial address %ld index %d\n", 
               ScratchVar, (var * MoL_Num_Scratch_Levels + 
                            scratchindex) * totalsize);
        */
      }
      else
      {
        ScratchVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1, 
                                                  EvolvedVariableIndex[var]);
      }
      
      if ( (alpha > MoL_Tiny)||(alpha < -MoL_Tiny) )
      {
        for (index = 0; index < totalsize; index++)
        {
          UpdateVar[index] += alpha * ScratchVar[index];
          /*    
          printf("Variable: %d. Index: %d. step: %d. alpha: %f. Scratch: %f. q: %f.\n",
            var, index, (*MoL_Intermediate_Step), alpha, ScratchVar[index], UpdateVar[index]);
          */
        }
      }
      
    }
    
  }

  if (*MoL_Intermediate_Step > 1)
  {
    for (var = 0; var < MoLNumEvolvedVariables; var++)
    {
      UpdateVar = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                EvolvedVariableIndex[var]);
      ScratchVar = &ScratchSpace[(var * MoL_Num_Scratch_Levels + 
                                 MoL_Intermediate_Steps -
                                 (*MoL_Intermediate_Step)) * totalsize];
      /*
      printf("Writing to scratch space, initial address %ld, index %d \n", 
             ScratchVar, (var * MoL_Num_Scratch_Levels + 
                                 MoL_Intermediate_Steps -
                                 (*MoL_Intermediate_Step)) * totalsize);
      */  
      for (index = 0; index < totalsize; index++)
      {
        ScratchVar[index] = UpdateVar[index];
        /*
        printf("Variable: %d. Index: %d. step: %d. Scratch: %f.\n",
               var, index, (*MoL_Intermediate_Step), ScratchVar[index]);
        */
      }
    }
  }
        
  return;
}

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

CCTK_INT AlphaIndex(CCTK_INT Step_Number, CCTK_INT Scratch_Level)
{
  DECLARE_CCTK_PARAMETERS;
  
  return (MoL_Intermediate_Steps - Step_Number) * MoL_Intermediate_Steps + 
    Scratch_Level;
}