aboutsummaryrefslogtreecommitdiff
path: root/src/RK2.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/RK2.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/RK2.c')
-rw-r--r--src/RK2.c140
1 files changed, 140 insertions, 0 deletions
diff --git a/src/RK2.c b/src/RK2.c
new file mode 100644
index 0000000..34aa068
--- /dev/null
+++ b/src/RK2.c
@@ -0,0 +1,140 @@
+ /*@@
+ @file RK2.c
+ @date Sun May 26 04:13:45 2002
+ @author Ian Hawke
+ @desc
+ A specialized second order Runge-Kutta time integrator. This is
+ the integrator that Shu refers to as the optimal TVD second
+ order method (see reference in documentation). It is equivalent
+ to Heun's predictor-corrector method, or the MacCormack method.
+ @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_RK2_c);
+
+/********************************************************************
+ ********************* Local Data Types ***********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ***************** Scheduled Routine Prototypes *********************
+ ********************************************************************/
+
+void MoL_RK2Add(CCTK_ARGUMENTS);
+
+/********************************************************************
+ ********************* Other Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Data *****************************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* External Routines **********************
+ ********************************************************************/
+
+ /*@@
+ @routine MoL_RK2Add
+ @date Sun May 26 04:17:23 2002
+ @author Ian Hawke
+ @desc
+ Performs second order Runge-Kutta time integration.
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+void MoL_RK2Add(CCTK_ARGUMENTS)
+{
+
+ DECLARE_CCTK_ARGUMENTS;
+ DECLARE_CCTK_PARAMETERS;
+
+ CCTK_INT index, var;
+ CCTK_INT totalsize;
+ CCTK_REAL *OldVar;
+ CCTK_REAL *UpdateVar;
+ CCTK_REAL *RHSVar;
+
+ /*
+ printf("Inside RK2.\nStep %d.\nRefinement %d.\nTimestep %g.\nSpacestep %g.\nTime %g\n",
+ MoL_Intermediate_Steps - *MoL_Intermediate_Step + 1,
+ *cctk_levfac,
+ cctk_delta_time,
+ CCTK_DELTA_SPACE(0),
+ cctk_time);
+ */
+
+ totalsize = cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2];
+
+ switch (*MoL_Intermediate_Step)
+ {
+
+ case 2:
+ {
+ for (var = 0; var < MoLNumEvolvedVariables; var++)
+ {
+ UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ EvolvedVariableIndex[var]);
+ RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ RHSVariableIndex[var]);
+
+ for (index = 0; index < totalsize; index++)
+ {
+ UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index];
+ }
+ }
+ break;
+ }
+ case 1:
+ {
+ for (var = 0; var < MoLNumEvolvedVariables; var++)
+ {
+ OldVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1,
+ EvolvedVariableIndex[var]);
+ UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ EvolvedVariableIndex[var]);
+ RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
+ RHSVariableIndex[var]);
+
+ for (index = 0; index < totalsize; index++)
+ {
+ UpdateVar[index] = 0.5 * (OldVar[index] + UpdateVar[index]) +
+ CCTK_DELTA_TIME * RHSVar[index];
+ }
+ }
+ break;
+ }
+ default:
+ {
+ CCTK_WARN(0, "RK2 expects MoL_Intermediate_Step to be in [1,2]. This should be caught at ParamCheck - bug Ian!");
+ break;
+ }
+
+ }
+
+ return;
+
+}
+
+/********************************************************************
+ ********************* Local Routines *************************
+ ********************************************************************/