aboutsummaryrefslogtreecommitdiff
path: root/src/GenericRK.c
diff options
context:
space:
mode:
authorrhaas <rhaas@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2014-03-12 05:13:28 +0000
committerrhaas <rhaas@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2014-03-12 05:13:28 +0000
commit9558b41ae2b68ba801fdf7f1cdd9b9b542de9f0f (patch)
tree5d3414e2d3503294e5db48c180a75c8f17543301 /src/GenericRK.c
parentd7bb127b4e5101e8f8d8c644d9ae7229effc3ea6 (diff)
automatically count number of evolved variables
this patch obviates the MoL_Num_Evolved_Variables etc accumulator parameters. Instead MoL counts the number of registered variables and allocates memory for scratch space etc by adjusting the number of time levels. No change of user code is required, however users no longer have to adjust MoL_Num_Evolved_Variables in their parfiles. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@215 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/GenericRK.c')
-rw-r--r--src/GenericRK.c39
1 files changed, 17 insertions, 22 deletions
diff --git a/src/GenericRK.c b/src/GenericRK.c
index 6c88ca6..f723e19 100644
--- a/src/GenericRK.c
+++ b/src/GenericRK.c
@@ -87,6 +87,7 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
#endif
static CCTK_INT scratchspace_firstindex = -99;
+ static CCTK_INT complexscratchspace_firstindex = -99;
CCTK_INT index, var, scratchstep, alphaindex, scratchindex;
CCTK_INT totalsize;
CCTK_REAL alpha, beta;
@@ -105,6 +106,7 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
if (scratchspace_firstindex == -99)
{
scratchspace_firstindex = CCTK_FirstVarIndex("MOL::SCRATCHSPACE");
+ complexscratchspace_firstindex = CCTK_FirstVarIndex("MOL::COMPLEXSCRATCHSPACE");
}
beta = RKBetaCoefficients[MoL_Intermediate_Steps -
@@ -181,16 +183,14 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
ScratchVar = &ScratchSpace[(var * MoL_Num_Scratch_Levels +
scratchindex) * totalsize];
*/
- ScratchVar = CCTK_VarDataPtrI(cctkGH, 0,
+ ScratchVar = CCTK_VarDataPtrI(cctkGH, var,
scratchspace_firstindex
- + var
- + MoL_Num_Evolved_Vars * scratchindex);
+ + scratchindex);
#ifdef MOLDEBUG
if (CCTK_EQUALS(verbose,"extreme"))
{
- printf("Reading from scratch space, initial address %ld index %d\n",
- ScratchVar, (var * MoL_Num_Scratch_Levels +
- scratchindex) * totalsize);
+ printf("Reading from scratch space var %d, initial index %d\n",
+ var, scratchindex);
}
#endif
}
@@ -237,16 +237,12 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
MoL_Intermediate_Steps -
(*MoL_Intermediate_Step)) * totalsize];
*/
- ScratchVar = CCTK_VarDataPtrI(cctkGH, 0,
+ ScratchVar = CCTK_VarDataPtrI(cctkGH, var,
scratchspace_firstindex
- + var
- + MoL_Num_Evolved_Vars * (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));
+ + (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));
#ifdef MOLDEBUG
- printf("Writing to scratch space, initial address %ld, "
- "index %d, totalsize %d \n",
- ScratchVar, (var * MoL_Num_Scratch_Levels +
- MoL_Intermediate_Steps -
- (*MoL_Intermediate_Step)) * totalsize, totalsize);
+ printf("Writing to scratch space, var %d, index %d\n",
+ var, MoL_Intermediate_Steps - (*MoL_Intermediate_Step));
#endif
#pragma omp parallel for
for (index = 0; index < totalsize; index++)
@@ -300,14 +296,14 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
if (scratchstep)
{
- ScratchComplexVar = &ComplexScratchSpace[(var *
- MoL_Num_Scratch_Levels +
- scratchindex) * totalsize];
+ ScratchComplexVar = CCTK_VarDataPtrI(cctkGH, var,
+ complexscratchspace_firstindex
+ + scratchindex);
}
else
{
ScratchComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedComplexVariableIndex[var]);
+ EvolvedComplexVariableIndex[var]);
}
if ( (alpha > MoL_Tiny)||(alpha < -MoL_Tiny) )
@@ -332,10 +328,9 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
{
UpdateComplexVar = (CCTK_COMPLEX *)CCTK_VarDataPtrI(cctkGH, 0,
EvolvedComplexVariableIndex[var]);
- ScratchComplexVar = &ComplexScratchSpace[(var * MoL_Num_Scratch_Levels +
- MoL_Intermediate_Steps -
- (*MoL_Intermediate_Step)) *
- totalsize];
+ ScratchComplexVar = CCTK_VarDataPtrI(cctkGH, var,
+ complexscratchspace_firstindex
+ + (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));
#pragma omp parallel for
for (index = 0; index < totalsize; index++)
{