aboutsummaryrefslogtreecommitdiff
path: root/src/ICN.c
blob: 589e4d966c78c85984cd335b4478c9adba864ba6 (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
 /*@@
   @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;
  
  CCTK_INT index, var;
  CCTK_INT totalsize;
  CCTK_REAL *OldVar;
  CCTK_REAL *UpdateVar;
  CCTK_REAL *RHSVar;
  
  /*
  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);
  */

  totalsize = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
  /*
  printf("MoL: the ICN routine says dt = %f.\n", CCTK_DELTA_TIME);
  */
  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];
      /*
      if (CCTK_VarIndex("wavetoymol::phi") == EvolvedVariableIndex[var]){
        CCTK_VWarn(1,__LINE__,__FILE__,"MoL",
                   "ICN: proc %d index %d, dt %g, oldvar %g, rhs %g, new %g\n",
                   CCTK_MyProc(cctkGH), index,CCTK_DELTA_TIME, 
                   OldVar[index], RHSVar[index], UpdateVar[index]);
                   }*/
    }
  }
  
  return;

}

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