aboutsummaryrefslogtreecommitdiff
path: root/src/ICN.c
blob: b32bf15857312c2354646f6d044356319de3535a (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
 /*@@
   @file      ICN.c
   @date      Sun May 26 04:29:07 2002
   @author    Ian Hawke
   @desc 
   This implements the more efficient Iterative Crank Nicholson integrator.
   This follows the implementation of ICN in all AEI codes and is 
   equivalent to (but hopefully more efficient than) the generic ICN
   integrator also implemented by MoL.
   @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_ICN_c);

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

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

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

void MoL_ICNAdd(CCTK_ARGUMENTS);

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

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

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

 /*@@
   @routine    MoL_ICNAdd
   @date       Sun May 26 04:17:23 2002
   @author     Ian Hawke
   @desc 
   Performs Iterative Crank Nicholson time integration. The number of
   steps is arbitrary.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

void MoL_ICNAdd(CCTK_ARGUMENTS)
{

  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  cGroupDynamicData arraydata;
  CCTK_INT groupindex, ierr;
  CCTK_INT arraytotalsize, arraydim;

  CCTK_INT index, var;
  CCTK_INT totalsize;
  CCTK_REAL *OldVar;
  CCTK_REAL *UpdateVar;
  CCTK_REAL *RHSVar;
  
  /* FIXME */

#ifdef MOLDOESCOMPLEX

  CCTK_COMPLEX *OldComplexVar;
  CCTK_COMPLEX *UpdateComplexVar;
  CCTK_COMPLEX *RHSComplexVar;
  CCTK_COMPLEX Complex_Delta_Time = CCTK_Cmplx(CCTK_DELTA_TIME, 0);

#endif
  
#ifdef MOLDEBUG
  printf("Inside ICN.\nProcessor %d.\nStep %d.\nRefinement %d.\nTimestep %g.\nSpacestep %g.\nTime %g\n",
         CCTK_MyProc(cctkGH),
         MoL_Intermediate_Steps - *MoL_Intermediate_Step + 1,
         *cctk_levfac,
         cctk_delta_time,
         CCTK_DELTA_SPACE(0),
         cctk_time);
#endif  

  totalsize = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
#ifdef MOLDEBUG
  printf("MoL: the ICN routine says dt = %f.\n", CCTK_DELTA_TIME);
#endif
  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] = OldVar[index] + CCTK_DELTA_TIME * RHSVar[index];
    }
  }
  
  for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
  {
    OldVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
                                             EvolvedArrayVariableIndex[var]);
    UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
                                             EvolvedArrayVariableIndex[var]);
    RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, 
                                          RHSArrayVariableIndex[var]);
    
    groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]);
    ierr = CCTK_GroupDynamicData(cctkGH, groupindex,
                                 &arraydata);
    if (ierr)
    {
      CCTK_VWarn(0, __LINE__, __FILE__, "MoL", 
                 "The driver does not return group information for group '%s'.", 
                 CCTK_GroupName(groupindex));
    }
    arraytotalsize = 1;
    for (arraydim = 0; arraydim < arraydata.dim; arraydim++)
    {
      arraytotalsize *= arraydata.lsh[arraydim];
    }

    for (index = 0; index < arraytotalsize; index++)
    {
      UpdateVar[index] = OldVar[index] + CCTK_DELTA_TIME * RHSVar[index];
    }
  }

  /* FIXME */

#ifdef MOLDOESCOMPLEX
  
  for (var = 0; var < MoLNumEvolvedComplexVariables; var++)
  {
    OldComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 1,
                                                    EvolvedComplexVariableIndex[var]);
    UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
                                                       EvolvedComplexVariableIndex[var]);
    RHSComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0, 
                                                    RHSComplexVariableIndex[var]);
    
    for (index = 0; index < totalsize; index++)
    {
      UpdateComplexVar[index] = CCTK_CmplxAdd(OldComplexVar[index],
                                              CCTK_CmplxMul(Complex_Delta_Time,
                                                            RHSComplexVar[index]));
    }
  }

#endif
  
  return;

}

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