aboutsummaryrefslogtreecommitdiff
path: root/src/ChangeType.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/ChangeType.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/ChangeType.c')
-rw-r--r--src/ChangeType.c572
1 files changed, 572 insertions, 0 deletions
diff --git a/src/ChangeType.c b/src/ChangeType.c
new file mode 100644
index 0000000..0894361
--- /dev/null
+++ b/src/ChangeType.c
@@ -0,0 +1,572 @@
+ /*@@
+ @file ChangeType.c
+ @date Thu May 30 16:16:40 2002
+ @author Ian Hawke
+ @desc
+ The external functions called (via function aliasing) by physics
+ thorns to tell MoL that they want these GFs to be treated as a
+ different type to the original declaration
+ @enddesc
+ @version $Header$
+ @@*/
+
+#include "cctk.h"
+
+#include "ExternalVariables.h"
+
+static const char *rcsid = "$Header$";
+
+CCTK_FILEVERSION(CactusMoL2_MoL2_ChangeType_c);
+
+/********************************************************************
+ ********************* Local Data Types ***********************
+ ********************************************************************/
+
+#define MOL_UNKNOWN_VARTYPE 0
+#define MOL_EVOLVED_VARTYPE 1
+#define MOL_CONSTRAINED_VARTYPE 2
+#define MOL_SANDR_VARTYPE 3
+
+/********************************************************************
+ ********************* Local Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ***************** Scheduled Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Other Routine Prototypes *********************
+ ********************************************************************/
+
+CCTK_INT MoL_ChangeToEvolved(CCTK_INT EvolvedIndex, CCTK_INT RHSIndex);
+
+CCTK_INT MoL_ChangeToConstrained(CCTK_INT ConstrainedIndex);
+
+CCTK_INT MoL_ChangeToSaveAndRestore(CCTK_INT SandRIndex);
+
+CCTK_INT MoL_ChangeToNone(CCTK_INT RemoveIndex);
+
+/********************************************************************
+ ********************* Local Data *****************************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* External Routines **********************
+ ********************************************************************/
+
+ /*@@
+ @routine MoL_ChangeToEvolved
+ @date Thu May 30 16:45:30 2002
+ @author Ian Hawke
+ @desc
+ Changes a variable to evolved type. Checks to see which type it was
+ before and does the bookkeeping on the index arrays.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+CCTK_INT MoL_ChangeToEvolved(CCTK_INT EvolvedIndex, CCTK_INT RHSIndex)
+{
+
+ CCTK_INT index, usedindex;
+ CCTK_INT vartype; /* See the defines at the top of file */
+ CCTK_INT timelevs;
+
+ vartype = 0;
+ usedindex = -1;
+
+ for (index = 0; (index < MoLNumEvolvedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_EVOLVED_VARTYPE *
+ (EvolvedVariableIndex[index] == EvolvedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumConstrainedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_CONSTRAINED_VARTYPE *
+ (ConstrainedVariableIndex[index] == EvolvedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumSandRVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_SANDR_VARTYPE *
+ (SandRVariableIndex[index] == EvolvedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ switch (vartype)
+ {
+
+ case MOL_UNKNOWN_VARTYPE:
+
+ {
+ timelevs = CCTK_NumTimeLevelsFromVarI(EvolvedIndex);
+ if (timelevs < 1)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved.", EvolvedIndex);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ else if (timelevs == 1)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved.", EvolvedIndex);
+ CCTK_WARN(0, "The index passed only has a single timelevel.");
+ }
+ timelevs = CCTK_NumTimeLevelsFromVarI(RHSIndex);
+ if (timelevs < 1) {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved (RHS GF).", RHSIndex);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ EvolvedVariableIndex[MoLNumEvolvedVariables] = EvolvedIndex;
+ RHSVariableIndex[MoLNumEvolvedVariables] = RHSIndex;
+ MoLNumEvolvedVariables++;
+ break;
+ }
+
+ case MOL_EVOLVED_VARTYPE:
+
+ {
+ RHSVariableIndex[usedindex] = RHSIndex;
+ break;
+ }
+
+ case MOL_CONSTRAINED_VARTYPE:
+
+ {
+ timelevs = CCTK_NumTimeLevelsFromVarI(RHSIndex);
+ if (timelevs < 1) {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved type from constrained.", RHSIndex);
+ CCTK_WARN(0, "The RHS index passed does not correspond to a GF.");
+ }
+ for (index = usedindex; index < MoLNumConstrainedVariables - 1; index++)
+ {
+ ConstrainedVariableIndex[index] = ConstrainedVariableIndex[index+1];
+ }
+ MoLNumConstrainedVariables--;
+ EvolvedVariableIndex[index] = EvolvedIndex;
+ RHSVariableIndex[index] = RHSIndex;
+ MoLNumEvolvedVariables++;
+ break;
+ }
+
+ case MOL_SANDR_VARTYPE:
+
+ {
+ timelevs = CCTK_NumTimeLevelsFromVarI(RHSIndex);
+ if (timelevs < 1) {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved type from save and restore.", RHSIndex);
+ CCTK_WARN(0, "The RHS index passed does not correspond to a GF.");
+ }
+ for (index = usedindex; index < MoLNumSandRVariables - 1; index++)
+ {
+ SandRVariableIndex[index] = SandRVariableIndex[index+1];
+ }
+ MoLNumSandRVariables--;
+ EvolvedVariableIndex[index] = EvolvedIndex;
+ RHSVariableIndex[index] = RHSIndex;
+ MoLNumEvolvedVariables++;
+ break;
+ }
+
+ default:
+
+ {
+ CCTK_WARN(0, "Something is seriously wrong in ChangeType.c! Case out of range in switch statement.");
+ }
+
+ }
+
+ return 0;
+
+}
+
+ /*@@
+ @routine MoL_ChangeToConstrained
+ @date Thu May 30 16:47:08 2002
+ @author Ian Hawke
+ @desc
+ Changes a variable to be constrained.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+CCTK_INT MoL_ChangeToConstrained(CCTK_INT ConstrainedIndex)
+{
+
+ CCTK_INT index, usedindex;
+ CCTK_INT vartype; /* See the defines at the top of file */
+ CCTK_INT timelevs;
+
+ vartype = 0;
+ usedindex = -1;
+
+ for (index = 0; (index < MoLNumEvolvedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_EVOLVED_VARTYPE *
+ (EvolvedVariableIndex[index] == ConstrainedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumConstrainedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_CONSTRAINED_VARTYPE *
+ (ConstrainedVariableIndex[index] == ConstrainedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumSandRVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_SANDR_VARTYPE *
+ (SandRVariableIndex[index] == ConstrainedIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ switch (vartype)
+ {
+
+ case MOL_UNKNOWN_VARTYPE:
+
+ {
+ timelevs = CCTK_NumTimeLevelsFromVarI(ConstrainedIndex);
+ if (timelevs < 1)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved.", ConstrainedIndex);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ else if (timelevs == 1)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to evolved.", ConstrainedIndex);
+ CCTK_WARN(0, "The index passed only has a single timelevel.");
+ }
+ ConstrainedVariableIndex[MoLNumConstrainedVariables] = ConstrainedIndex;
+ MoLNumConstrainedVariables++;
+ break;
+ }
+
+ case MOL_EVOLVED_VARTYPE:
+
+ {
+
+ break;
+ }
+
+ case MOL_CONSTRAINED_VARTYPE:
+
+ {
+ break;
+ }
+
+ case MOL_SANDR_VARTYPE:
+
+ {
+ for (index = usedindex; index < MoLNumSandRVariables - 1; index++)
+ {
+ SandRVariableIndex[index] = SandRVariableIndex[index+1];
+ }
+ MoLNumSandRVariables--;
+ ConstrainedVariableIndex[index] = ConstrainedIndex;
+ MoLNumConstrainedVariables++;
+ break;
+ }
+
+ default:
+
+ {
+ CCTK_WARN(0, "Something is seriously wrong in ChangeType.c! Case out of range in switch statement.");
+ }
+
+ }
+
+ return 0;
+
+}
+
+ /*@@
+ @routine MoL_ChangeToSaveAndRestore
+ @date Thu May 30 16:50:36 2002
+ @author Ian Hawke
+ @desc
+ Changes the variable type to save and restore.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+CCTK_INT MoL_ChangeToSaveAndRestore(CCTK_INT SandRIndex)
+{
+
+ CCTK_INT index, usedindex;
+ CCTK_INT vartype; /* See the defines at the top of file */
+ CCTK_INT timelevs;
+
+ vartype = 0;
+ usedindex = -1;
+
+ for (index = 0; (index < MoLNumEvolvedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_EVOLVED_VARTYPE *
+ (EvolvedVariableIndex[index] == SandRIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumConstrainedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_CONSTRAINED_VARTYPE *
+ (ConstrainedVariableIndex[index] == SandRIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumSandRVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_SANDR_VARTYPE *
+ (SandRVariableIndex[index] == SandRIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ switch (vartype)
+ {
+
+ case MOL_UNKNOWN_VARTYPE:
+
+ {
+ timelevs = CCTK_NumTimeLevelsFromVarI(SandRIndex);
+ if (timelevs < 1)
+ {
+ CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning whilst trying to change variable index %i to save and restore.", SandRIndex);
+ CCTK_WARN(0, "The index passed does not correspond to a GF.");
+ }
+ SandRVariableIndex[MoLNumSandRVariables] = SandRIndex;
+ MoLNumSandRVariables++;
+ break;
+ }
+
+ case MOL_EVOLVED_VARTYPE:
+
+ {
+
+ break;
+ }
+
+ case MOL_CONSTRAINED_VARTYPE:
+
+ {
+ for (index = usedindex; index < MoLNumConstrainedVariables - 1; index++)
+ {
+ ConstrainedVariableIndex[index] = ConstrainedVariableIndex[index+1];
+ }
+ MoLNumConstrainedVariables--;
+ SandRVariableIndex[index] = SandRIndex;
+ MoLNumSandRVariables++;
+ break;
+ }
+
+ case MOL_SANDR_VARTYPE:
+
+ {
+ break;
+ }
+
+ default:
+
+ {
+ CCTK_WARN(0, "Something is seriously wrong in ChangeType.c! Case out of range in switch statement.");
+ }
+
+ }
+
+ return 0;
+
+}
+
+ /*@@
+ @routine MoL_ChangeToNone
+ @date Thu May 30 16:53:46 2002
+ @author Ian Hawke
+ @desc
+ Changes a variable type to none (i.e., MoL no longer considers it).
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+CCTK_INT MoL_ChangeToNone(CCTK_INT RemoveIndex)
+{
+
+ CCTK_INT index, usedindex;
+ CCTK_INT vartype; /* See the defines at the top of file */
+ CCTK_INT timelevs;
+
+ vartype = 0;
+ usedindex = -1;
+
+ for (index = 0; (index < MoLNumEvolvedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_EVOLVED_VARTYPE *
+ (EvolvedVariableIndex[index] == RemoveIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumConstrainedVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_CONSTRAINED_VARTYPE *
+ (ConstrainedVariableIndex[index] == RemoveIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ for (index = 0; (index < MoLNumSandRVariables)&&(!vartype); index++)
+ {
+ vartype = MOL_SANDR_VARTYPE *
+ (SandRVariableIndex[index] == RemoveIndex);
+ if (vartype)
+ {
+ usedindex = index;
+ }
+ }
+
+ switch (vartype)
+ {
+
+ case MOL_UNKNOWN_VARTYPE:
+
+ {
+ break;
+ }
+
+ case MOL_EVOLVED_VARTYPE:
+
+ {
+ for (index = usedindex; index < MoLNumEvolvedVariables - 1; index++)
+ {
+ EvolvedVariableIndex[index] = EvolvedVariableIndex[index+1];
+ RHSVariableIndex[index] = RHSVariableIndex[index+1];
+ }
+ MoLNumEvolvedVariables--;
+ break;
+ }
+
+ case MOL_CONSTRAINED_VARTYPE:
+
+ {
+ for (index = usedindex; index < MoLNumConstrainedVariables - 1; index++)
+ {
+ ConstrainedVariableIndex[index] = ConstrainedVariableIndex[index+1];
+ }
+ MoLNumConstrainedVariables--;
+ break;
+ }
+
+ case MOL_SANDR_VARTYPE:
+
+ {
+ for (index = usedindex; index < MoLNumSandRVariables - 1; index++)
+ {
+ SandRVariableIndex[index] = SandRVariableIndex[index+1];
+ }
+ MoLNumSandRVariables--;
+ break;
+ }
+
+ default:
+
+ {
+ CCTK_WARN(0, "Something is seriously wrong in ChangeType.c! Case out of range in switch statement.");
+ }
+
+ }
+
+ return 0;
+
+}
+
+
+
+/*
+ Fortran wrappers for the above functions.
+ Should be replaced by using function aliasing eventually.
+ */
+
+void CCTK_FCALL CCTK_FNAME(MoL_ChangeToEvolved)(int *ierr,
+ CCTK_INT *EvolvedIndex,
+ CCTK_INT *RHSIndex)
+{
+ *ierr = MoL_ChangeToEvolved(*EvolvedIndex, *RHSIndex);
+ return;
+}
+
+void CCTK_FCALL CCTK_FNAME(MoL_ChangeToConstrained)(int *ierr,
+ CCTK_INT *ConstrainedIndex)
+{
+ *ierr = MoL_ChangeToConstrained(*ConstrainedIndex);
+ return;
+}
+
+void CCTK_FCALL CCTK_FNAME(MoL_ChangeToSaveAndRestore)(int *ierr,
+ CCTK_INT *EvolvedIndex)
+{
+ *ierr = MoL_ChangeToSaveAndRestore(*EvolvedIndex);
+ return;
+}
+
+void CCTK_FCALL CCTK_FNAME(MoL_ChangeToNone)(int *ierr,
+ CCTK_INT *EvolvedIndex)
+{
+ *ierr = MoL_ChangeToNone(*EvolvedIndex);
+ return;
+}
+
+/********************************************************************
+ ********************* Local Routines *************************
+ ********************************************************************/