aboutsummaryrefslogtreecommitdiff
path: root/src/InitialCopy.c
diff options
context:
space:
mode:
authoreschnett <eschnett@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2013-01-22 21:00:28 +0000
committereschnett <eschnett@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2013-01-22 21:00:28 +0000
commita4ed7cfd73db6ef8418ffb47f4237cc529f8c710 (patch)
treebfd603e9eb86fc9ad01ee2b841013903752d1859 /src/InitialCopy.c
parent15f9511f821e350bc8c0c850b6f0a65ff786df87 (diff)
MoL Update
New integrator Euler. This is an explicit, first-order method, mostly useful for debugging. New integrators AB (Adams-Bashforth) with various orders. These are explicit integrators using several past timelevels to provide higher-order integration with a single RHS evaluation each. Introduce LinearCombination, a generic routine to calculate linear combinations of grid functions. This simplifies existing code, and can be overloaded if MoL should run on a device (e.g. with OpenCL). Replace cctk_lsh with cctk_ash where necessary. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@190 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/InitialCopy.c')
-rw-r--r--src/InitialCopy.c50
1 files changed, 14 insertions, 36 deletions
diff --git a/src/InitialCopy.c b/src/InitialCopy.c
index b421908..5cef13c 100644
--- a/src/InitialCopy.c
+++ b/src/InitialCopy.c
@@ -20,6 +20,7 @@
#include "cctk_Parameters.h"
#include "ExternalVariables.h"
+#include "Operators.h"
static const char *rcsid = "$Header$";
@@ -104,11 +105,15 @@ void MoL_InitialCopy(CCTK_ARGUMENTS)
totalsize = 1;
for (arraydim = 0; arraydim < cctk_dim; arraydim++)
{
- totalsize *= cctk_lsh[arraydim];
+ totalsize *= cctk_ash[arraydim];
}
for (var = 0; var < MoLNumEvolvedVariables; var++)
{
+ const int nsrc = 1;
+ const int srcs[1] = {EvolvedVariableIndex[var]};
+ const int tls[1] = {1};
+ const CCTK_REAL facts[1] = {1.0};
StorageOn = CCTK_QueryGroupStorageI(cctkGH,
CCTK_GroupIndexFromVarI(EvolvedVariableIndex[var]));
@@ -130,20 +135,8 @@ void MoL_InitialCopy(CCTK_ARGUMENTS)
CCTK_WARN(0, "The grid function does not have storage assigned.");
}
- PreviousVar = (CCTK_REAL const*)CCTK_VarDataPtrI(cctkGH, 1,
- EvolvedVariableIndex[var]);
- CurrentVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- EvolvedVariableIndex[var]);
- if (PreviousVar && CurrentVar)
- {
- memcpy(CurrentVar, PreviousVar, totalsize * sizeof(CCTK_REAL));
- }
- else
- {
- CCTK_VWarn(0,__LINE__,__FILE__,CCTK_THORNSTRING,"Null pointer for variable %s",
- CCTK_VarName(EvolvedVariableIndex[var]));
- }
-
+ MoL_LinearCombination(cctkGH, EvolvedVariableIndex[var], 0.0,
+ srcs, tls, facts, nsrc);
}
/* Set up the array sizes */
@@ -186,7 +179,7 @@ void MoL_InitialCopy(CCTK_ARGUMENTS)
arraytotalsize = 1;
for (arraydim = 0; arraydim < arraydata.dim; arraydim++)
{
- arraytotalsize *= arraydata.lsh[arraydim];
+ arraytotalsize *= arraydata.ash[arraydim];
}
ArrayScratchSizes[var] = arraytotalsize;
@@ -526,12 +519,11 @@ void MoL_InitRHS(CCTK_ARGUMENTS)
totalsize = 1;
for (arraydim = 0; arraydim < cctk_dim; arraydim++)
{
- totalsize *= cctk_lsh[arraydim];
+ totalsize *= cctk_ash[arraydim];
}
for (var = 0; var < MoLNumEvolvedVariables; var++)
{
-
StorageOn = CCTK_QueryGroupStorageI(cctkGH,
CCTK_GroupIndexFromVarI(RHSVariableIndex[var]));
@@ -552,22 +544,8 @@ void MoL_InitRHS(CCTK_ARGUMENTS)
CCTK_WARN(0, "The grid function does not have storage assigned.");
}
- RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0,
- RHSVariableIndex[var]);
- if (RHSVar)
- {
-#pragma omp parallel for
- for (index = 0; index < totalsize; index++)
- {
- RHSVar[index] = 0;
- }
- }
- else
- {
- CCTK_VWarn(0,__LINE__,__FILE__,CCTK_THORNSTRING,"Null pointer for variable %s",
- CCTK_VarName(RHSVariableIndex[var]));
- }
-
+ MoL_LinearCombination(cctkGH, RHSVariableIndex[var], 0.0,
+ NULL, NULL, NULL, 0);
}
for (var = 0; var < MoLNumEvolvedArrayVariables; var++)
@@ -588,7 +566,7 @@ void MoL_InitRHS(CCTK_ARGUMENTS)
arraytotalsize = 1;
for (arraydim = 0; arraydim < arraydata.dim; arraydim++)
{
- arraytotalsize *= arraydata.lsh[arraydim];
+ arraytotalsize *= arraydata.ash[arraydim];
}
if (arraytotalsize)
@@ -690,7 +668,7 @@ void MoL_FillAllLevels(CCTK_ARGUMENTS)
totalsize = 1;
for (arraydim = 0; arraydim < cctk_dim; arraydim++)
{
- totalsize *= cctk_lsh[arraydim];
+ totalsize *= cctk_ash[arraydim];
}
for (var = 0; var < MoLNumEvolvedVariables; var++)