From: Lei Qiao Date: Thu, 17 Sep 2015 22:18:13 +0000 (-0500) Subject: use same counting function for master and slave node in parallel::distributed::GridRe... X-Git-Tag: v8.4.0-rc2~391^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=d69f5cac382fa5963182a4c1c9852c5acc80ea0d;p=dealii.git use same counting function for master and slave node in parallel::distributed::GridRefinement::refine_and_coarsen_fixed_fraction() --- diff --git a/source/distributed/grid_refinement.cc b/source/distributed/grid_refinement.cc index 8377382b5f..7325bc9e0e 100644 --- a/source/distributed/grid_refinement.cc +++ b/source/distributed/grid_refinement.cc @@ -352,22 +352,23 @@ namespace */ template number - master_compute_threshold (const Vector &criteria, - const std::pair global_min_and_max, - const double target_error, - MPI_Comm mpi_communicator) + compute_threshold (const Vector &criteria, + const std::pair global_min_and_max, + const double target_error, + MPI_Comm mpi_communicator) { double interesting_range[2] = { global_min_and_max.first, global_min_and_max.second }; adjust_interesting_range (interesting_range); + const unsigned int master_mpi_rank = 0; unsigned int iteration = 0; do { MPI_Bcast (&interesting_range[0], 2, MPI_DOUBLE, - 0, mpi_communicator); + master_mpi_rank, mpi_communicator); if (interesting_range[0] == interesting_range[1]) { @@ -382,7 +383,7 @@ namespace double final_threshold = std::min (interesting_range[0], global_min_and_max.second); MPI_Bcast (&final_threshold, 1, MPI_DOUBLE, - 0, mpi_communicator); + master_mpi_rank, mpi_communicator); return final_threshold; } @@ -390,8 +391,7 @@ namespace const double test_threshold = (interesting_range[0] > 0 ? - std::exp((std::log(interesting_range[0]) + - std::log(interesting_range[1])) / 2) + std::sqrt(interesting_range[0] * interesting_range[1]) : (interesting_range[0] + interesting_range[1]) / 2); @@ -405,12 +405,15 @@ namespace double total_error; MPI_Reduce (&my_error, &total_error, 1, MPI_DOUBLE, - MPI_SUM, 0, mpi_communicator); + MPI_SUM, master_mpi_rank, mpi_communicator); // now adjust the range. if we have to many cells, we take // the upper half of the previous range, otherwise the lower // half. if we have hit the right number, then set the range - // to the exact value + // to the exact value. + // slave nodes also update their own interesting_range, however + // their results are not significant since the values will be + // overwritten by MPI_Bcast from the master node in next loop. if (total_error > target_error) interesting_range[0] = test_threshold; else if (total_error < target_error) @@ -437,60 +440,6 @@ namespace Assert (false, ExcInternalError()); return -1; } - - - /** - * The corresponding function to - * the one above, to be run on the - * slaves. - */ - template - number - slave_compute_threshold (const Vector &criteria, - MPI_Comm mpi_communicator) - { - do - { - double interesting_range[2] = { -1, -1 }; - MPI_Bcast (&interesting_range[0], 2, MPI_DOUBLE, - 0, mpi_communicator); - - if (interesting_range[0] == interesting_range[1]) - { - // we got the final value. the master adjusts the - // value one more time, so receive it from there - double final_threshold = -1; - MPI_Bcast (&final_threshold, 1, MPI_DOUBLE, - 0, mpi_communicator); - - return final_threshold; - } - - // count how many elements - // there are that are bigger - // than the following trial - // threshold - const double test_threshold - = (interesting_range[0] > 0 - ? - std::exp((std::log(interesting_range[0]) + - std::log(interesting_range[1])) / 2) - : - (interesting_range[0] + interesting_range[1]) / 2); - - double my_error = 0; - for (unsigned int i=0; i test_threshold) - my_error += criteria(i); - - MPI_Reduce (&my_error, 0, 1, MPI_DOUBLE, - MPI_SUM, 0, mpi_communicator); - } - while (true); - - Assert (false, ExcInternalError()); - return -1; - } } } @@ -633,72 +582,33 @@ namespace parallel const double total_error = compute_global_sum (locally_owned_indicators, mpi_communicator); - - // from here on designate a - // master and slaves double top_threshold, bottom_threshold; - if (Utilities::MPI::this_mpi_process (mpi_communicator) == 0) - { - // this is the master - // processor - top_threshold - = - RefineAndCoarsenFixedFraction:: - master_compute_threshold (locally_owned_indicators, - global_min_and_max, - top_fraction_of_error * - total_error, - mpi_communicator); - // compute bottom - // threshold only if - // necessary. otherwise - // use a threshold lower - // than the smallest - // value we have locally - if (bottom_fraction_of_error > 0) - bottom_threshold - = - RefineAndCoarsenFixedFraction:: - master_compute_threshold (locally_owned_indicators, - global_min_and_max, - (1-bottom_fraction_of_error) * - total_error, - mpi_communicator); - else - { - bottom_threshold = *std::min_element (criteria.begin(), - criteria.end()); - bottom_threshold -= std::fabs(bottom_threshold); - } - } + top_threshold = + RefineAndCoarsenFixedFraction:: + compute_threshold (locally_owned_indicators, + global_min_and_max, + top_fraction_of_error * + total_error, + mpi_communicator); + // compute bottom + // threshold only if + // necessary. otherwise + // use a threshold lower + // than the smallest + // value we have locally + if (bottom_fraction_of_error > 0) + bottom_threshold = + RefineAndCoarsenFixedFraction:: + compute_threshold (locally_owned_indicators, + global_min_and_max, + (1-bottom_fraction_of_error) * + total_error, + mpi_communicator); else { - // this is a slave - // processor - top_threshold - = - RefineAndCoarsenFixedFraction:: - slave_compute_threshold (locally_owned_indicators, - mpi_communicator); - - // compute bottom - // threshold only if - // necessary. otherwise - // use a threshold lower - // than the smallest - // value we have locally - if (bottom_fraction_of_error > 0) - bottom_threshold - = - RefineAndCoarsenFixedFraction:: - slave_compute_threshold (locally_owned_indicators, - mpi_communicator); - else - { - bottom_threshold = *std::min_element (criteria.begin(), - criteria.end()); - bottom_threshold -= std::fabs(bottom_threshold); - } + bottom_threshold = *std::min_element (criteria.begin(), + criteria.end()); + bottom_threshold -= std::fabs(bottom_threshold); } // now refine the mesh