aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorRoland Haas <rhaas@tapir.caltech.edu>2013-09-14 14:53:58 -0700
committerRoland Haas <rhaas@tapir.caltech.edu>2013-09-27 15:10:18 -0700
commit911ab59919a935b1d25146830673b2d200c74fed (patch)
tree39716a6f8e8dae7eade077a95500c6875627c407
parentd8cc74e37df4ff3750ccaed11b67e152cc0295f3 (diff)
Carpet: add new group tag no_split_directions
-rw-r--r--Carpet/Carpet/src/Recompose.cc10
-rw-r--r--Carpet/Carpet/src/SetupGH.cc11
-rw-r--r--Carpet/Carpet/src/functions.hh6
3 files changed, 20 insertions, 7 deletions
diff --git a/Carpet/Carpet/src/Recompose.cc b/Carpet/Carpet/src/Recompose.cc
index 9654ecefd..6f18fb657 100644
--- a/Carpet/Carpet/src/Recompose.cc
+++ b/Carpet/Carpet/src/Recompose.cc
@@ -1165,12 +1165,13 @@ namespace Carpet {
void
SplitRegions_Automatic (cGH const * const cctkGH,
vector<region_t> & superregs,
- vector<region_t> & regs)
+ vector<region_t> & regs,
+ bvect const & no_split_dims)
{
assert (regs.empty());
vector<vector<region_t> > superregss (1, superregs);
vector<vector<region_t> > regss (1);
- SplitRegionsMaps_Automatic (cctkGH, superregss, regss);
+ SplitRegionsMaps_Automatic (cctkGH, superregss, regss, no_split_dims);
assert (superregss.size() == 1);
superregs = superregss.AT(0);
assert (regss.size() == 1);
@@ -1595,7 +1596,8 @@ namespace Carpet {
void
SplitRegionsMaps_Automatic (cGH const * const cctkGH,
vector<vector<region_t> > & superregss,
- vector<vector<region_t> > & regss)
+ vector<vector<region_t> > & regss,
+ bvect const & no_split_dims)
{
DECLARE_CCTK_PARAMETERS;
@@ -1710,7 +1712,7 @@ namespace Carpet {
for (int r=0, p=0; r<nregs; p+=mynprocs.AT(r), ++r) {
if (recompose_verbose) cout << "SRMA superreg[" << r << "] " << superregs.AT(r) << endl;
// bvect const dims = false;
- bvect const dims = not split_components;
+ bvect const dims = not split_components || no_split_dims;
SplitRegionsMaps_Automatic_Recursively
(dims, p, mynprocs.AT(r), superregs.AT(r), newregs);
} // for r
diff --git a/Carpet/Carpet/src/SetupGH.cc b/Carpet/Carpet/src/SetupGH.cc
index f24d3f6d4..80bb952a5 100644
--- a/Carpet/Carpet/src/SetupGH.cc
+++ b/Carpet/Carpet/src/SetupGH.cc
@@ -1249,7 +1249,16 @@ namespace Carpet {
// Split it into components, one for each processor
switch (gdata.disttype) {
case CCTK_DISTRIB_DEFAULT: {
- SplitRegions_Automatic (cctkGH, superregs, regs);
+ CCTK_INT no_split_directions[dim];
+ bvect no_split_dims(false);
+ int const nvals =
+ Util_TableGetIntArray(gdata.tagstable, dim, no_split_directions,
+ "no_split_directions");
+ assert((0 <= nvals && nvals <= dim) || nvals == UTIL_ERROR_TABLE_NO_SUCH_KEY);
+ for(int i=0;i<nvals;++i) {
+ no_split_dims[i] = true;
+ }
+ SplitRegions_Automatic (cctkGH, superregs, regs, no_split_dims);
break;
}
case CCTK_DISTRIB_CONSTANT: {
diff --git a/Carpet/Carpet/src/functions.hh b/Carpet/Carpet/src/functions.hh
index 5a4c4fea1..c9eabbe75 100644
--- a/Carpet/Carpet/src/functions.hh
+++ b/Carpet/Carpet/src/functions.hh
@@ -181,7 +181,8 @@ namespace Carpet {
void
SplitRegions_Automatic (cGH const * cctkGH,
vector<region_t> & superregs,
- vector<region_t> & regs);
+ vector<region_t> & regs,
+ bvect const & no_split_dims = false);
void
SplitRegionsMaps_Recursively (cGH const * cctkGH,
vector<vector<region_t> > & superregss,
@@ -194,7 +195,8 @@ namespace Carpet {
void
SplitRegionsMaps_Automatic (cGH const * cctkGH,
vector<vector<region_t> > & superregss,
- vector<vector<region_t> > & regss);
+ vector<vector<region_t> > & regss,
+ bvect const & no_split_dims = false);
void
SplitRegionsMaps_Recursively (cGH const * cctkGH,
vector<vector<region_t> > & superregss,