aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorknarf <knarf@b6729ddc-ac74-4bd1-908c-9dc7244c52a1>2010-01-12 17:18:27 +0000
committerknarf <knarf@b6729ddc-ac74-4bd1-908c-9dc7244c52a1>2010-01-12 17:18:27 +0000
commit5e56e8927d0e6772996b6073e03b04965d2496b4 (patch)
tree02bd293052794bcf803e590202bb24184dd714fa
parente32ca78bb29e9f7b8188a8233305cd9dab60a4a6 (diff)
fix bug
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/Hydro_Analysis/trunk@33 b6729ddc-ac74-4bd1-908c-9dc7244c52a1
-rw-r--r--src/Hydro_Analysis.c43
1 files changed, 43 insertions, 0 deletions
diff --git a/src/Hydro_Analysis.c b/src/Hydro_Analysis.c
index 7d182ea..3f1e3f2 100644
--- a/src/Hydro_Analysis.c
+++ b/src/Hydro_Analysis.c
@@ -77,3 +77,46 @@ void Hydro_Analysis_LocationSearch(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS
DECLARE_CCTK_PARAMETERS
+
+ /* Initialized MPI-local quantities */
+ CCTK_REAL local_rho_max_loc[3];
+ local_rho_max_loc[0] = 0.0;
+ local_rho_max_loc[1] = 0.0;
+ local_rho_max_loc[2] = 0.0;
+
+ /* Look for the location of the global maximum.
+ * This algorithm will have problems when that occurs at more than
+ * one location. */
+ #pragma omp parallel
+ {
+ LC_LOOP3(location_reduction, i,j,k,
+ cctk_nghostzones[0], cctk_nghostzones[1], cctk_nghostzones[2],
+ cctk_lsh[0] - cctk_nghostzones[0], cctk_lsh[1] - cctk_nghostzones[1],
+ cctk_lsh[2] - cctk_nghostzones[2], cctk_lsh[0], cctk_lsh[1], cctk_lsh[2])
+ {
+ CCTK_INT i3D = CCTK_GFINDEX3D(cctkGH, i, j, k);
+ if (rho[i3D] == *Hydro_Analysis_rho_max &&
+ (!Hydro_Analysis_rho_max_loc_only_positive_x || x[i3D] >= 0.0))
+ {
+ #pragma omp critical
+ {
+ local_rho_max_loc[0] = x[i3D];
+ local_rho_max_loc[1] = y[i3D];
+ local_rho_max_loc[2] = z[i3D];
+ }
+ }
+ }
+ LC_ENDLOOP3(location_reduction);
+ }
+ /* After this at least one MPI process should have the location set.
+ * The others should have a value of 0.0, so that the sum reduction
+ * gives the value on each MPI process. This however is not true sometimes
+ * and this is the point where the algorithm can fail. */
+ CCTK_INT handle_sum = CCTK_LocalArrayReductionHandle("sum");
+ if (handle_sum < 0)
+ CCTK_WARN(0, "Unable to get reduction handle 'sum'.");
+ if (CCTK_ReduceLocArrayToArray1D(cctkGH, -1, handle_sum, &local_rho_max_loc,
+ Hydro_Analysis_rho_max_loc, 3, CCTK_VARIABLE_REAL))
+ CCTK_WARN(0, "Error while reducing local_rho_max_loc");
+}
+