aboutsummaryrefslogtreecommitdiff
path: root/Carpet/Carpet/src
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@aei.mpg.de>2005-02-01 22:58:00 +0000
committerErik Schnetter <schnetter@aei.mpg.de>2005-02-01 22:58:00 +0000
commit0914bc88c7aea61eedf470a903c8fa252bdc97dd (patch)
treec11153db9c2272c28d1861a42e543942c49a7dd1 /Carpet/Carpet/src
parentf9fe6d4b5a573027170f45784dae4b094160c546 (diff)
global: Change the way in which the grid hierarchy is stored
Change the way in which the grid hierarchy is stored. The new hierarchy is map mglevel reflevel component timelevel i.e., mglevel moved from the bottom to almost the top. This is because mglevel used to be a true multigrid level, but is now meant to be a convergence level. Do not allocate all storage all the time. Allow storage to be switched on an off per refinement level (and for a single mglevel, which prompted the change above). Handle storage management with CCTK_{In,De}creaseGroupStorage instead of CCTK_{En,Dis}ableGroupStorage. darcs-hash:20050201225827-891bb-eae3b6bd092ae8d6b5e49be84c6f09f0e882933e.gz
Diffstat (limited to 'Carpet/Carpet/src')
-rw-r--r--Carpet/Carpet/src/Recompose.cc128
-rw-r--r--Carpet/Carpet/src/SetupGH.cc58
-rw-r--r--Carpet/Carpet/src/Storage.cc161
-rw-r--r--Carpet/Carpet/src/functions.hh10
-rw-r--r--Carpet/Carpet/src/modes.cc16
-rw-r--r--Carpet/Carpet/src/variables.hh1
6 files changed, 209 insertions, 165 deletions
diff --git a/Carpet/Carpet/src/Recompose.cc b/Carpet/Carpet/src/Recompose.cc
index 8e4618c5a..59b5ecc73 100644
--- a/Carpet/Carpet/src/Recompose.cc
+++ b/Carpet/Carpet/src/Recompose.cc
@@ -65,47 +65,50 @@ namespace Carpet {
- void CheckRegions (const gh::rexts & bbsss,
+ void CheckRegions (const gh::mexts & bbsss,
const gh::rbnds & obss,
const gh::rprocs& pss)
{
- // At least one level
+ // At least one multigrid level
if (bbsss.size() == 0) {
- CCTK_WARN (0, "I cannot set up a grid hierarchy with zero refinement levels.");
+ CCTK_WARN (0, "I cannot set up a grid hierarchy with zero multigrid levels.");
}
assert (bbsss.size() > 0);
+ // At least one refinement level
+ if (bbsss.at(0).size() == 0) {
+ CCTK_WARN (0, "I cannot set up a grid hierarchy with zero refinement levels.");
+ }
+ assert (bbsss.at(0).size() > 0);
// At most maxreflevels levels
- if ((int)bbsss.size() > maxreflevels) {
+ if ((int)bbsss.at(0).size() > maxreflevels) {
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"I cannot set up a grid hierarchy with more than Carpet::max_refinement_levels refinement levels. I found Carpet::max_refinement_levels=%d, while %d levels were requested.",
- (int)maxreflevels, (int)bbsss.size());
- }
- assert ((int)bbsss.size() <= maxreflevels);
- for (int rl=0; rl<(int)bbsss.size(); ++rl) {
- // No empty levels
- assert (bbsss.at(rl).size() > 0);
- for (int c=0; c<(int)bbsss.at(rl).size(); ++c) {
- // At least one multigrid level
- assert (bbsss.at(rl).at(c).size() > 0);
- for (int ml=0; ml<(int)bbsss.at(rl).at(c).size(); ++ml) {
+ (int)maxreflevels, (int)bbsss.at(0).size());
+ }
+ assert ((int)bbsss.at(0).size() <= maxreflevels);
+ for (int ml=0; ml<(int)bbsss.size(); ++ml) {
+ for (int rl=0; rl<(int)bbsss.at(0).size(); ++rl) {
+ // No empty levels
+ assert (bbsss.at(ml).at(rl).size() > 0);
+ for (int c=0; c<(int)bbsss.at(ml).at(rl).size(); ++c) {
// Check sizes
// Do allow processors with zero grid points
// assert (all(bbsss.at(rl).at(c).at(ml).lower() <= bbsssi.at(rl).at(c).at(ml).upper()));
// Check strides
const int str = ipow(reffact, maxreflevels-rl-1) * ipow(mgfact, ml);
- assert (all(bbsss.at(rl).at(c).at(ml).stride() == str));
+ assert (all(bbsss.at(ml).at(rl).at(c).stride() == str));
// Check alignments
- assert (all(bbsss.at(rl).at(c).at(ml).lower() % str == 0));
- assert (all(bbsss.at(rl).at(c).at(ml).upper() % str == 0));
+ assert (all(bbsss.at(ml).at(rl).at(c).lower() % str == 0));
+ assert (all(bbsss.at(ml).at(rl).at(c).upper() % str == 0));
}
}
}
- assert (pss.size() == bbsss.size());
- assert (obss.size() == bbsss.size());
- for (int rl=0; rl<(int)bbsss.size(); ++rl) {
- assert (obss.at(rl).size() == bbsss.at(rl).size());
- assert (pss.at(rl).size() == bbsss.at(rl).size());
+ assert (pss.size() == bbsss.at(0).size());
+ assert (obss.size() == bbsss.at(0).size());
+ for (int rl=0; rl<(int)bbsss.at(0).size(); ++rl) {
+ assert (obss.at(rl).size() == bbsss.at(0).at(rl).size());
+ assert (pss.at(rl).size() == bbsss.at(0).at(rl).size());
}
}
@@ -128,7 +131,7 @@ namespace Carpet {
bool did_change = false;
BEGIN_MAP_LOOP(cgh, CCTK_GF) {
- gh::rexts bbsss = vhh.at(map)->extents();
+ gh::mexts bbsss = vhh.at(map)->extents();
gh::rbnds obss = vhh.at(map)->outer_boundaries();
gh::rprocs pss = vhh.at(map)->processors();
@@ -177,13 +180,13 @@ namespace Carpet {
{
CCTK_INFO ("New grid structure (grid points):");
cout << " Refinement level " << reflevel << ", map " << map << endl;
- for (int rl=0; rl<hh.reflevels(); ++rl) {
- for (int c=0; c<hh.components(rl); ++c) {
- for (int ml=0; ml<hh.mglevels(rl,c); ++ml) {
+ for (int ml=0; ml<hh.mglevels(); ++ml) {
+ for (int rl=0; rl<hh.reflevels(); ++rl) {
+ for (int c=0; c<hh.components(rl); ++c) {
const int convfact = ipow(mgfact, ml);
const int levfact = ipow(reffact, rl);
- const ivect lower = hh.extents().at(rl).at(c).at(ml).lower();
- const ivect upper = hh.extents().at(rl).at(c).at(ml).upper();
+ const ivect lower = hh.extents().at(ml).at(rl).at(c).lower();
+ const ivect upper = hh.extents().at(ml).at(rl).at(c).upper();
assert (all(lower * levfact % maxreflevelfact == 0));
assert (all(upper * levfact % maxreflevelfact == 0));
assert (all(((upper - lower) * levfact / maxreflevelfact)
@@ -205,13 +208,13 @@ namespace Carpet {
}
CCTK_INFO ("New grid structure (coordinates):");
- for (int rl=0; rl<hh.reflevels(); ++rl) {
- for (int c=0; c<hh.components(rl); ++c) {
- for (int ml=0; ml<hh.mglevels(rl,c); ++ml) {
+ for (int ml=0; ml<hh.mglevels(); ++ml) {
+ for (int rl=0; rl<hh.reflevels(); ++rl) {
+ for (int c=0; c<hh.components(rl); ++c) {
const rvect origin = origin_space.at(0);
const rvect delta = delta_space;
- const ivect lower = hh.extents().at(rl).at(c).at(ml).lower();
- const ivect upper = hh.extents().at(rl).at(c).at(ml).upper();
+ const ivect lower = hh.extents().at(ml).at(rl).at(c).lower();
+ const ivect upper = hh.extents().at(ml).at(rl).at(c).upper();
const int convfact = ipow(mgfact, ml);
const int levfact = ipow(reffact, rl);
cout << " [" << ml << "][" << rl << "][" << m << "][" << c << "]"
@@ -230,7 +233,7 @@ namespace Carpet {
void OutputGridStructure (const cGH * const cgh,
const int m,
- const gh::rexts & bbsss,
+ const gh::mexts & bbsss,
const gh::rbnds & obss,
const gh::rprocs& pss)
{
@@ -274,13 +277,13 @@ namespace Carpet {
file << "iteration " << cgh->cctk_iteration << endl;
file << "maps " << maps << endl;
- file << m << " reflevels " << bbsss.size() << endl;
- for (int rl=0; rl<(int)bbsss.size(); ++rl) {
- file << m << " " << rl << " components " << bbsss.at(rl).size() << endl;
- for (int c=0; c<(int)bbsss.at(rl).size(); ++c) {
- file << m << " " << rl << " " << c << " mglevels " << bbsss.at(rl).at(c).size() << endl;
- for (int ml=0; ml<(int)bbsss.at(rl).at(c).size(); ++ml) {
- file << m << " " << rl << " " << c << " " << ml << " " << pss.at(rl).at(c) << " " << bbsss.at(rl).at(c).at(ml) << obss.at(rl).at(c) << endl;
+ file << m << " mglevels " << bbsss.size() << endl;
+ for (int ml=0; ml<(int)bbsss.size(); ++ml) {
+ file << m << " " << ml << " reflevels " << bbsss.at(ml).size() << endl;
+ for (int rl=0; rl<(int)bbsss.at(ml).size(); ++rl) {
+ file << m << " " << ml << " " << rl << " components " << bbsss.at(ml).at(rl).size() << endl;
+ for (int c=0; c<(int)bbsss.at(ml).at(rl).size(); ++c) {
+ file << m << " " << ml << " " << rl << " " << c << " " << pss.at(rl).at(c) << " " << bbsss.at(ml).at(rl).at(c) << obss.at(rl).at(c) << endl;
}
}
}
@@ -829,19 +832,42 @@ namespace Carpet {
}
void MakeMultigridBoxes (const cGH* cgh,
- vector<ibbox> const & bbs,
- vector<bbvect> const & obs,
- vector<vector<ibbox> >& bbss)
+ vector<vector<ibbox> > const & bbss,
+ vector<vector<bbvect> > const & obss,
+ vector<vector<vector<ibbox> > > & bbsss)
{
- assert (bbs.size() == obs.size());
- ibbox base;
- for (int c=0; c<(int)bbs.size(); ++c) {
- base = base.expanded_containing(bbs.at(c));
+ assert (bbss.size() == obss.size());
+ for (int rl=0; rl<(int)bbss.size(); ++rl) {
+ assert (bbss.at(rl).size() == obss.at(rl).size());
}
- bbss.resize(bbs.size());
- for (int c=0; c<(int)bbs.size(); ++c) {
- MakeMultigridBoxes (cgh, base, bbs.at(c), obs.at(c), bbss.at(c));
+
+ bbsss.resize (mglevels);
+ for (int ml=0; ml<mglevels; ++ml) {
+ bbsss.at(ml).resize (bbss.size());
+ for (int rl=0; rl<(int)bbss.size(); ++rl) {
+ bbsss.at(ml).at(rl).resize (bbss.at(rl).size());
+ }
}
+
+ for (int rl=0; rl<(int)bbss.size(); ++rl) {
+
+ ibbox base;
+ for (int c=0; c<(int)bbss.at(rl).size(); ++c) {
+ base = base.expanded_containing(bbss.at(rl).at(c));
+ }
+
+ for (int c=0; c<(int)bbss.at(rl).size(); ++c) {
+ vector<ibbox> mg_bbs;
+ MakeMultigridBoxes
+ (cgh, base, bbss.at(rl).at(c), obss.at(rl).at(c), mg_bbs);
+ assert ((int)mg_bbs.size() == mglevels);
+ for (int ml=0; ml<mglevels; ++ml) {
+ bbsss.at(ml).at(rl).at(c) = mg_bbs.at(ml);
+ }
+
+ } // for c
+
+ } // for rl
}
} // namespace Carpet
diff --git a/Carpet/Carpet/src/SetupGH.cc b/Carpet/Carpet/src/SetupGH.cc
index e97e8fd83..e7f06fb54 100644
--- a/Carpet/Carpet/src/SetupGH.cc
+++ b/Carpet/Carpet/src/SetupGH.cc
@@ -596,18 +596,18 @@ namespace Carpet {
vector<int> ps;
SplitRegions (cgh, bbs, obs, ps);
- // Create all multigrid levels
- vector<vector<ibbox> > bbss;
- MakeMultigridBoxes (cgh, bbs, obs, bbss);
-
// Only one refinement level
- vector<vector<vector<ibbox> > > bbsss(1);
+ vector<vector<ibbox> > bbss(1);
vector<vector<bbvect> > obss(1);
vector<vector<int> > pss(1);
- bbsss.at(0) = bbss;
+ bbss.at(0) = bbs;
obss.at(0) = obs;
pss.at(0) = ps;
+ // Create all multigrid levels
+ vector<vector<vector<ibbox> > > bbsss;
+ MakeMultigridBoxes (cgh, bbss, obss, bbsss);
+
// Check the regions
CheckRegions (bbsss, obss, pss);
@@ -769,11 +769,11 @@ namespace Carpet {
}
void print_grid_structure (vector<gh*> & vhh, int m) {
- const int rl = 0;
- for (int c=0; c<vhh.at(m)->components(rl); ++c) {
- for (int ml=0; ml<vhh.at(m)->mglevels(rl,c); ++ml) {
- const ivect lower = vhh.at(m)->extents().at(rl).at(c).at(ml).lower();
- const ivect upper = vhh.at(m)->extents().at(rl).at(c).at(ml).upper();
+ for (int ml=0; ml<vhh.at(m)->mglevels(); ++ml) {
+ const int rl = 0;
+ for (int c=0; c<vhh.at(m)->components(rl); ++c) {
+ const ivect lower = vhh.at(m)->extents().at(ml).at(rl).at(c).lower();
+ const ivect upper = vhh.at(m)->extents().at(ml).at(rl).at(c).upper();
const int convfact = ipow(mgfact, ml);
assert (all(lower % maxreflevelfact == 0));
assert (all(upper % maxreflevelfact == 0));
@@ -952,9 +952,17 @@ namespace Carpet {
} else {
SplitRegions_Automatic (cgh, bbs, obs, ps);
}
+
+ // Only one refinement level
+ vector<vector<ibbox> > bbss(1);
+ vector<vector<bbvect> > obss(1);
+ vector<vector<int> > pss(1);
+ bbss.at(0) = bbs;
+ obss.at(0) = obs;
+ pss.at(0) = ps;
// Create all multigrid levels
- vector<vector<ibbox> > bbss (bbs.size());
+ vector<vector<vector<ibbox> > > bbsss (mglevels);
ivect amgfact;
iivect aoffset;
for (int d=0; d<dim; ++d) {
@@ -962,17 +970,17 @@ namespace Carpet {
aoffset[d][0] = 0;
aoffset[d][1] = convoffsets[d];
}
- for (size_t c=0; c<bbs.size(); ++c) {
- bbss.at(c).resize (mglevels);
- bbss.at(c).at(0) = bbs.at(c);
- for (int ml=1; ml<mglevels; ++ml) {
+ bbsss.at(0) = bbss;
+ for (int ml=1; ml<mglevels; ++ml) {
+ const int rl = 0;
+ for (int c=0; c<(int)bbss.at(rl).size(); ++c) {
// this base
ivect const baselo = ivect(0);
ivect const baselo1 = baselo;
// next finer grid
- ivect const flo = bbss.at(c).at(ml-1).lower();
- ivect const fhi = bbss.at(c).at(ml-1).upper();
- ivect const fstr = bbss.at(c).at(ml-1).stride();
+ ivect const flo = bbsss.at(ml-1).at(rl).at(c).lower();
+ ivect const fhi = bbsss.at(ml-1).at(rl).at(c).upper();
+ ivect const fstr = bbsss.at(ml-1).at(rl).at(c).stride();
// this grid
ivect const str = fstr * amgfact;
ivect const lo = flo + xpose(obs.at(c))[0].ifthen
@@ -983,17 +991,9 @@ namespace Carpet {
ivect(0));
ivect const lo1 = baselo1 + (lo - baselo1 + str - 1) / str * str;
ivect const hi1 = lo1 + (hi - lo1) / str * str;
- bbss.at(c).at(ml) = ibbox(lo1, hi1, str);
+ bbsss.at(ml).at(rl).at(c) = ibbox(lo1, hi1, str);
}
}
-
- // Only one refinement level
- vector<vector<vector<ibbox> > > bbsss(1);
- vector<vector<bbvect> > obss(1);
- vector<vector<int> > pss(1);
- bbsss.at(0) = bbss;
- obss.at(0) = obs;
- pss.at(0) = ps;
// Recompose for this map
char * const groupname = CCTK_GroupName (group);
@@ -1027,7 +1027,7 @@ namespace Carpet {
arrdata.at(group).at(m).data.resize(CCTK_NumVarsInGroupI(group));
for (int var=0; var<(int)arrdata.at(group).at(m).data.size(); ++var) {
- arrdata.at(group).at(m).data.at(var) = 0;
+ arrdata.at(group).at(m).data.at(var) = NULL;
}
}
diff --git a/Carpet/Carpet/src/Storage.cc b/Carpet/Carpet/src/Storage.cc
index 0e606d9bd..c0cccde03 100644
--- a/Carpet/Carpet/src/Storage.cc
+++ b/Carpet/Carpet/src/Storage.cc
@@ -35,14 +35,19 @@ namespace Carpet {
assert (timelevels);
for (int n=0; n<n_groups; ++n) {
assert (groups[n] >= 0 and groups[n] < CCTK_NumGroups());
-#if 0
- for (int nn=0; nn<n; ++nn) {
- assert (groups[nn] != groups[n]);
- }
-#endif
+ // TODO: timelevels[n] can also be -1; in that case, all time
+ // levels should be activated / deactivated
assert (timelevels[n] >= 0);
}
+ bool const can_do = is_meta_mode() or is_global_mode() or is_level_mode();
+ bool const all_ml = is_meta_mode();
+ int const min_ml = all_ml ? 0 : mglevel;
+ int const max_ml = all_ml ? mglevels : mglevel+1;
+ bool const all_rl = is_meta_mode() or is_global_mode();
+ int const min_rl = all_rl ? 0 : reflevel;
+ int const max_rl = all_rl ? reflevels : reflevel+1;
+
int total_num_timelevels = 0;
for (int n=0; n<n_groups; ++n) {
@@ -68,76 +73,93 @@ namespace Carpet {
// Record previous number of allocated time levels
if (status) {
- status[n] = groupdata.at(group).activetimelevels;
+ status[n] = groupdata.at(group).info.activetimelevels;
}
// Only do something if the number of time levels actually needs
// to be changed -- do nothing otherwise
- if ((inc and timelevels[n] > groupdata.at(group).activetimelevels)
- or
- (! inc and timelevels[n] < groupdata.at(group).activetimelevels))
- {
+
+ const bool do_increase
+ = inc and timelevels[n] > groupdata.at(group).info.activetimelevels;
+ const bool do_decrease
+ = ! inc and timelevels[n] < groupdata.at(group).info.activetimelevels;
+ if (do_increase or do_decrease) {
- // Set the new number of active time levels
- groupdata.at(group).activetimelevels = timelevels[n];
+ if (! can_do) {
+ char * const groupname = CCTK_GroupName (group);
+ char const * const modestring
+ = (is_meta_mode() ? "meta"
+ : is_global_mode() ? "global"
+ : is_level_mode() ? "level"
+ : is_singlemap_mode() ? "singlemap"
+ : is_local_mode() ? "local"
+ : NULL);
+ CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Cannot change storage for group \"%s\" in %s mode",
+ groupname, modestring);
+ free (groupname);
+ }
+ assert (can_do);
- // There is a difference between the Cactus time levels and
- // the Carpet time levels. If there are n time levels, then
- // the Cactus time levels are numbered 0 ... n-1, with the
- // current time level being 0. In Carpet, the time levels are
- // numbered -(n-1) ... 0, where the current time level is also
- // 0.
- const int tmin = - timelevels[n] + 1;
- const int tmax = 0;
+ // Set the new number of active time levels
+ groupdata.at(group).info.activetimelevels = timelevels[n];
// Allocate the time levels as well
- for (int m=0; m<(int)arrdata.at(group).size(); ++m) {
- for (int var=0; var<gp.numvars; ++var) {
- const int vectorindex
- = (gp.contiguous
- ? var
- : gp.vectorgroup ? var % gp.vectorlength : 0);
- const int vectorlength
- = (gp.contiguous
- ? gp.numvars
- : gp.vectorgroup ? gp.vectorlength : 1);
- assert (vectorindex>=0 and vectorindex<gp.numvars);
- assert (vectorlength>0 and vectorlength<=gp.numvars);
- ggf* vectorleader
- = (vectorindex>0
- ? arrdata.at(group).at(m).data.at(var - vectorindex)
- : NULL);
- const int varindex = firstvarindex + var;
- switch (gp.vartype) {
+ for (int ml=min_ml; ml<max_ml; ++ml) {
+ for (int rl=min_rl; rl<max_rl; ++rl) {
+ for (int m=0; m<(int)arrdata.at(group).size(); ++m) {
+ for (int var=0; var<gp.numvars; ++var) {
+ const int vectorindex
+ = gp.vectorgroup ? var % gp.vectorlength : 0;
+ const int vectorlength
+ = gp.vectorgroup ? gp.vectorlength : 1;
+ assert (vectorindex>=0 and vectorindex<gp.numvars);
+ assert (vectorlength>0 and vectorlength<=gp.numvars);
+ ggf* const vectorleader
+ = (vectorindex>0
+ ? arrdata.at(group).at(m).data.at(var - vectorindex)
+ : NULL);
+ const int varindex = firstvarindex + var;
+#warning "TODO: allocate these in SetupGH, and after recomposing"
+ if (! arrdata.at(group).at(m).data.at(var)) {
+ switch (gp.vartype) {
#define TYPECASE(N,T) \
- case N: \
- arrdata.at(group).at(m).data.at(var) = new gf<T> \
- (varindex, groupdata.at(group).transport_operator, \
- *arrdata.at(group).at(m).tt, *arrdata.at(group).at(m).dd, \
- tmin, tmax, prolongation_order_time, \
- vectorlength, vectorindex, (gf<T>*)vectorleader); \
- break;
+ case N: \
+ arrdata.at(group).at(m).data.at(var) = new gf<T> \
+ (varindex, \
+ groupdata.at(group).transport_operator, \
+ *arrdata.at(group).at(m).tt, \
+ *arrdata.at(group).at(m).dd, \
+ prolongation_order_time, \
+ vectorlength, vectorindex, (gf<T>*)vectorleader); \
+ break;
#include "typecase"
#undef TYPECASE
- default:
- UnsupportedVarType (varindex);
- } // switch gp.vartype
-
- // Set the data pointers for grid arrays
- if (gp.grouptype != CCTK_GF) {
- assert (m == 0);
- int const c = CCTK_MyProc(cgh);
- for (int tl=0; tl<gp.numtimelevels; ++tl) {
- cgh->data[varindex][tl]
- = (tl < groupdata.at(group).activetimelevels
- ? ((*arrdata.at(group).at(m).data.at(var))
- (-tl, 0, c, 0)->storage())
- : NULL);
- }
- }
-
- } // for var
- } // for m
+ default:
+ UnsupportedVarType (varindex);
+ } // switch gp.vartype
+ } // if ! allocated
+
+ arrdata.at(group).at(m).data.at(var)->set_timelevels
+ (ml, rl, timelevels[n]);
+
+ // Set the data pointers for grid arrays
+ if (gp.grouptype != CCTK_GF) {
+ assert (rl==0 and m==0);
+ int const c = CCTK_MyProc(cgh);
+ for (int tl=0; tl<gp.numtimelevels; ++tl) {
+ cgh->data[varindex][tl]
+ = (tl < groupdata.at(group).info.activetimelevels
+ ? ((*arrdata.at(group).at(m).data.at(var))
+ (tl, 0, c, 0)->storage())
+ : NULL);
+ }
+ } // if grouptype != GF
+
+ } // for var
+ } // for m
+ } // for rl
+ } // for ml
} // if really change the number of active time levels
@@ -146,8 +168,8 @@ namespace Carpet {
if (gp.grouptype == CCTK_GF) {
if (max_refinement_levels > 1) {
if (groupdata.at(group).transport_operator != op_none) {
- if (groupdata.at(group).activetimelevels != 0
- and (groupdata.at(group).activetimelevels
+ if (groupdata.at(group).info.activetimelevels != 0
+ and (groupdata.at(group).info.activetimelevels
<= prolongation_order_time))
{
char * const groupname = CCTK_GroupName (group);
@@ -163,7 +185,7 @@ namespace Carpet {
#endif
// Record current number of time levels
- total_num_timelevels += groupdata.at(group).activetimelevels;
+ total_num_timelevels += groupdata.at(group).info.activetimelevels;
} // for n
@@ -229,11 +251,8 @@ namespace Carpet {
group = CCTK_GroupIndex(groupname);
}
assert (group>=0 and group<CCTK_NumGroups());
- const int timelevels = 0;
- int status;
- GroupStorageIncrease (cgh, 1, &group, &timelevels, &status);
// Return whether storage is allocated
- return status;
+ return groupdata.at(group).info.activetimelevels > 0;
}
diff --git a/Carpet/Carpet/src/functions.hh b/Carpet/Carpet/src/functions.hh
index e5fdf4ed2..14ba8c50d 100644
--- a/Carpet/Carpet/src/functions.hh
+++ b/Carpet/Carpet/src/functions.hh
@@ -40,7 +40,7 @@ namespace Carpet {
// Helpers for recomposing the grid hierarchy
- void CheckRegions (const gh::rexts & bbsss,
+ void CheckRegions (const gh::mexts & bbsss,
const gh::rbnds & obss,
const gh::rprocs& pss);
@@ -48,7 +48,7 @@ namespace Carpet {
void OutputGridStructure (const cGH *cgh,
const int m,
- const gh::rexts & bbsss,
+ const gh::mexts & bbsss,
const gh::rbnds & obss,
const gh::rprocs& pss);
@@ -66,9 +66,9 @@ namespace Carpet {
vector<bbvect>& obs, vector<int>& ps);
void MakeMultigridBoxes (const cGH* cgh,
- vector<ibbox> const & bbs,
- vector<bbvect> const & obs,
- vector<vector<ibbox> > & bbss);
+ gh::rexts const & bbss,
+ gh::rbnds const & obss,
+ gh::mexts & bbsss);
} // namespace Carpet
diff --git a/Carpet/Carpet/src/modes.cc b/Carpet/Carpet/src/modes.cc
index d489cf88b..10237e2df 100644
--- a/Carpet/Carpet/src/modes.cc
+++ b/Carpet/Carpet/src/modes.cc
@@ -90,9 +90,9 @@ namespace Carpet {
const int m = 0;
const int c = CCTK_MyProc(cgh);
- const ibbox& base = arrdata.at(group).at(m).hh->bases().at(rl).at(ml);
+ const ibbox& base = arrdata.at(group).at(m).hh->bases().at(ml).at(rl);
const bbvect& obnds = arrdata.at(group).at(m).hh->outer_boundaries().at(rl).at(c);
- const ibbox& ext = arrdata.at(group).at(m).dd->boxes.at(rl).at(c).at(ml).exterior;
+ const ibbox& ext = arrdata.at(group).at(m).dd->boxes.at(ml).at(rl).at(c).exterior;
ivect::ref(const_cast<int*>(groupdata.at(group).info.nghostzones))
= arrdata.at(group).at(m).dd->lghosts;
@@ -136,7 +136,7 @@ namespace Carpet {
for (int tl=0; tl<num_tl; ++tl) {
ggf * const ff = arrdata.at(group).at(m).data.at(var);
if (ff) {
- gdata * const data = (*ff) (-tl, rl, c, ml);
+ gdata * const data = (*ff) (tl, rl, c, ml);
assert (data);
cgh->data[firstvar+var][tl] = data->storage();
} else {
@@ -288,8 +288,8 @@ namespace Carpet {
carpetGH.map = map = m;
// Set grid shape
- const ibbox& coarseext = vdd.at(map)->bases.at(0 ).at(mglevel).exterior;
- const ibbox& baseext = vdd.at(map)->bases.at(reflevel).at(mglevel).exterior;
+ const ibbox& coarseext = vdd.at(map)->bases.at(mglevel).at(0 ).exterior;
+ const ibbox& baseext = vdd.at(map)->bases.at(mglevel).at(reflevel).exterior;
assert (all (baseext.lower() % baseext.stride() == 0));
assert (all ((baseext.lower() - coarseext.lower()) % baseext.stride() == 0));
ivect::ref(cgh->cctk_levoff) = (baseext.lower() - coarseext.lower()) / baseext.stride();
@@ -353,9 +353,9 @@ namespace Carpet {
component = c;
// Set cGH fields
- const ibbox& baseext = vdd.at(map)->bases.at(reflevel).at(mglevel).exterior;
+ const ibbox& baseext = vdd.at(map)->bases.at(mglevel).at(reflevel).exterior;
const bbvect& obnds = vhh.at(map)->outer_boundaries().at(reflevel).at(component);
- const ibbox& ext = vdd.at(map)->boxes.at(reflevel).at(component).at(mglevel).exterior;
+ const ibbox& ext = vdd.at(map)->boxes.at(mglevel).at(reflevel).at(component).exterior;
ivect::ref(cgh->cctk_lsh) = ext.shape() / ext.stride();
ivect::ref(cgh->cctk_lbnd)
@@ -416,7 +416,7 @@ namespace Carpet {
for (int tl=0; tl<num_tl; ++tl) {
ggf * const ff = arrdata.at(group).at(map).data.at(var);
if (ff) {
- gdata * const data = (*ff) (-tl, reflevel, component, mglevel);
+ gdata * const data = (*ff) (tl, reflevel, component, mglevel);
assert (data);
cgh->data[firstvar+var][tl] = data->storage();
} else {
diff --git a/Carpet/Carpet/src/variables.hh b/Carpet/Carpet/src/variables.hh
index 821b77ae3..ae5ade6df 100644
--- a/Carpet/Carpet/src/variables.hh
+++ b/Carpet/Carpet/src/variables.hh
@@ -115,7 +115,6 @@ namespace Carpet {
// Data for the groups
struct groupdesc {
cGroupDynamicData info;
- int activetimelevels;
operator_type transport_operator; // prolongation and restriction
};
extern vector<groupdesc> groupdata; // [group]