aboutsummaryrefslogtreecommitdiff
path: root/src/RK3.c
diff options
context:
space:
mode:
Diffstat (limited to 'src/RK3.c')
-rw-r--r--src/RK3.c250
1 files changed, 73 insertions, 177 deletions
diff --git a/src/RK3.c b/src/RK3.c
index f3b5e9c..38a55ab 100644
--- a/src/RK3.c
+++ b/src/RK3.c
@@ -15,6 +15,7 @@
#include "cctk_Parameters.h"
#include "ExternalVariables.h"
+#include "Operators.h"
static const char *rcsid = "$Header$";
@@ -67,15 +68,10 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
- cGroupDynamicData arraydata;
- CCTK_INT groupindex, ierr;
- CCTK_INT arraytotalsize, arraydim;
+ CCTK_INT arraydim;
- CCTK_INT index, var;
+ CCTK_INT var;
CCTK_INT totalsize;
- CCTK_REAL const * restrict OldVar;
- CCTK_REAL * restrict UpdateVar;
- CCTK_REAL const * restrict RHSVar;
/* FIXME */
@@ -106,7 +102,7 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
totalsize = 1;
for (arraydim = 0; arraydim < cctk_dim; arraydim++)
{
- totalsize *= cctk_lsh[arraydim];
+ totalsize *= cctk_ash[arraydim];
}
switch (*MoL_Intermediate_Step)
@@ -116,45 +112,24 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
{
for (var = 0; var < MoLNumEvolvedVariables; var++)
{
- UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedVariableIndex[var]);
- RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 1;
+ CCTK_INT const srcs[] = {RHSVariableIndex[var]};
+ CCTK_INT const tls[] = {0};
+ CCTK_REAL const facts[] = {CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedVariableIndex[var], 1.0,
+ srcs, tls, facts, nsrcs);
}
for (var = 0; var < MoLNumEvolvedArrayVariables; 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__, CCTK_THORNSTRING,
- "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];
- }
-
-#pragma omp parallel for
- for (index = 0; index < arraytotalsize; index++)
- {
- UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 1;
+ CCTK_INT const srcs[] = {RHSArrayVariableIndex[var]};
+ CCTK_INT const tls[] = {0};
+ CCTK_REAL const facts[] = {CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedArrayVariableIndex[var], 1.0,
+ srcs, tls, facts, nsrcs);
}
/* FIXME */
@@ -163,18 +138,13 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
for (var = 0; var < MoLNumEvolvedComplexVariables; var++)
{
- UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedComplexVariableIndex[var]);
- RHSComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSComplexVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateComplexVar[index] = CCTK_CmplxAdd(UpdateComplexVar[index],
- CCTK_CmplxMul(Complex_Delta_Time,
- RHSComplexVar[index]));
- }
+ int const nsrcs = 1;
+ CCTK_INT const srcs[] = {RHSComplexVariableIndex[var]};
+ CCTK_INT const tls[] = {0};
+ CCTK_REAL const facts[] = {CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedComplexVariableIndex[var], 1.0,
+ srcs, tls, facts, nsrcs);
}
#endif
@@ -186,52 +156,26 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
{
for (var = 0; var < MoLNumEvolvedVariables; var++)
{
- OldVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedVariableIndex[var]);
- UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedVariableIndex[var]);
- RHSVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateVar[index] = 0.25 * (3*OldVar[index] +
- UpdateVar[index]) +
- CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedVariableIndex[var], RHSVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {0.75, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedVariableIndex[var], 0.25,
+ srcs, tls, facts, nsrcs);
}
for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
{
- OldVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedArrayVariableIndex[var]);
- UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedArrayVariableIndex[var]);
- RHSVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSArrayVariableIndex[var]);
- groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]);
- ierr = CCTK_GroupDynamicData(cctkGH, groupindex,
- &arraydata);
- if (ierr)
- {
- CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
- "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];
- }
-
-#pragma omp parallel for
- for (index = 0; index < arraytotalsize; index++)
- {
- UpdateVar[index] = 0.25*(3*OldVar[index] + UpdateVar[index]) +
- CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedArrayVariableIndex[var], RHSArrayVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {0.75, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedArrayVariableIndex[var], 0.25,
+ srcs, tls, facts, nsrcs);
}
/* FIXME */
@@ -240,24 +184,14 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
for (var = 0; var < MoLNumEvolvedComplexVariables; var++)
{
- OldComplexVar = (CCTK_COMPLEX const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedComplexVariableIndex[var]);
- UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedComplexVariableIndex[var]);
- RHSComplexVar = (CCTK_COMPLEX const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSComplexVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateComplexVar[index] =
- CCTK_CmplxAdd(CCTK_CmplxAdd(CCTK_CmplxMul(Complex_ThreeQuarter,
- OldComplexVar[index]),
- CCTK_CmplxMul(ComplexQuarter,
- UpdateComplexVar[index])),
- CCTK_CmplxMul(CCTK_CmplxMul(Complex_Delta_Time,
- RHSComplexVar[index])));
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedComplexVariableIndex[var], RHSComplexVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {0.75, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedComplexVariableIndex[var], 0.25,
+ srcs, tls, facts, nsrcs);
}
#endif
@@ -266,56 +200,28 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
}
case 1:
{
- const CCTK_REAL one_third = 1.0 / 3.0;
-
for (var = 0; var < MoLNumEvolvedVariables; var++)
{
- OldVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedVariableIndex[var]);
- UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedVariableIndex[var]);
- RHSVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateVar[index] = (OldVar[index] + 2*UpdateVar[index]) * one_third
- + CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedVariableIndex[var], RHSVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {1.0/3.0, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedVariableIndex[var], 2.0/3.0,
+ srcs, tls, facts, nsrcs);
}
for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
{
- OldVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedArrayVariableIndex[var]);
- UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedArrayVariableIndex[var]);
- RHSVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSArrayVariableIndex[var]);
-
- groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]);
- ierr = CCTK_GroupDynamicData(cctkGH, groupindex,
- &arraydata);
- if (ierr)
- {
- CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
- "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];
- }
-
-#pragma omp parallel for
- for (index = 0; index < arraytotalsize; index++)
- {
- UpdateVar[index] = (OldVar[index] + 2*UpdateVar[index]) * one_third
- + CCTK_DELTA_TIME * RHSVar[index];
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedArrayVariableIndex[var], RHSArrayVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {1.0/3.0, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedArrayVariableIndex[var], 2.0/3.0,
+ srcs, tls, facts, nsrcs);
}
/* FIXME */
@@ -324,24 +230,14 @@ void MoL_RK3Add(CCTK_ARGUMENTS)
for (var = 0; var < MoLNumEvolvedComplexVariables; var++)
{
- OldComplexVar = (CCTK_COMPLEX const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedComplexVariableIndex[var]);
- UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedComplexVariableIndex[var]);
- RHSComplexVar = (CCTK_COMPLEX const*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSComplexVariableIndex[var]);
-
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- UpdateComplexVar[index] =
- CCTK_CmplxAdd(CCTK_CmplxAdd(CCTK_CmplxMul(Complex_Third,
- OldComplexVar[index]),
- CCTK_CmplxMul(ComplexTwoThird,
- UpdateComplexVar[index])),
- CCTK_CmplxMul(Complex_Delta_Time,
- RHSComplexVar[index]));
- }
+ int const nsrcs = 2;
+ CCTK_INT const srcs[] =
+ {EvolvedComplexVariableIndex[var], RHSComplexVariableIndex[var]};
+ CCTK_INT const tls[] = {1, 0};
+ CCTK_REAL const facts[] = {1.0/3.0, CCTK_DELTA_TIME};
+ MoL_LinearCombination(cctkGH,
+ EvolvedComplexVariableIndex[var], 2.0/3.0,
+ srcs, tls, facts, nsrcs);
}
#endif