aboutsummaryrefslogtreecommitdiff
path: root/Carpet
diff options
context:
space:
mode:
authorhawke <>2004-01-13 12:50:00 +0000
committerhawke <>2004-01-13 12:50:00 +0000
commit299eaed96c11757c8dbff59e8a9d256eea29071b (patch)
tree09cd61ea5c9e2dd523c75c1b976bcc5b49e8de36 /Carpet
parentd409df948e48cccc2145830bfad683cb226c09a2 (diff)
Changes to CarpetRegrid.
Changes to CarpetRegrid. - Does progressive MR correctly again. - Will take a criteria for the number of levels to activate from an aliased function if requested, allowing for adaptive refinement in time. darcs-hash:20040113125005-58737-15008754eaae7b272c5a144333241be8ab6bbde2.gz
Diffstat (limited to 'Carpet')
-rw-r--r--Carpet/CarpetRegrid/interface.ccl7
-rw-r--r--Carpet/CarpetRegrid/param.ccl9
-rw-r--r--Carpet/CarpetRegrid/src/regrid.cc59
3 files changed, 64 insertions, 11 deletions
diff --git a/Carpet/CarpetRegrid/interface.ccl b/Carpet/CarpetRegrid/interface.ccl
index 4c2504112..7fc64b1cc 100644
--- a/Carpet/CarpetRegrid/interface.ccl
+++ b/Carpet/CarpetRegrid/interface.ccl
@@ -1,5 +1,5 @@
# Interface definition for thorn CarpetRegrid
-# $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/interface.ccl,v 1.3 2003/09/04 16:23:22 tradke Exp $
+# $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/interface.ccl,v 1.4 2004/01/13 13:50:05 hawke Exp $
implements: CarpetRegrid
@@ -13,3 +13,8 @@ uses include header: vect.hh
uses include header: gf.hh
uses include header: gh.hh
+
+CCTK_INT FUNCTION RegridLevel (CCTK_POINTER_TO_CONST IN cctkGH, \
+ CCTK_INT IN current_reflevel, \
+ CCTK_INT IN current_max_reflevel)
+USES FUNCTION RegridLevel
diff --git a/Carpet/CarpetRegrid/param.ccl b/Carpet/CarpetRegrid/param.ccl
index 63255ac6a..9f2dc9bf9 100644
--- a/Carpet/CarpetRegrid/param.ccl
+++ b/Carpet/CarpetRegrid/param.ccl
@@ -1,5 +1,5 @@
# Parameter definitions for thorn CarpetRegrid
-# $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/param.ccl,v 1.12 2003/12/19 16:01:55 schnetter Exp $
+# $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/param.ccl,v 1.13 2004/01/13 13:50:05 hawke Exp $
@@ -268,3 +268,10 @@ CCTK_STRING errorvar "Name of grid function that contains the error" STEERABLE=a
{
".*" :: "must be the name of a grid function"
} ""
+
+
+# Should we set the regridding criteria from an aliased function?
+
+BOOLEAN regrid_from_function "Regridding criteria from the aliased function?"
+{
+} "no"
diff --git a/Carpet/CarpetRegrid/src/regrid.cc b/Carpet/CarpetRegrid/src/regrid.cc
index 3144197ca..9d1731862 100644
--- a/Carpet/CarpetRegrid/src/regrid.cc
+++ b/Carpet/CarpetRegrid/src/regrid.cc
@@ -24,7 +24,7 @@
#include "regrid.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.cc,v 1.31 2003/11/21 12:51:56 schnetter Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.cc,v 1.32 2004/01/13 13:49:38 hawke Exp $";
CCTK_FILEVERSION(Carpet_CarpetRegrid_regrid_cc);
}
@@ -39,7 +39,15 @@ namespace CarpetRegrid {
using namespace std;
using namespace Carpet;
-
+ int initial_refinement_levels = 0; // Store the initial number of
+ // levels to refine for correct
+ // progressive MR.
+ // Note that this does not need
+ // to be checkpointed as it is
+ // reset at the beginning of every
+ // iteration.
+ int last_regridded_iteration = -1; // The last iteration where any
+ // regridding was done.
int CarpetRegridStartup ()
{
@@ -71,14 +79,47 @@ namespace CarpetRegrid {
// Return if we want to regrid during initial data only, and this
// is not the time for initial data
if (regrid_every == 0 && cctkGH->cctk_iteration != 0) return 0;
-
+
// Return if we want to regrid regularly, but not at this time
if (regrid_every > 0 && cctkGH->cctk_iteration != 0
&& cctkGH->cctk_iteration % regrid_every != 0) {
return 0;
}
- int newnumlevels = refinement_levels + activate_newlevels_on_regrid;
+ // Find the number of levels that should be active.
+ int newnumlevels;
+
+ // Is this the first time that we have regridded this iteration?
+ if (cctkGH->cctk_iteration > last_regridded_iteration) {
+ last_regridded_iteration = cctkGH->cctk_iteration;
+ initial_refinement_levels = refinement_levels;
+ }
+
+ if (cctkGH->cctk_iteration == 0) {
+ newnumlevels = refinement_levels;
+ } else {
+ if (regrid_from_function) {
+ if (CCTK_IsFunctionAliased("RegridLevel")) {
+ int tempnewnumlevels = 0;
+ newnumlevels = 0;
+ BEGIN_MGLEVEL_LOOP(cctkGH) {
+ tempnewnumlevels =
+ RegridLevel(cctkGH, reflevel, refinement_levels);
+ newnumlevels = max(tempnewnumlevels, newnumlevels);
+ } END_MGLEVEL_LOOP;
+ } else {
+ CCTK_WARN(1, "No thorn has provided the function "
+ "\"RegridLevel\". Regridding will not be done.");
+ }
+ } else {
+ // The standard progressive MR number of levels.
+ newnumlevels = initial_refinement_levels +
+ activate_newlevels_on_regrid;
+ }
+ }
+
+ newnumlevels = min(newnumlevels, maxreflevels);
+
if ( ( newnumlevels >= 1) && (newnumlevels <= maxreflevels )) {
char numlevelstring[10];
sprintf(numlevelstring,"%d",newnumlevels);
@@ -94,7 +135,7 @@ namespace CarpetRegrid {
} else if (CCTK_EQUALS(refined_regions, "centre")) {
- MakeRegions_RefineCentre (cctkGH, refinement_levels, bbl, obl);
+ MakeRegions_RefineCentre (cctkGH, newnumlevels, bbl, obl);
} else if (CCTK_EQUALS(refined_regions, "manual-gridpoints")) {
@@ -109,7 +150,7 @@ namespace CarpetRegrid {
upper[1] = ivect (l2ixmax, l2iymax, l2izmax);
lower[2] = ivect (l3ixmin, l3iymin, l3izmin);
upper[2] = ivect (l3ixmax, l3iymax, l3izmax);
- MakeRegions_AsSpecified (cctkGH, refinement_levels, lower, upper,
+ MakeRegions_AsSpecified (cctkGH, newnumlevels, lower, upper,
bbl, obl);
} else if (CCTK_EQUALS(refined_regions, "manual-coordinates")) {
@@ -125,7 +166,7 @@ namespace CarpetRegrid {
upper[1] = rvect (l2xmax, l2ymax, l2zmax);
lower[2] = rvect (l3xmin, l3ymin, l3zmin);
upper[2] = rvect (l3xmax, l3ymax, l3zmax);
- MakeRegions_AsSpecified (cctkGH, refinement_levels, lower, upper,
+ MakeRegions_AsSpecified (cctkGH, newnumlevels, lower, upper,
bbl, obl);
} else if (CCTK_EQUALS(refined_regions, "manual-gridpoint-list")) {
@@ -161,7 +202,7 @@ namespace CarpetRegrid {
}
}
- MakeRegions_AsSpecified (cctkGH, refinement_levels, bbss, obss,
+ MakeRegions_AsSpecified (cctkGH, newnumlevels, bbss, obss,
bbl, obl);
} else if (CCTK_EQUALS(refined_regions, "manual-coordinate-list")) {
@@ -197,7 +238,7 @@ namespace CarpetRegrid {
}
}
- MakeRegions_AsSpecified (cctkGH, refinement_levels, bbss, obss,
+ MakeRegions_AsSpecified (cctkGH, newnumlevels, bbss, obss,
bbl, obl);
} else if (CCTK_EQUALS(refined_regions, "automatic")) {