From d5bc7cf9626df7800644836c81712ace1e92ed0f Mon Sep 17 00:00:00 2001 From: eschnett Date: Tue, 18 Feb 2014 18:52:34 +0000 Subject: Support multiple refinement and time levels in LinearCombination API git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@210 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b --- interface.ccl | 7 +++++++ src/Euler.c | 10 ++++++++-- src/InitialCopy.c | 16 +++++++++++++-- src/Operators.c | 60 ++++++++++++++++++++++++++++++++++++------------------- src/Operators.h | 2 ++ src/RK2.c | 18 +++++++++++------ src/RK3.c | 24 +++++++++++++--------- src/RK4.c | 14 +++++++++---- src/RK87.c | 14 +++++++++---- 9 files changed, 118 insertions(+), 47 deletions(-) diff --git a/interface.ccl b/interface.ccl index 51d9647..f12606b 100644 --- a/interface.ccl +++ b/interface.ccl @@ -211,11 +211,17 @@ PROVIDES FUNCTION MoLNumIntegratorSubsteps WITH MoL_NumIntegratorSubsteps \ ### to override low-level grid variable operations. ### ############################################################## +CCTK_INT FUNCTION GetRefinementLevel \ + (CCTK_POINTER_TO_CONST IN cctkGH) +USES FUNCTION GetRefinementLevel + # Computes: # var = scale * var + \sum_i^nsrcs facts[i] * scrcs[i][tls[i]] CCTK_INT FUNCTION LinearCombination \ (CCTK_POINTER_TO_CONST IN cctkGH, \ CCTK_INT IN var, \ + CCTK_INT IN rl , \ + CCTK_INT IN tl , \ CCTK_REAL IN scale, \ CCTK_INT ARRAY IN srcs, \ CCTK_INT ARRAY IN tls, \ @@ -226,6 +232,7 @@ USES FUNCTION LinearCombination CCTK_INT FUNCTION Accelerator_NotifyDataModified \ (CCTK_POINTER_TO_CONST IN cctkGH, \ CCTK_INT ARRAY IN variables, \ + CCTK_INT ARRAY IN rls, \ CCTK_INT ARRAY IN tls, \ CCTK_INT IN nvariables, \ CCTK_INT IN on_device) diff --git a/src/Euler.c b/src/Euler.c index 0033cd1..294b8eb 100644 --- a/src/Euler.c +++ b/src/Euler.c @@ -72,6 +72,12 @@ void MoL_EulerAdd(CCTK_ARGUMENTS) totalsize *= cctk_ash[arraydim]; } + CCTK_INT rl = 0; + if (CCTK_IsFunctionAliased("GetRefinementLevel")) { + rl = GetRefinementLevel(cctkGH); + } + CCTK_INT tl = 0; + switch (*MoL_Intermediate_Step) { case 1: @@ -84,7 +90,7 @@ void MoL_EulerAdd(CCTK_ARGUMENTS) CCTK_INT const tls[] = {1, 0}; CCTK_REAL const facts[] = {1.0, CCTK_DELTA_TIME}; MoL_LinearCombination(cctkGH, - EvolvedVariableIndex[var], 0.0, + EvolvedVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrcs); } @@ -96,7 +102,7 @@ void MoL_EulerAdd(CCTK_ARGUMENTS) CCTK_INT const tls[] = {1, 0}; CCTK_REAL const facts[] = {1.0, CCTK_DELTA_TIME}; MoL_LinearCombination(cctkGH, - EvolvedArrayVariableIndex[var], 0.0, + EvolvedArrayVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrcs); } diff --git a/src/InitialCopy.c b/src/InitialCopy.c index 5cef13c..5a8979d 100644 --- a/src/InitialCopy.c +++ b/src/InitialCopy.c @@ -108,6 +108,12 @@ void MoL_InitialCopy(CCTK_ARGUMENTS) totalsize *= cctk_ash[arraydim]; } + int rl = 0; + if (CCTK_IsFunctionAliased("GetRefinementLevel")) { + rl = GetRefinementLevel(cctkGH); + } + int tl = 0; + for (var = 0; var < MoLNumEvolvedVariables; var++) { const int nsrc = 1; @@ -135,7 +141,7 @@ void MoL_InitialCopy(CCTK_ARGUMENTS) CCTK_WARN(0, "The grid function does not have storage assigned."); } - MoL_LinearCombination(cctkGH, EvolvedVariableIndex[var], 0.0, + MoL_LinearCombination(cctkGH, EvolvedVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrc); } @@ -522,6 +528,12 @@ void MoL_InitRHS(CCTK_ARGUMENTS) totalsize *= cctk_ash[arraydim]; } + CCTK_INT rl = 0; + if (CCTK_IsFunctionAliased("GetRefinementLevel")) { + rl = GetRefinementLevel(cctkGH); + } + CCTK_INT tl = 0; + for (var = 0; var < MoLNumEvolvedVariables; var++) { StorageOn = CCTK_QueryGroupStorageI(cctkGH, @@ -544,7 +556,7 @@ void MoL_InitRHS(CCTK_ARGUMENTS) CCTK_WARN(0, "The grid function does not have storage assigned."); } - MoL_LinearCombination(cctkGH, RHSVariableIndex[var], 0.0, + MoL_LinearCombination(cctkGH, RHSVariableIndex[var], rl, tl, 0.0, NULL, NULL, NULL, 0); } diff --git a/src/Operators.c b/src/Operators.c index 21ff8b2..2cbde16 100644 --- a/src/Operators.c +++ b/src/Operators.c @@ -1,5 +1,10 @@ #include "Operators.h" + +#include +#include + #include +#include #include #include @@ -13,11 +18,12 @@ static void -error_no_storage(int const var, int const tl) +error_no_storage(int const var, int const rl, int const tl) { char *const fullname = CCTK_FullName(var); - CCTK_VWarn(CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING, - "Variable %s, timelevel %d has no storage", fullname, tl); + CCTK_VError(__LINE__, __FILE__, CCTK_THORNSTRING, + "Variable %s, refinement level %d, timelevel %d has no storage", + fullname, rl, tl); free(fullname); return; } @@ -164,6 +170,16 @@ op_real_update_3(CCTK_REAL *restrict const varptr, @vtype int @vio inout @endvar + @var rl + @vdesc refinement level of all variables + @vtype int + @vio in + @endvar + @var tl + @vdesc time level of target variable + @vtype int + @vio in + @endvar @var scale @vdesc scale target by this @vtype CCTK_REAL @@ -200,19 +216,24 @@ op_real_update_3(CCTK_REAL *restrict const varptr, CCTK_INT MoL_LinearCombination(cGH const *const cctkGH, CCTK_INT const var, + CCTK_INT const rl, + CCTK_INT const tl, CCTK_REAL const scale, CCTK_INT const srcs[], CCTK_INT const tls[], CCTK_REAL const facts[], CCTK_INT const nsrcs) { + DECLARE_CCTK_PARAMETERS; + // Forward call to aliased function, if it is defined static int is_aliased = -1; if (is_aliased < 0) { is_aliased = CCTK_IsFunctionAliased("LinearCombination"); } if (is_aliased) { - return LinearCombination(cctkGH, var, scale, srcs, tls, facts, nsrcs); + return + LinearCombination(cctkGH, var, rl, tl, scale, srcs, tls, facts, nsrcs); } // Determine grid variable size @@ -231,12 +252,13 @@ MoL_LinearCombination(cGH const *const cctkGH, case CCTK_VARIABLE_REAL: { // Obtain pointer to variable data // TODO: check that all variable types are CCTK_REAL - CCTK_REAL *restrict const varptr = CCTK_VarDataPtrI(cctkGH, 0, var); - if (!varptr) error_no_storage(var, 0); + CCTK_REAL *restrict const varptr = CCTK_VarDataPtrI(cctkGH, tl, var); + if (!varptr) error_no_storage(var, rl, tl); CCTK_REAL const *restrict srcptrs[nsrcs]; for (int n=0; ncctk_timefac * beta}; MoL_LinearCombination(cctkGH, - EvolvedVariableIndex[var], 0.0, + EvolvedVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrcs); time = facts[0] * old_time + facts[1] * time_rhs; } @@ -169,7 +175,7 @@ CCTK_WARN(0, "not implemented"); if ((*MoL_Intermediate_Step) == MoL_Intermediate_Steps) { MoL_LinearCombination(cctkGH, - scratchspace_firstindex + var, 0.0, + scratchspace_firstindex + var, rl, tl, 0.0, NULL, NULL, NULL, 0); scratch_time = 0.0; } @@ -181,7 +187,7 @@ CCTK_WARN(0, "not implemented"); CCTK_INT const tls[] = {0}; CCTK_REAL const facts[] = {alpha}; MoL_LinearCombination(cctkGH, - scratchspace_firstindex + var, 1.0, + scratchspace_firstindex + var, rl, tl, 1.0, srcs, tls, facts, nsrcs); scratch_time += facts[0] * time; } @@ -193,7 +199,7 @@ CCTK_WARN(0, "not implemented"); CCTK_INT const tls[] = {0, 1}; CCTK_REAL const facts[] = {1.0, -4.0/3.0}; MoL_LinearCombination(cctkGH, - EvolvedVariableIndex[var], 1.0, + EvolvedVariableIndex[var], rl, tl, 1.0, srcs, tls, facts, nsrcs); time += facts[0] * scratch_time + facts[1] * old_time; } diff --git a/src/RK87.c b/src/RK87.c index 21b00ce..fb16b37 100644 --- a/src/RK87.c +++ b/src/RK87.c @@ -136,6 +136,12 @@ void MoL_RK87Add(CCTK_ARGUMENTS) /* First store (dt times) the rhs in the scratch array. */ + CCTK_INT rl = 0; + if (CCTK_IsFunctionAliased("GetRefinementLevel")) { + rl = GetRefinementLevel(cctkGH); + } + CCTK_INT tl = 0; + for (var = 0; var < MoLNumEvolvedVariables; var++) { int const step = MoL_Intermediate_Steps - (*MoL_Intermediate_Step); @@ -147,7 +153,7 @@ void MoL_RK87Add(CCTK_ARGUMENTS) CCTK_INT const tls[] = {0}; CCTK_REAL const facts[] = {(*Original_Delta_Time) / cctkGH->cctk_timefac}; MoL_LinearCombination(cctkGH, - scratchvarindex, 1.0, + scratchvarindex, rl, tl, 1.0, srcs, tls, facts, nsrcs); } @@ -176,7 +182,7 @@ void MoL_RK87Add(CCTK_ARGUMENTS) } MoL_LinearCombination(cctkGH, - EvolvedVariableIndex[var], 0.0, + EvolvedVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrcs); } else @@ -200,7 +206,7 @@ void MoL_RK87Add(CCTK_ARGUMENTS) } MoL_LinearCombination(cctkGH, - EvolvedVariableIndex[var], 0.0, + EvolvedVariableIndex[var], rl, tl, 0.0, srcs, tls, facts, nsrcs); } @@ -221,7 +227,7 @@ void MoL_RK87Add(CCTK_ARGUMENTS) } MoL_LinearCombination(cctkGH, - error_firstindex + var, 0.0, + error_firstindex + var, rl, tl, 0.0, srcs, tls, facts, nsrcs); } -- cgit v1.2.3