aboutsummaryrefslogtreecommitdiff
path: root/src/SandR.c
blob: 577a8fe6e581909e91d2e546e351e67b546ff45b (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
 /*@@
   @file      SandR.c
   @date      Sun May 26 03:35:58 2002
   @author    Ian Hawke
   @desc 
   Restores the Save and Restore variables to their original positions.
   @enddesc 
   @version   $Header$
 @@*/

#include <string.h>

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

#include "ExternalVariables.h"

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

CCTK_FILEVERSION(CactusBase_MoL_SandR_c);

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

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

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

void MoL_RestoreSandR(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_RestoreSandR
   @date       Sun May 26 03:39:02 2002
   @author     Ian Hawke
   @desc 
   Save and Restore variables are those that the physics thorn may 
   need to know to calculate the RHS, but which may be evolved by
   something other than MoL. In order to get the timelevels correct,
   the previous timelevel is copied to the current before the MoL step.
   As we do not know whether the variable will be evolved before or
   after MoL we must save the data that was in the current timelevel, 
   and then restore it at the end of the MoL timestep. This routine
   restores the variables.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_RestoreSandR(CCTK_ARGUMENTS) 
{
  
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INT /*  index, */ var;
  CCTK_INT totalsize, arraydim;
  CCTK_REAL *SandRDataArray;
  CCTK_REAL *ScratchVar;

  /* FIXME */

#ifdef MOLDOESCOMPLEX

  CCTK_COMPLEX *SandRComplexDataArray;
  CCTK_COMPLEX *ComplexScratchVar;

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

  for (var = 0; var < MoLNumSandRVariables; var++) 
  {
    
    SandRDataArray = (CCTK_REAL *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                   SandRVariableIndex[var]);

    ScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                  CCTK_FirstVarIndex("MOL::SANDRSCRATCHSPACE")
                                  + var);

#ifdef MOLDEBUG    
    printf("Restore:Variable %s, first entry %g, scratch %g\n",
           CCTK_VarName(SandRVariableIndex[var]), SandRDataArray[0], 
           ScratchVar[0]);
#endif

    memcpy(SandRDataArray, ScratchVar, 
           totalsize * sizeof(CCTK_REAL));
  }

  /* FIXME */

#ifdef MOLDOESCOMPLEX

  for (var = 0; var < MoLNumSandRComplexVariables; var++) 
  {
    
    SandRComplexDataArray = (CCTK_COMPLEX *)CCTK_VarDataPtrI(cctkGH, 0, 
                                                             SandRComplexVariableIndex[var]);
    
    ComplexScratchVar = CCTK_VarDataPtrI(cctkGH, 0, 
                                         CCTK_FirstVarIndex("MOL::COMPLEXSANDRSCRATCHSPACE")
                                         + var);

    memcpy(SandRComplexDataArray, ComplexScratchVar, 
           totalsize * sizeof(CCTK_COMPLEX));
  }

#endif
  
  return;
  
}

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