aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetInterp
diff options
context:
space:
mode:
Diffstat (limited to 'Carpet/CarpetInterp')
-rw-r--r--Carpet/CarpetInterp/param.ccl4
-rw-r--r--Carpet/CarpetInterp/src/interp.cc15
2 files changed, 19 insertions, 0 deletions
diff --git a/Carpet/CarpetInterp/param.ccl b/Carpet/CarpetInterp/param.ccl
index 37f4d791a..3f6b86379 100644
--- a/Carpet/CarpetInterp/param.ccl
+++ b/Carpet/CarpetInterp/param.ccl
@@ -1,5 +1,9 @@
# Parameter definitions for thorn CarpetInterp
+BOOLEAN barriers "Insert barriers at strategic places for debugging purposes (slows down execution)" STEERABLE=always
+{
+} "no"
+
SHARES: Cactus
USES CCTK_REAL cctk_initial_time
diff --git a/Carpet/CarpetInterp/src/interp.cc b/Carpet/CarpetInterp/src/interp.cc
index a2b0f29a0..49e5d1754 100644
--- a/Carpet/CarpetInterp/src/interp.cc
+++ b/Carpet/CarpetInterp/src/interp.cc
@@ -200,6 +200,8 @@ namespace CarpetInterp {
CCTK_INT const output_array_type_codes [],
CCTK_POINTER const output_arrays [])
{
+ DECLARE_CCTK_PARAMETERS;
+
cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
assert (cctkGH);
assert (0 <= N_dims and N_dims <= dim);
@@ -309,6 +311,19 @@ namespace CarpetInterp {
}
+ if (barriers) {
+ static unsigned int magic = 0x2ad072fbUL; // a random starting value
+ unsigned int recv = magic;
+ MPI_Bcast (& recv, 1, MPI_UNSIGNED, 0, dist::comm());
+ if (recv != magic) {
+ CCTK_WARN (CCTK_WARN_ABORT,
+ "Inconsistent communication schedule: not all processes call CarpetInterp");
+ }
+ ++ magic;
+ MPI_Barrier (dist::comm());
+ }
+
+
// Find range of refinement levels
assert (maps > 0);
for (int m=1; m<maps; ++m) {