aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2004-03-23 17:53:29 +0000
committerschnetter <schnetter@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2004-03-23 17:53:29 +0000
commitf7f56a61852ea7125d2d0b539a1cd5b6849e87b3 (patch)
treec9b6b412257f34be96cb6402f53cf85dcb84cd97
parentb6b896ef68bebf49b1b2cc241986da174b64e205 (diff)
Initialise the RHS with zero before scheduling the physics routines.
This is necessary because MoL switches off boundary prolongation when there is mesh refinement, leaving certain grid points undefined. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@54 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
-rw-r--r--schedule.ccl9
-rw-r--r--src/InitialCopy.c182
2 files changed, 191 insertions, 0 deletions
diff --git a/schedule.ccl b/schedule.ccl
index 83832c5..8d71041 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -444,6 +444,15 @@ if (CCTK_Equals(ODE_Method,"ICN-avg"))
} "Averages the time levels for the averaging ICN method"
}
+##############################################
+### The time step is initialised to zero ###
+##############################################
+
+schedule MoL_InitRHS IN MoL_Step BEFORE MoL_CalcRHS
+{
+ LANG: C
+} "Initialise the RHS functions"
+
#####################################################
### The group where all the physics takes place ###
#####################################################
diff --git a/src/InitialCopy.c b/src/InitialCopy.c
index 89c67ae..05a8dfb 100644
--- a/src/InitialCopy.c
+++ b/src/InitialCopy.c
@@ -37,6 +37,8 @@ CCTK_FILEVERSION(AlphaThorns_MoL_InitialCopy_c);
void MoL_InitialCopy(CCTK_ARGUMENTS);
+void MoL_InitRHS(CCTK_ARGUMENTS);
+
void MoL_FillAllLevels(CCTK_ARGUMENTS);
void MoL_ReportNumberVariables(CCTK_ARGUMENTS);
@@ -53,6 +55,21 @@ void MoL_ReportNumberVariables(CCTK_ARGUMENTS);
********************* External Routines **********************
********************************************************************/
+ /*@@
+ @routine MoL_InitialCopy
+ @date ???
+ @author Ian Hawke
+ @desc
+ Copy the previous time level to the current time level.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
void MoL_InitialCopy(CCTK_ARGUMENTS)
{
@@ -390,6 +407,171 @@ void MoL_InitialCopy(CCTK_ARGUMENTS)
}
/*@@
+ @routine MoL_InitRHS
+ @date Tue Mar 23 2004
+ @author Erik Schnetter
+ @desc
+ Initialise all RHS variables with zero.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+void MoL_InitRHS(CCTK_ARGUMENTS)
+{
+
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ cGroupDynamicData arraydata;
+ CCTK_INT groupindex, ierr;
+ CCTK_INT arraytotalsize, arraydim;
+
+ CCTK_INT var;
+ CCTK_INT index;
+/* CCTK_INT i,j,k; */
+ CCTK_INT totalsize;
+
+ CCTK_REAL *RHSVar;
+ CCTK_INT StorageOn;
+
+ /* FIXME */
+
+#ifdef MOLDOESCOMPLEX
+
+ CCTK_COMPLEX *RHSComplexVar;
+
+#endif
+
+ totalsize = cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2];
+
+ for (var = 0; var < MoLNumEvolvedVariables; var++)
+ {
+
+ StorageOn = CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(RHSVariableIndex[var]));
+
+ if (StorageOn < 0)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for index %i",
+ RHSVariableIndex[var]);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ else if (StorageOn == 0) {
+#ifdef MOLDEBUG
+ printf("Aargh! Vars %d var %d index %d name %s\n",
+ MoLNumEvolvedVariables, var, RHSVariableIndex[var],
+ CCTK_VarName(RHSVariableIndex[var]));
+#endif
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for GF %s",
+ CCTK_VarName(RHSVariableIndex[var]));
+ CCTK_WARN(0, "The grid function does not have storage assigned.");
+ }
+
+ RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ RHSVariableIndex[var]);
+ if (RHSVar)
+ {
+ for (index = 0; index < totalsize; index++)
+ {
+ RHSVar[index] = 0;
+ }
+ }
+ else
+ {
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(RHSVariableIndex[var]));
+ }
+
+ }
+
+ for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
+ {
+ RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ RHSArrayVariableIndex[var]);
+
+ groupindex = CCTK_GroupIndexFromVarI(RHSArrayVariableIndex[var]);
+ ierr = CCTK_GroupDynamicData(cctkGH, groupindex,
+ &arraydata);
+ if (ierr)
+ {
+ CCTK_VWarn(0, __LINE__, __FILE__, "MoL",
+ "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];
+ }
+
+ if (RHSVar)
+ {
+ for (index = 0; index < arraytotalsize; index++)
+ {
+ RHSVar[index] = 0;
+ }
+ }
+ else
+ {
+ printf("The pointer is %p (rhs)\n.",
+ RHSVar);
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(RHSArrayVariableIndex[var]));
+ }
+
+ }
+
+ /* FIXME */
+
+#ifdef MOLDOESCOMPLEX
+
+ for (var = 0; var < MoLNumEvolvedComplexVariables; var++)
+ {
+
+ StorageOn = CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(RHSComplexVariableIndex[var]));
+
+ if (StorageOn < 0)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for index %i",
+ RHSComplexVariableIndex[var]);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ else if (StorageOn == 0) {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for GF %s",
+ CCTK_VarName(RHSComplexVariableIndex[var]));
+ CCTK_WARN(0, "The grid function does not have storage assigned.");
+ }
+
+ RHSComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0,
+ RHSComplexVariableIndex[var]);
+ if (RHSVar)
+ {
+ for (index = 0; index < totalsize; index++)
+ {
+ RHSVar[index] = 0;
+ }
+ }
+ else
+ {
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(RHSComplexVariableIndex[var]));
+ }
+
+ }
+
+#endif
+
+ return;
+}
+
+ /*@@
@routine MoL_FillAllLevels
@date Fri Apr 25 16:11:18 2003
@author Ian Hawke