aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <>2002-01-11 16:37:00 +0000
committerschnetter <>2002-01-11 16:37:00 +0000
commit3dbc34d97aad4664f5006f25d230d8e0d8f5eb58 (patch)
treec82eda31644e44129a9fc9b5201ac3638d5fd8ec
parentf34dfb0a25cb9f18b87d3a89c45cbc099bcc7518 (diff)
Allow grid extents to be specified by integer grid point locations as
Allow grid extents to be specified by integer grid point locations as well as by physical coordinates. darcs-hash:20020111163713-07bb3-fa1a86617d4805de7faba76ee1612f36d3a41b6d.gz
-rw-r--r--Carpet/CarpetRegrid/param.ccl99
-rw-r--r--Carpet/CarpetRegrid/src/regrid.cc87
-rw-r--r--Carpet/CarpetRegrid/src/regrid.hh6
3 files changed, 165 insertions, 27 deletions
diff --git a/Carpet/CarpetRegrid/param.ccl b/Carpet/CarpetRegrid/param.ccl
index eb2d7300f..e6902bc80 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.2 2002/01/11 17:19:47 schnetter Exp $
+# $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/param.ccl,v 1.3 2002/01/11 17:37:13 schnetter Exp $
@@ -21,15 +21,102 @@ CCTK_INT regrid_every "Regrid every n time steps"
KEYWORD refined_regions "Regions where the grid is refined"
{
- "none" :: "Don't refine"
- "centre" :: "Refine (arbitrarily) the centre of the grid only"
- "manual" :: "Refine the regions specified by l[123][xyz]{min,max}"
- "automatic" :: "Refine automatically"
+ "none" :: "Don't refine"
+ "centre" :: "Refine around the centre of the grid only"
+ "manual-gridpoints" :: "Refine the regions specified by integer grid points l[123]i[xyz]{min,max}"
+ "manual-coordinates" :: "Refine the regions specified by coordinates l[123][xyz]{min,max}"
+ "automatic" :: "Refine automatically"
} "centre"
-# Region specifications for manual refinement
+# Region specifications for manual gridpoint refinement
+
+CCTK_INT l1ixmin "Lower boundary of level 1 box in x-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l1iymin "Lower boundary of level 1 box in y-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l1izmin "Lower boundary of level 1 box in z-direction"
+{
+ : :: ""
+} 0
+
+CCTK_INT l1ixmax "Upper boundary of level 1 box in x-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l1iymax "Upper boundary of level 1 box in y-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l1izmax "Upper boundary of level 1 box in z-direction"
+{
+ : :: ""
+} -1
+
+
+
+CCTK_INT l2ixmin "Lower boundary of level 2 box in x-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l2iymin "Lower boundary of level 2 box in y-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l2izmin "Lower boundary of level 2 box in z-direction"
+{
+ : :: ""
+} 0
+
+CCTK_INT l2ixmax "Upper boundary of level 2 box in x-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l2iymax "Upper boundary of level 2 box in y-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l2izmax "Upper boundary of level 2 box in z-direction"
+{
+ : :: ""
+} -1
+
+
+
+CCTK_INT l3ixmin "Lower boundary of level 3 box in x-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l3iymin "Lower boundary of level 3 box in y-direction"
+{
+ : :: ""
+} 0
+CCTK_INT l3izmin "Lower boundary of level 3 box in z-direction"
+{
+ : :: ""
+} 0
+
+CCTK_INT l3ixmax "Upper boundary of level 3 box in x-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l3iymax "Upper boundary of level 3 box in y-direction"
+{
+ : :: ""
+} -1
+CCTK_INT l3izmax "Upper boundary of level 3 box in z-direction"
+{
+ : :: ""
+} -1
+
+
+
+# Region specifications for manual coordinate refinement
CCTK_REAL l1xmin "Lower boundary of level 1 box in x-direction"
{
diff --git a/Carpet/CarpetRegrid/src/regrid.cc b/Carpet/CarpetRegrid/src/regrid.cc
index e8243fc4f..fc2a62d75 100644
--- a/Carpet/CarpetRegrid/src/regrid.cc
+++ b/Carpet/CarpetRegrid/src/regrid.cc
@@ -19,7 +19,7 @@
#include "carpet.hh"
#include "regrid.hh"
-static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.cc,v 1.5 2002/01/11 17:19:48 schnetter Exp $";
+static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.cc,v 1.6 2002/01/11 17:37:13 schnetter Exp $";
@@ -79,7 +79,22 @@ namespace CarpetRegrid {
MakeRegions_RefineCentre (cctkGH, refinement_levels, bbl);
- } else if (CCTK_EQUALS(refined_regions, "manual")) {
+ } else if (CCTK_EQUALS(refined_regions, "manual-gridpoints")) {
+
+ if (refinement_levels > 4) {
+ CCTK_WARN (0, "Cannot currently specify manual refinement regions for more than 4 refinement levels");
+ }
+ assert (refinement_levels<=4);
+ vector<vect<int,dim> > lower(3), upper(3);
+ lower[0] = vect<int,dim> (l1ixmin, l1iymin, l1izmin);
+ upper[0] = vect<int,dim> (l1ixmax, l1iymax, l1izmax);
+ lower[1] = vect<int,dim> (l2ixmin, l2iymin, l2izmin);
+ upper[1] = vect<int,dim> (l2ixmax, l2iymax, l2izmax);
+ lower[2] = vect<int,dim> (l3ixmin, l3iymin, l3izmin);
+ upper[2] = vect<int,dim> (l3ixmax, l3iymax, l3izmax);
+ MakeRegions_AsSpecified (cctkGH, refinement_levels, lower, upper, bbl);
+
+ } else if (CCTK_EQUALS(refined_regions, "manual-coordinates")) {
if (refinement_levels > 4) {
CCTK_WARN (0, "Cannot currently specify manual refinement regions for more than 4 refinement levels");
@@ -209,35 +224,22 @@ namespace CarpetRegrid {
- void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
- const vector<vect<CCTK_REAL,dim> > lower,
- const vector<vect<CCTK_REAL,dim> > upper,
- list<bbox<int,dim> >& bbl)
+ static void
+ MakeRegions_AsSpecified_OneLevel (const cGH* cctkGH, const int reflevels,
+ const vect<int,dim> ilower,
+ const vect<int,dim> iupper,
+ list<bbox<int,dim> >& bbl)
{
assert (bbl.empty());
if (reflevel+1 >= reflevels) return;
- vect<CCTK_REAL,dim> global_lower, global_upper;
- for (int d=0; d<dim; ++d) {
- const int ierr = CCTK_CoordRange
- (cctkGH, &global_lower[d], &global_upper[d], d+1, 0, "cart3d");
- if (ierr<0) {
- global_lower[d] = 0;
- global_upper[d] = 1;
- }
- }
- const vect<int,dim> global_extent (hh->baseextent.upper() - hh->baseextent.lower() + hh->baseextent.stride() * (dd->lghosts + dd->ughosts));
-
const vect<int,dim> rstr = hh->baseextent.stride();
const vect<int,dim> rlb = hh->baseextent.lower();
const vect<int,dim> rub = hh->baseextent.upper();
const int rl = reflevel+1;
- const vect<int,dim> ilower = map(floor, (lower[rl-1] - global_lower) / (global_upper - global_lower) * global_extent + 0.5);
- const vect<int,dim> iupper = map(floor, (upper[rl-1] - global_lower) / (global_upper - global_lower) * global_extent + 0.5);
-
const int levfac = ipow(hh->reffact, rl);
assert (all (rstr % levfac == 0));
const vect<int,dim> str (rstr / levfac);
@@ -264,6 +266,51 @@ namespace CarpetRegrid {
+ void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
+ const vector<vect<int,dim> > lower,
+ const vector<vect<int,dim> > upper,
+ list<bbox<int,dim> >& bbl)
+ {
+ if (reflevel+1 >= reflevels) return;
+
+ const int rl = reflevel+1;
+
+ const vect<int,dim> ilower = lower[rl-1];
+ const vect<int,dim> iupper = upper[rl-1];
+
+ MakeRegions_AsSpecified_OneLevel (cctkGH, reflevels, ilower, iupper, bbl);
+ }
+
+
+
+ void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
+ const vector<vect<CCTK_REAL,dim> > lower,
+ const vector<vect<CCTK_REAL,dim> > upper,
+ list<bbox<int,dim> >& bbl)
+ {
+ if (reflevel+1 >= reflevels) return;
+
+ vect<CCTK_REAL,dim> global_lower, global_upper;
+ for (int d=0; d<dim; ++d) {
+ const int ierr = CCTK_CoordRange
+ (cctkGH, &global_lower[d], &global_upper[d], d+1, 0, "cart3d");
+ if (ierr<0) {
+ global_lower[d] = 0;
+ global_upper[d] = 1;
+ }
+ }
+ const vect<int,dim> global_extent (hh->baseextent.upper() - hh->baseextent.lower() + hh->baseextent.stride() * (dd->lghosts + dd->ughosts));
+
+ const int rl = reflevel+1;
+
+ const vect<int,dim> ilower = map(floor, (lower[rl-1] - global_lower) / (global_upper - global_lower) * global_extent + 0.5);
+ const vect<int,dim> iupper = map(floor, (upper[rl-1] - global_lower) / (global_upper - global_lower) * global_extent + 0.5);
+
+ MakeRegions_AsSpecified_OneLevel (cctkGH, reflevels, ilower, iupper, bbl);
+ }
+
+
+
static void
MakeRegions_Adaptively_Recombine (list<bbox<int,dim> >& bbl1,
list<bbox<int,dim> >& bbl2,
diff --git a/Carpet/CarpetRegrid/src/regrid.hh b/Carpet/CarpetRegrid/src/regrid.hh
index 0491816b8..cef62a2fb 100644
--- a/Carpet/CarpetRegrid/src/regrid.hh
+++ b/Carpet/CarpetRegrid/src/regrid.hh
@@ -1,4 +1,4 @@
-// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.hh,v 1.2 2002/01/11 17:19:49 schnetter Exp $
+// $Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetRegrid/src/regrid.hh,v 1.3 2002/01/11 17:37:14 schnetter Exp $
#ifndef REGRID_HH
#define REGRID_HH
@@ -39,6 +39,10 @@ namespace CarpetRegrid {
list<bbox<int,dim> >& bbl);
void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
+ const vector<vect<int,dim> > lower,
+ const vector<vect<int,dim> > upper,
+ list<bbox<int,dim> >& bbl);
+ void MakeRegions_AsSpecified (const cGH* cctkGH, const int reflevels,
const vector<vect<CCTK_REAL,dim> > lower,
const vector<vect<CCTK_REAL,dim> > upper,
list<bbox<int,dim> >& bbl);