aboutsummaryrefslogtreecommitdiff
path: root/src/InitialCopy.c
diff options
context:
space:
mode:
authorhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2003-04-23 14:59:00 +0000
committerhawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2003-04-23 14:59:00 +0000
commitfd1b7d175f18cb5f1d168c8a9e95508b41ae39e9 (patch)
tree438f83a13ed6c0dcd80b495222d57b5ee59c898c /src/InitialCopy.c
parent1100eb9f3017973b6cfa5b6bff3306773ff3e532 (diff)
The Method of Lines thorn (version 2 - see below).
MoL provides generic integration methods for multiple thorns simultaneously. By providing a layer between the driver and evolution thorns, this should mean that some technical issues to do with mesh refinement can be ignored. It also allows you to choose different evolution methods (in time, at least). But the primary purpose is to unambiguously evolve models in different thorns at the same time. This is version 2 - the one that will work with mesh refinement. It's a straight copy of HawkeCVS/Public/CactusMoL2/MoL2 and I haven't checked that it will work "as is" with the new name (one of the reasons it goes into Alpha). At the moment the only evolution method guaranteed to work with mesh refinement is ICN, although RK2 should. The generic RK methods will do something that will be subtly wrong... Note that the "old" way of registering variables (through the functions declared in header files) is still there. As soon as function aliasing settles down, this will be removed. Also to be done: Better documentation Tidy up the code (especially the debugging statements) Optimize Add various useful time evolution methods. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@2 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/InitialCopy.c')
-rw-r--r--src/InitialCopy.c312
1 files changed, 312 insertions, 0 deletions
diff --git a/src/InitialCopy.c b/src/InitialCopy.c
new file mode 100644
index 0000000..9d3c19a
--- /dev/null
+++ b/src/InitialCopy.c
@@ -0,0 +1,312 @@
+ /*@@
+ @file InitialCopy.c
+ @date Sun May 26 04:43:06 2002
+ @author Ian Hawke
+ @desc
+ Performs the initial copy from the previous timelevel to the
+ current. This is required because the driver has rotated the
+ timelevels, but the physics thorns are expecting data in the
+ current.
+ @enddesc
+ @version $Header$
+ @@*/
+
+#include "cctk.h"
+#include "cctk_Arguments.h"
+#include "cctk_Parameters.h"
+
+#include "ExternalVariables.h"
+
+static const char *rcsid = "$Header$";
+
+CCTK_FILEVERSION(CactusMoL2_MoL2_InitialCopy_c);
+
+/********************************************************************
+ ********************* Local Data Types ***********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ***************** Scheduled Routine Prototypes *********************
+ ********************************************************************/
+
+void MoL_InitialCopy(CCTK_ARGUMENTS);
+
+/********************************************************************
+ ********************* Other Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Data *****************************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* External Routines **********************
+ ********************************************************************/
+
+void MoL_InitialCopy(CCTK_ARGUMENTS)
+{
+
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+
+ CCTK_INT var;
+ CCTK_INT index;
+ CCTK_INT i,j,k;
+ CCTK_INT totalsize;
+
+ CCTK_REAL *CurrentVar;
+ CCTK_REAL *PreviousVar;
+ CCTK_REAL *ScratchVar;
+ CCTK_INT StorageOn;
+
+ totalsize = cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2];
+
+ for (var = 0; var < MoLNumEvolvedVariables; var++)
+ {
+
+ StorageOn = CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(EvolvedVariableIndex[var]));
+
+ if (StorageOn < 0)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for index %i",
+ EvolvedVariableIndex[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(EvolvedVariableIndex[var]));
+ CCTK_WARN(0, "The grid function does not have storage assigned.");
+ }
+
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
+ EvolvedVariableIndex[var]);
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ EvolvedVariableIndex[var]);
+ /*
+ if (EvolvedVariableIndex[var] == CCTK_VarIndex("adm_bssn::ADM_BS_gxx"))
+ {
+ printf("Initial copy: Pointer to time levels %ld %ld\n",PreviousVar,
+ CurrentVar);
+ printf("Level %d:First point of first two levels is %g and %g.\n",
+ cctk_levfac[0],PreviousVar[0],CurrentVar[0]);
+ }
+ */
+ if (PreviousVar && CurrentVar)
+ {
+ memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
+ }
+ else
+ {
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(EvolvedVariableIndex[var]));
+ }
+
+ /*
+ // if (EvolvedVariableIndex[var] == CCTK_VarIndex("adm_bssn::adm_bs_gxx"))
+ if (EvolvedVariableIndex[var] == CCTK_VarIndex("wavetoymol::phi"))
+ {
+ printf("Init copy: Lev %d Current %g.\n",cctk_levfac[0],
+ CurrentVar[0]);
+ }
+ */
+ }
+
+ /*
+ Now the Save and Restore variables. Shift the data in the
+ current level to the scratch space, then do the copy
+ */
+
+ for (var = 0; var < MoLNumSandRVariables; var++)
+ {
+
+ StorageOn = CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(SandRVariableIndex[var]));
+
+ if (StorageOn < 0)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for index %i",
+ SandRVariableIndex[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(SandRVariableIndex[var]));
+ CCTK_WARN(0, "The grid function does not have storage assigned.");
+ }
+
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
+ SandRVariableIndex[var]);
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ SandRVariableIndex[var]);
+ /* ScratchVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ SandRScratchSpace[var]);
+ */
+ ScratchVar = &SandRScratchSpace[var*totalsize];
+ /*
+ printf("Pointers for the SandR vars are to %ld, %ld and %ld.\n",
+ PreviousVar, CurrentVar, ScratchVar);
+
+ printf("Init1:Variable %s, current %g, previous %g, scratch %g\n",
+ CCTK_VarName(SandRVariableIndex[var]), CurrentVar[0],
+ PreviousVar[0], ScratchVar[0]);
+ */
+ if (PreviousVar && CurrentVar && ScratchVar)
+ {
+ memcpy(ScratchVar, CurrentVar, totalsize * sizeof(CCTK_REAL));
+ memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
+ }
+ else
+ {
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(SandRScratchSpace[var]));
+ }
+ /*
+ printf("Init2:Variable %s, current %g, previous %g, scratch %g\n",
+ CCTK_VarName(SandRVariableIndex[var]), CurrentVar[0],
+ PreviousVar[0], ScratchVar[0]);
+ */
+ }
+
+ /*
+ Now do the constrained variables.
+ */
+
+ for (var = 0; var < MoLNumConstrainedVariables; var++)
+ {
+
+ StorageOn = CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(ConstrainedVariableIndex[var]));
+
+ if (StorageOn < 0)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for index %i",
+ ConstrainedVariableIndex[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(ConstrainedVariableIndex[var]));
+ CCTK_WARN(0, "The grid function does not have storage assigned.");
+ }
+
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
+ ConstrainedVariableIndex[var]);
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ ConstrainedVariableIndex[var]);
+
+ if (PreviousVar && CurrentVar)
+ {
+ memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
+ }
+ else
+ {
+ CCTK_VWarn(0,__LINE__,__FILE__,"MoL","Null pointer for variable %s",
+ CCTK_VarName(ConstrainedVariableIndex[var]));
+ }
+
+ }
+
+ return;
+}
+
+/********************************************************************
+ ********************* Local Routines *************************
+ ********************************************************************/
+
+/*
+Hack
+*/
+
+void MoL_FillAllLevels(CCTK_ARGUMENTS)
+{
+
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+
+ CCTK_INT var, level;
+ CCTK_INT totalsize;
+
+ CCTK_REAL *CurrentVar;
+ CCTK_REAL *PreviousVar;
+
+ totalsize = cctk_lsh[0]*cctk_lsh[1]*cctk_lsh[2];
+
+// printf("\n\n Fill all levels \n\n");
+
+
+ for (var = 0; var < MoLNumEvolvedVariables; var++)
+ {
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ EvolvedVariableIndex[var]);
+ for (level = 1; level < CCTK_NumTimeLevelsFromVarI(EvolvedVariableIndex[var]); level++)
+ {
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, level,
+ EvolvedVariableIndex[var]);
+ if (PreviousVar)
+ {
+ memcpy(PreviousVar, CurrentVar, totalsize * sizeof(CCTK_REAL));
+ }
+ /*
+ if (EvolvedVariableIndex[var] == CCTK_VarIndex("wavetoymol::phi"))
+ {
+ printf("All levels copy: Lev %d tl %d Current %g.\n",cctk_levfac[0],
+ level,PreviousVar[0]);
+ }
+ */
+ }
+ }
+
+
+ for (var = 0; var < MoLNumConstrainedVariables; var++)
+ {
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ ConstrainedVariableIndex[var]);
+ for (level = 1; level < CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(ConstrainedVariableIndex[var])); level++)
+ {
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, level,
+ ConstrainedVariableIndex[var]);
+ if (PreviousVar)
+ {
+ memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
+ }
+ }
+ }
+
+
+ for (var = 0; var < MoLNumSandRVariables; var++)
+ {
+ CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ SandRVariableIndex[var]);
+ for (level = 1; level < CCTK_QueryGroupStorage(cctkGH,
+ CCTK_GroupNameFromVarI(SandRVariableIndex[var])); level++)
+ {
+ PreviousVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, level,
+ SandRVariableIndex[var]);
+ if (PreviousVar)
+ {
+ memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
+ }
+ }
+ }
+
+ CCTK_VInfo(CCTK_THORNSTRING,
+ "The maximum number of evolved variables is %d. %d are registered",
+ MoL_Num_Evolved_Vars,MoLNumEvolvedVariables);
+
+ CCTK_VInfo(CCTK_THORNSTRING,
+ "The maximum number of constrained variables is %d. %d are registered",
+ MoL_Num_Constrained_Vars,MoLNumConstrainedVariables);
+
+ CCTK_VInfo(CCTK_THORNSTRING,
+ "The maximum number of SandR variables is %d. %d are registered",
+ MoL_Num_SaveAndRestore_Vars,MoLNumSandRVariables);
+
+
+ return;
+}