aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@gmail.com>2012-07-03 21:33:41 -0400
committerErik Schnetter <schnetter@gmail.com>2012-07-03 21:33:41 -0400
commit6679e32674b71ed65cc2dbbb97e114f1699fee11 (patch)
treec55ff175168495cfa9d6b45e4eb172fef6b89493
parentd826d9a2ffdbf30b8babe8ba5e39074fc1a659da (diff)
Introduce cctk_ash, retire cctk_lssh
Introduce cctk_ash, describing the process-local array shape that has been allocated. This may be larger than cctk_lsh, the process-local shape that should be used. Retire cctk_lssh and related infrastructure to handle staggered grid functions.
-rw-r--r--Carpet/Carpet/src/SetupGH.cc13
-rw-r--r--Carpet/Carpet/src/modes.cc97
-rw-r--r--Carpet/CarpetInterp/src/interp.cc5
-rw-r--r--Carpet/CarpetLib/src/restrict_3d_cc_o5_rf2.cc29
-rw-r--r--Carpet/CarpetRegrid2/src/amr.cc6
5 files changed, 41 insertions, 109 deletions
diff --git a/Carpet/Carpet/src/SetupGH.cc b/Carpet/Carpet/src/SetupGH.cc
index 59ae79550..2e814628f 100644
--- a/Carpet/Carpet/src/SetupGH.cc
+++ b/Carpet/Carpet/src/SetupGH.cc
@@ -1306,9 +1306,7 @@ namespace Carpet {
groupdata.AT(group).info.dim = gdata.dim;
groupdata.AT(group).info.gsh = new int [dim];
groupdata.AT(group).info.lsh = new int [dim];
-#ifdef CCTK_GROUPDYNAMICDATA_HAS_LSSH
- groupdata.AT(group).info.lssh = new int [CCTK_NSTAGGER*dim];
-#endif
+ groupdata.AT(group).info.ash = new int [dim];
groupdata.AT(group).info.lbnd = new int [dim];
groupdata.AT(group).info.ubnd = new int [dim];
groupdata.AT(group).info.bbox = new int [2*dim];
@@ -2737,15 +2735,6 @@ namespace Carpet {
free (groupname);
}
#endif
-
- // Staggered groups are not supported
- if (gdata.stagtype != 0) {
- char * const groupname = CCTK_GroupName (group);
- CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
- "The group \"%s\" is staggered. Staggered groups are not yet supported",
- groupname);
- free (groupname);
- }
}
diff --git a/Carpet/Carpet/src/modes.cc b/Carpet/Carpet/src/modes.cc
index b5c9a9fb6..387d3a3ff 100644
--- a/Carpet/Carpet/src/modes.cc
+++ b/Carpet/Carpet/src/modes.cc
@@ -128,30 +128,17 @@ namespace Carpet {
= (ext.lower() - baseext.lower()) / ext.stride();
ivect_ref(info.ubnd)
= (ext.upper() - baseext.lower()) / ext.stride();
- ivect_ref(info.lsh) = gdata::allocated_memory_shape (ext);
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- for (int d=0; d<dim; ++d) {
- // TODO: support staggering
- const_cast<int*>(info.lssh)[CCTK_LSSH_IDX(stg,d)]
- = (ext.shape() / ext.stride())[d];
- }
- }
-#endif
+ ivect_ref(info.lsh) = ext.shape() / ext.stride();
+ ivect_ref(info.ash) = gdata::allocated_memory_shape (ext);
if (gp.disttype == CCTK_DISTRIB_CONSTANT) {
int const dir = gp.dim==0 ? 0 : gp.dim-1;
ivect & gsh = ivect_ref(info.gsh);
- ivect lssh = ivect_ref(info.lsh);
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- for (int d=0; d<dim; ++d) {
- lssh[d] = info.lssh[CCTK_LSSH_IDX(0,d)];
- }
-#endif
ivect & lbnd = ivect_ref(info.lbnd);
ivect & ubnd = ivect_ref(info.ubnd);
- gsh[dir] = lssh[dir];
+ ivect const & lsh = ivect_ref(info.lsh);
+ gsh[dir] = lsh[dir];
lbnd[dir] = 0;
- ubnd[dir] = lssh[dir] - 1;
+ ubnd[dir] = lsh[dir] - 1;
}
for (int d=0; d<dim; ++d) {
const_cast<int*>(info.bbox)[2*d ] = obnds[0][d];
@@ -162,22 +149,13 @@ namespace Carpet {
for (int d=0; d<dim; ++d) {
assert (info.lsh[d]>=0);
- //assert (info.lsh[d]<=info.gsh[d]);
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- assert (info.lssh[CCTK_LSSH_IDX(0,d)]>=0);
- assert (info.lssh[CCTK_LSSH_IDX(0,d)]<=info.gsh[d]);
- }
-#endif
+ assert (info.lsh[d]<=info.gsh[d]);
assert (info.lbnd[d]>=0);
assert (info.lbnd[d]<=info.ubnd[d]+1);
assert (info.ubnd[d]<info.gsh[d]);
- //assert (info.lbnd[d] + info.lsh[d] - 1 == info.ubnd[d]);
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- assert (info.lbnd[d] + info.lssh[CCTK_LSSH_IDX(0,d)] - 1
- == info.ubnd[d]);
-#endif
+ assert (info.lbnd[d] + info.lsh[d] - 1 == info.ubnd[d]);
assert (info.lbnd[d]<=info.ubnd[d]+1);
+ assert (info.lsh[d]<=info.ash[d]);
}
const int numvars = CCTK_NumVarsInGroupI (group);
@@ -256,13 +234,7 @@ namespace Carpet {
const_cast<int*>(info.bbox)[2*d ] = deadbeef;
const_cast<int*>(info.bbox)[2*d+1] = deadbeef;
}
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- for (int d=0; d<dim; ++d) {
- const_cast<int*>(info.lssh)[CCTK_LSSH_IDX(stg,d)] = deadbeef;
- }
- }
-#endif
+ ivect_ref(info.ash) = deadbeef;
info.activetimelevels = deadbeef;
const int numvars = CCTK_NumVarsInGroupI (group);
@@ -551,37 +523,28 @@ namespace Carpet {
(ext.lower() - baseext.lower()) / ext.stride();
ivect_ref(cctkGH->cctk_ubnd) =
(ext.upper() - baseext.lower()) / ext.stride();
- ivect_ref(cctkGH->cctk_lsh) = gdata::allocated_memory_shape (ext);
+ ivect_ref(cctkGH->cctk_lsh) = ext.shape() / ext.stride();
+ ivect_ref(cctkGH->cctk_ash) = gdata::allocated_memory_shape (ext);
for (int d=0; d<dim; ++d) {
cctkGH->cctk_bbox[2*d ] = obnds[0][d];
cctkGH->cctk_bbox[2*d+1] = obnds[1][d];
}
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- for (int d=0; d<dim; ++d) {
- // TODO: support staggering
- cctkGH->cctk_lssh[CCTK_LSSH_IDX(stg,d)]
- = (ext.shape() / ext.stride())[d];
- }
- }
for (int d=0; d<dim; ++d) {
cctkGH->cctk_from[d] = 0;
- cctkGH->cctk_to [d] = cctkGH->cctk_lssh[CCTK_LSSH_IDX(0,d)];
+ cctkGH->cctk_to [d] = cctkGH->cctk_lsh[d];
}
for (int d=0; d<dim; ++d) {
- //assert (cctkGH->cctk_lsh[d] >= 0);
- //assert (cctkGH->cctk_lsh[d] <= cctkGH->cctk_gsh[d]);
- assert (cctkGH->cctk_lssh[CCTK_LSSH_IDX(0,d)] >= 0);
- assert (cctkGH->cctk_lssh[CCTK_LSSH_IDX(0,d)] <= cctkGH->cctk_gsh[d]);
+ assert (cctkGH->cctk_lsh[d] >= 0);
+ assert (cctkGH->cctk_lsh[d] <= cctkGH->cctk_gsh[d]);
assert (cctkGH->cctk_lbnd[d] >= 0);
assert (cctkGH->cctk_lbnd[d] <= cctkGH->cctk_ubnd[d] + 1);
assert (cctkGH->cctk_ubnd[d] < cctkGH->cctk_gsh[d]);
- //assert (cctkGH->cctk_lbnd[d] + cctkGH->cctk_lsh[d] - 1 == cctkGH->cctk_ubnd[d]);
- assert (cctkGH->cctk_lbnd[d] + cctkGH->cctk_lssh[CCTK_LSSH_IDX(0,d)] - 1 == cctkGH->cctk_ubnd[d]);
+ assert (cctkGH->cctk_lbnd[d] + cctkGH->cctk_lsh[d] - 1 == cctkGH->cctk_ubnd[d]);
assert (cctkGH->cctk_lbnd[d] <= cctkGH->cctk_ubnd[d]+1);
- assert (cctkGH->cctk_lssh[CCTK_LSSH_IDX(0,d)] <= cctkGH->cctk_lsh[d]);
+ assert (cctkGH->cctk_lsh[d] <= cctkGH->cctk_ash[d]);
assert (cctkGH->cctk_from[d] >= 0);
assert (cctkGH->cctk_from[d] <= cctkGH->cctk_to[d]);
assert (cctkGH->cctk_to[d] <= cctkGH->cctk_lsh[d]);
@@ -594,20 +557,13 @@ namespace Carpet {
ivect_ref(info.lbnd) = ivect_ref(cctkGH->cctk_lbnd);
ivect_ref(info.ubnd) = ivect_ref(cctkGH->cctk_ubnd);
-
ivect_ref(info.lsh) = ivect_ref(cctkGH->cctk_lsh);
+ ivect_ref(info.ash) = ivect_ref(cctkGH->cctk_ash);
for (int d=0; d<dim; ++d) {
const_cast<int*>(info.bbox)[2*d ] = cctkGH->cctk_bbox[2*d ];
const_cast<int*>(info.bbox)[2*d+1] = cctkGH->cctk_bbox[2*d+1];
}
-
- // for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- // for (int d=0; d<dim; ++d) {
- // const_cast<int*>(info.lssh)[CCTK_LSSH_IDX(stg,d)]
- // = cctkGH->cctk_lssh[CCTK_LSSH_IDX(stg,d)];
- // }
- // }
if (local_component != -1) {
const int numvars = CCTK_NumVarsInGroupI (group);
@@ -690,19 +646,13 @@ namespace Carpet {
ivect_ref(cctkGH->cctk_from) = -deadbeef;
ivect_ref(cctkGH->cctk_to ) = deadbeef;
ivect_ref(cctkGH->cctk_lsh ) = deadbeef;
+ ivect_ref(cctkGH->cctk_ash ) = deadbeef;
for (int d=0; d<dim; ++d) {
cctkGH->cctk_bbox[2*d ] = deadbeef;
cctkGH->cctk_bbox[2*d+1] = deadbeef;
}
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- for (int d=0; d<dim; ++d) {
- // TODO: support staggering
- cctkGH->cctk_lssh[CCTK_LSSH_IDX(stg,d)] = cctkGH->cctk_lsh[d];
- }
- }
-
for (int group=0; group<CCTK_NumGroups(); ++group) {
if (CCTK_GroupTypeI(group) == CCTK_GF) {
@@ -710,23 +660,14 @@ namespace Carpet {
ivect_ref(info.lbnd) = ivect_ref(cctkGH->cctk_lbnd);
ivect_ref(info.ubnd) = ivect_ref(cctkGH->cctk_ubnd);
-
ivect_ref(info.lsh) = ivect_ref(cctkGH->cctk_lsh);
+ ivect_ref(info.ash) = ivect_ref(cctkGH->cctk_ash);
for (int d=0; d<dim; ++d) {
const_cast<int*>(info.bbox)[2*d ] = cctkGH->cctk_bbox[2*d ];
const_cast<int*>(info.bbox)[2*d+1] = cctkGH->cctk_bbox[2*d+1];
}
-#ifdef CCTK_HAVE_CGROUPDYNAMICDATA_LSSH
- for (int stg=0; stg<CCTK_NSTAGGER; ++stg) {
- for (int d=0; d<dim; ++d) {
- const_cast<int*>(info.lssh)[CCTK_LSSH_IDX(stg,d)]
- = cctkGH->cctk_lssh[CCTK_LSSH_IDX(stg,d)];
- }
- }
-#endif
-
if (local_component != -1) {
const int numvars = CCTK_NumVarsInGroupI (group);
if (numvars>0) {
diff --git a/Carpet/CarpetInterp/src/interp.cc b/Carpet/CarpetInterp/src/interp.cc
index 77e4fb4cb..3418a932c 100644
--- a/Carpet/CarpetInterp/src/interp.cc
+++ b/Carpet/CarpetInterp/src/interp.cc
@@ -250,11 +250,6 @@ namespace CarpetInterp {
"input array variable %d is not of type CCTK_GF or "
"CCTK_ARRY", n);
}
- if (CCTK_GroupStaggerIndexGI (group) != 0) {
- CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
- "interpolation of staggered input array variable %d "
- "is not supported", n);
- }
if (coord_group < 0) {
coord_group = group;
CCTK_GroupDynamicData (cctkGH, coord_group, &coord_group_data);
diff --git a/Carpet/CarpetLib/src/restrict_3d_cc_o5_rf2.cc b/Carpet/CarpetLib/src/restrict_3d_cc_o5_rf2.cc
index c9067c9d9..a90ebb868 100644
--- a/Carpet/CarpetLib/src/restrict_3d_cc_o5_rf2.cc
+++ b/Carpet/CarpetLib/src/restrict_3d_cc_o5_rf2.cc
@@ -24,11 +24,12 @@ namespace CarpetLib {
index3 (dstioff + (i), dstjoff + (j), dstkoff + (k), \
dstiext, dstjext, dstkext)
#define SRCOFF3(i,j,k) \
- offset3 (srcioff + (i), srcjoff + (j), srckoff + (k), \
- srciext, srcjext, srckext)
+ offset3 (srcioff + (i), srcjoff + (j), srckoff + (k), \
+ srciext, srcjext, srckext)
#define DSTOFF3(i,j,k) \
- offset3 (dstioff + (i), dstjoff + (j), dstkoff + (k), \
- dstiext, dstjext, dstkext)
+ offset3 (dstioff + (i), dstjoff + (j), dstkoff + (k), \
+ dstiext, dstjext, dstkext)
+
// This operator offers fifth-order accurate restriction operators for cell
@@ -43,10 +44,12 @@ namespace CarpetLib {
// TODO: use prolongate_3d_rf2 instead of writing new operators, its the same
// interpolation.
+ namespace {
+
// 1D restriction
template <typename T>
- static inline T restrict1 (T const * restrict const p,
- size_t const d1)
+ inline T restrict1 (T const * restrict const p,
+ size_t const d1)
{
typedef typename typeprops<T>::real RT;
RT const den = 256;
@@ -63,8 +66,8 @@ namespace CarpetLib {
// 2D restriction
template <typename T>
- static inline T restrict2 (T const * restrict const p,
- size_t const d1, size_t const d2)
+ inline T restrict2 (T const * restrict const p,
+ size_t const d1, size_t const d2)
{
typedef typename typeprops<T>::real RT;
RT const den = 256;
@@ -81,8 +84,8 @@ namespace CarpetLib {
// 3D restriction
template <typename T>
- static inline T restrict3 (T const * restrict const p,
- size_t const d1, size_t const d2, size_t const d3)
+ inline T restrict3 (T const * restrict const p,
+ size_t const d1, size_t const d2, size_t const d3)
{
typedef typename typeprops<T>::real RT;
RT const den = 256;
@@ -97,6 +100,10 @@ namespace CarpetLib {
return res;
}
+ } // namespace
+
+
+
template <typename T>
void
restrict_3d_cc_o5_rf2 (T const * restrict const src,
@@ -221,7 +228,7 @@ namespace CarpetLib {
}
#endif
dst [DSTIND3(i, j, k)] =
- restrict3
+ CarpetLib::restrict3
(& src[SRCIND3(2*i, 2*j, 2*k)], srcdi, srcdj, srcdk);
}
diff --git a/Carpet/CarpetRegrid2/src/amr.cc b/Carpet/CarpetRegrid2/src/amr.cc
index 6d23dd9c0..26f78717c 100644
--- a/Carpet/CarpetRegrid2/src/amr.cc
+++ b/Carpet/CarpetRegrid2/src/amr.cc
@@ -102,12 +102,12 @@ namespace CarpetRegrid2 {
}
ivect const lbnd = ivect::ref(cctk_lbnd);
- ivect const lssh (CCTK_LSSH(0,0), CCTK_LSSH(0,1), CCTK_LSSH(0,2));
+ ivect const lsh = ivect::ref(cctk_lbnd);
ivect const bboxlo (cctk_bbox[0], cctk_bbox[2], cctk_bbox[4]);
ivect const bboxhi (cctk_bbox[1], cctk_bbox[3], cctk_bbox[5]);
- ivect const imin = 0 + either(bboxlo, 0, nghostzones);
- ivect const imax = lssh - either(bboxhi, 0, nghostzones);
+ ivect const imin = 0 + either(bboxlo, 0, nghostzones);
+ ivect const imax = lsh - either(bboxhi, 0, nghostzones);
ivect const bmin =
max (0,