aboutsummaryrefslogtreecommitdiff
path: root/src/RK2.c
blob: 72565c50467386e83a11fe84e6b490ffebb6812d (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
 /*@@
   @file      RK2.c
   @date      Sun May 26 04:13:45 2002
   @author    Ian Hawke
   @desc 
   A specialized second order Runge-Kutta time integrator. This is
   the integrator that Shu refers to as the optimal TVD second 
   order method (see reference in documentation). It is equivalent
   to Heun's predictor-corrector method, or the MacCormack method.
   @enddesc 
   @version   $Header$
 @@*/

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

#include "ExternalVariables.h"

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

CCTK_FILEVERSION(AlphaThorns_MoL_RK2_c);

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

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

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

void MoL_RK2Add(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_RK2Add
   @date       Sun May 26 04:17:23 2002
   @author     Ian Hawke
   @desc 
   Performs second order Runge-Kutta time integration.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_RK2Add(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INT index, var;
  CCTK_INT totalsize;
  CCTK_REAL *OldVar;
  CCTK_REAL *UpdateVar;
  CCTK_REAL *RHSVar;

  /*
  printf("Inside RK2.\nStep %d.\nRefinement %d.\nTimestep %g.\nSpacestep %g.\nTime %g\n",
         MoL_Intermediate_Steps - *MoL_Intermediate_Step + 1,
         *cctk_levfac,
         cctk_delta_time,
         CCTK_DELTA_SPACE(0),
         cctk_time);
  */

  totalsize = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];

  switch (*MoL_Intermediate_Step)
  {
  
    case 2:
      {
        for (var = 0; var < MoLNumEvolvedVariables; var++)
        {
          UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
                                                   EvolvedVariableIndex[var]);
          RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, 
                                                RHSVariableIndex[var]);
          
          for (index = 0; index < totalsize; index++)
          {
            UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index];
          }
        }
        break;
      }
    case 1:
      {
        for (var = 0; var < MoLNumEvolvedVariables; var++)
        {
          OldVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
                                                EvolvedVariableIndex[var]);
          UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
                                                   EvolvedVariableIndex[var]);
          RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, 
                                                RHSVariableIndex[var]);
          
          for (index = 0; index < totalsize; index++)
          {
            UpdateVar[index] = 0.5 * (OldVar[index] + UpdateVar[index]) + 
                                      CCTK_DELTA_TIME * RHSVar[index];
          }
        }
        break;
      }
    default:
      {
        CCTK_WARN(0, "RK2 expects MoL_Intermediate_Step to be in [1,2]. This should be caught at ParamCheck - bug Ian!");
        break;
      }
      
  }

  return;

}

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