From: Wolfgang Bangerth Date: Mon, 19 May 2014 13:00:52 +0000 (+0000) Subject: Patch by Michal Wichrowksi: run compress() also if there is nothing to do for one... X-Git-Tag: v8.2.0-rc1~478 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=13b573746c73c87efa1294cda30cc5eaca1a5ef5;p=dealii.git Patch by Michal Wichrowksi: run compress() also if there is nothing to do for one processor. git-svn-id: https://svn.dealii.org/trunk@32936 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 180f5f88d7..b4fc63b2a1 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -148,16 +148,23 @@ inconvenience this causes.

Specific improvements

    -
  1. New: AnyData::try_read() is a function that allows users to check -whether an entry exists and get a pointer to it without throwing an -exception in case of failure. -
    -(Guido Kanschat, 2014/05/16) -
  2. - -
  3. New: The GMRES solver of deal.II can now write an estimate of +
  4. Fixed: The MatrixTools::apply_boundary_values() variant that works + on PETSc matrices could produce a deadlock in parallel if one processor + had no boundary values to apply. This is now fixed. +
    + (Michal Wichrowski, 2014/05/19) +
  5. + +
  6. New: AnyData::try_read() is a function that allows users to check + whether an entry exists and get a pointer to it without throwing an + exception in case of failure. +
    + (Guido Kanschat, 2014/05/16) +
  7. + +
  8. New: The GMRES solver of deal.II can now write an estimate of eigenvalues to the log file, in analogy to the CG solver. This is enabled - by the flag SolverGMRES<>::AdditionalData::compute_eigenvalues. + by the flag SolverGMRES::AdditionalData::compute_eigenvalues.
    (Martin Kronbichler, 2014/05/11)
  9. diff --git a/deal.II/source/numerics/matrix_tools.cc b/deal.II/source/numerics/matrix_tools.cc index 6a726729c1..063602a74e 100644 --- a/deal.II/source/numerics/matrix_tools.cc +++ b/deal.II/source/numerics/matrix_tools.cc @@ -2334,74 +2334,75 @@ namespace MatrixTools Assert (matrix.n() == solution.size(), ExcDimensionMismatch(matrix.n(), solution.size())); - // if no boundary values are to be applied - // simply return - if (boundary_values.size() == 0) - return; - - const std::pair local_range - = matrix.local_range(); - Assert (local_range == right_hand_side.local_range(), - ExcInternalError()); - Assert (local_range == solution.local_range(), - ExcInternalError()); - - // determine the first nonzero diagonal - // entry from within the part of the - // matrix that we can see. if we can't - // find such an entry, take one - PetscScalar average_nonzero_diagonal_entry = 1; - for (types::global_dof_index i=local_range.first; i constrained_rows; - for (std::map::const_iterator - dof = boundary_values.begin(); - dof != boundary_values.end(); - ++dof) - if ((dof->first >= local_range.first) && - (dof->first < local_range.second)) - constrained_rows.push_back (dof->first); + // if no boundary values are to be applied, then + // jump straight to the compress() calls that we still have + // to perform because they are collective operations + if (boundary_values.size() > 0) + { + const std::pair local_range + = matrix.local_range(); + Assert (local_range == right_hand_side.local_range(), + ExcInternalError()); + Assert (local_range == solution.local_range(), + ExcInternalError()); + + // determine the first nonzero diagonal + // entry from within the part of the + // matrix that we can see. if we can't + // find such an entry, take one + PetscScalar average_nonzero_diagonal_entry = 1; + for (types::global_dof_index i=local_range.first; i constrained_rows; + for (std::map::const_iterator + dof = boundary_values.begin(); + dof != boundary_values.end(); + ++dof) + if ((dof->first >= local_range.first) && + (dof->first < local_range.second)) + constrained_rows.push_back (dof->first); - std::vector indices; - std::vector solution_values; - for (std::map::const_iterator - dof = boundary_values.begin(); - dof != boundary_values.end(); - ++dof) - if ((dof->first >= local_range.first) && - (dof->first < local_range.second)) - { - indices.push_back (dof->first); - solution_values.push_back (dof->second); - } - solution.set (indices, solution_values); + // then eliminate these rows and set + // their diagonal entry to what we have + // determined above. note that for petsc + // matrices interleaving read with write + // operations is very expensive. thus, we + // here always replace the diagonal + // element, rather than first checking + // whether it is nonzero and in that case + // preserving it. this is different from + // the case of deal.II sparse matrices + // treated in the other functions. + matrix.clear_rows (constrained_rows, average_nonzero_diagonal_entry); + + std::vector indices; + std::vector solution_values; + for (std::map::const_iterator + dof = boundary_values.begin(); + dof != boundary_values.end(); + ++dof) + if ((dof->first >= local_range.first) && + (dof->first < local_range.second)) + { + indices.push_back (dof->first); + solution_values.push_back (dof->second); + } + solution.set (indices, solution_values); - // now also set appropriate values for - // the rhs - for (unsigned int i=0; i