diff options
author | Erik Schnetter <schnetter@cct.lsu.edu> | 2008-02-19 05:24:00 +0000 |
---|---|---|
committer | Erik Schnetter <schnetter@cct.lsu.edu> | 2008-02-19 05:24:00 +0000 |
commit | 679ccfcacb45d91bc12f34241708350714906c19 (patch) | |
tree | 9d62a1ad9dc99130dc5b489acd63d9c5a04fff48 /Carpet/CarpetInterp | |
parent | 014e48022202a124db2288f1c21c1c438987af40 (diff) |
CarpetInterp: Add parameter "barriers" to catch calling errors
Add a new parameter CarpetInterp::barriers. If set, CarpetInterp
introduces MPI_Barrier statements. These statements catch calling
errors where not all processors call CarpetInterp at the same time.
darcs-hash:20080219052452-dae7b-49900543c4360e82f94fe5b5d50c66d594d77ad3.gz
Diffstat (limited to 'Carpet/CarpetInterp')
-rw-r--r-- | Carpet/CarpetInterp/param.ccl | 4 | ||||
-rw-r--r-- | Carpet/CarpetInterp/src/interp.cc | 15 |
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) { |