aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2005-12-11 18:53:12 +0000
committerschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2005-12-11 18:53:12 +0000
commit0863f0decdf6024f548f6178ce6221bc1a6fa722 (patch)
treed7a06b78aa6ce4adbf3def63596d74029868bb46 /src
parentd2581020e698490172cee2cdea2e2deb43bf9251 (diff)
Remember the variable index of the scratch variables.
This saves calling CCTK_FirstVarIndex during each MoL iteration. This function is expensive, as it performs case insensitive string comparisons, and it showed up high on a profile of some of the Mexico tests, which integrate small domains over long times. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@105 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src')
-rw-r--r--src/GenericRK.c11
1 files changed, 8 insertions, 3 deletions
diff --git a/src/GenericRK.c b/src/GenericRK.c
index 196c83e..60648f1 100644
--- a/src/GenericRK.c
+++ b/src/GenericRK.c
@@ -101,7 +101,12 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
{
totalsize *= cctk_lsh[arraydim];
}
-
+
+ if (scratchspace_firstindex == -99)
+ {
+ scratchspace_firstindex = CCTK_FirstVarIndex("MOL::SCRATCHSPACE");
+ }
+
beta = RKBetaCoefficients[MoL_Intermediate_Steps -
(*MoL_Intermediate_Step)];
@@ -176,7 +181,7 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
scratchindex) * totalsize];
*/
ScratchVar = CCTK_VarDataPtrI(cctkGH, 0,
- CCTK_FirstVarIndex("MOL::SCRATCHSPACE")
+ scratchspace_firstindex
+ var
+ MoL_Num_Evolved_Vars * scratchindex);
#ifdef MOLDEBUG
@@ -231,7 +236,7 @@ void MoL_GenericRKAdd(CCTK_ARGUMENTS)
(*MoL_Intermediate_Step)) * totalsize];
*/
ScratchVar = CCTK_VarDataPtrI(cctkGH, 0,
- CCTK_FirstVarIndex("MOL::SCRATCHSPACE")
+ scratchspace_firstindex
+ var
+ MoL_Num_Evolved_Vars * (MoL_Intermediate_Steps - (*MoL_Intermediate_Step)));
#ifdef MOLDEBUG