From 0e33900f3e7b79823a21d95f29cdbc3681924d45 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Tue, 30 Aug 2011 19:14:51 +0000 Subject: [PATCH] No longer try to call PETScWrappers::MatrixBase::compress() automatically, but instead only error out if the user tries to do something that won't work. git-svn-id: https://svn.dealii.org/trunk@24219 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/doc/news/changes.h | 12 ++--- .../include/deal.II/lac/petsc_matrix_base.h | 45 ++++++++++++------- deal.II/source/lac/petsc_matrix_base.cc | 2 + tests/mpi/petsc_02.cc | 17 ++++--- tests/mpi/petsc_02/ncpu_1/cmp/generic | 1 - tests/mpi/petsc_02/ncpu_10/cmp/generic | 1 - tests/mpi/petsc_02/ncpu_4/cmp/generic | 1 - 7 files changed, 45 insertions(+), 34 deletions(-) diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 2e10cc328d..5ca9275299 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -21,18 +21,20 @@ inconvenience this causes.

    -
  1. Changed: The PETScWrapper::VectorBase class tried to keep track of +
  2. Changed: The PETScWrapper::VectorBase and PETScWrapper::MatrixBase +classes tried to keep track of whether the last operation done on a vector was to add to an element or to write into one. If the previous such operation was of a different kind than the current one, we would flush buffers (see the description in @ref GlossCompress). However, trying to do this automatically turned out to be an endless source of hard-to-find bugs in %parallel programs. -The scheme has therefore now been changed to the following: the class -keeps track of the previous operation and if it differs from the +The scheme has therefore now been changed to the following: the classes +keep track of the previous operation and if it differs from the current one, reports an error stating that the user needs to call -PETScWrapper::VectorBase::compress() instead. +PETScWrapper::VectorBase::compress() or +PETScWrapper::MatrixBase::compress() instead.
    -(Wolfgang Bangerth, 2011/08/03) +(Wolfgang Bangerth, 2011/08/03, 2011/08/30)
  3. Changed: The classes Tensor, SymmetricTensor and Point now have an additional template argument for the number type. While a default template diff --git a/deal.II/include/deal.II/lac/petsc_matrix_base.h b/deal.II/include/deal.II/lac/petsc_matrix_base.h index 2d55c9dc83..2cb99594b0 100644 --- a/deal.II/include/deal.II/lac/petsc_matrix_base.h +++ b/deal.II/include/deal.II/lac/petsc_matrix_base.h @@ -1148,7 +1148,24 @@ namespace PETScWrappers * Exception */ DeclException0 (ExcSourceEqualsDestination); - + + /** + * Exception. + */ + DeclException2 (ExcWrongMode, + int, int, + << "You tried to do a " + << (arg1 == 1 ? + "'set'" : + (arg1 == 2 ? + "'add'" : "???")) + << " operation but the vector is currently in " + << (arg2 == 1 ? + "'set'" : + (arg2 == 2 ? + "'add'" : "???")) + << " mode. You first have to call 'compress()'."); + protected: /** * A generic matrix object in @@ -1190,9 +1207,10 @@ namespace PETScWrappers LastAction::Values last_action; /** - * Flush buffers on all CPUs when - * switching between inserting and - * adding to elements, no-op otherwise. + * Ensure that the add/set mode that + * is required for actions following + * this call is compatible with the + * current mode. * Should be called from all internal * functions accessing matrix elements. */ @@ -1433,7 +1451,6 @@ namespace PETScWrappers const unsigned int j, const PetscScalar value) { - Assert (numbers::is_finite(value), ExcNumberNotFinite()); set (i, 1, &j, &value, false); @@ -1793,18 +1810,12 @@ namespace PETScWrappers void MatrixBase::prepare_action(const LastAction::Values new_action) { - // flush PETSc buffers when switching - // actions, otherwise just return. - if (last_action == new_action) return; - - int ierr; - ierr = MatAssemblyBegin(matrix, MAT_FLUSH_ASSEMBLY); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - ierr = MatAssemblyEnd(matrix, MAT_FLUSH_ASSEMBLY); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - last_action = new_action; + if (last_action == new_action) + ; + else if (last_action == LastAction::none) + last_action = new_action; + else + Assert (false, ExcWrongMode (last_action, new_action)); } diff --git a/deal.II/source/lac/petsc_matrix_base.cc b/deal.II/source/lac/petsc_matrix_base.cc index b8195a4338..e07373baf2 100644 --- a/deal.II/source/lac/petsc_matrix_base.cc +++ b/deal.II/source/lac/petsc_matrix_base.cc @@ -247,6 +247,8 @@ namespace PETScWrappers ierr = MatAssemblyEnd (matrix,MAT_FINAL_ASSEMBLY); AssertThrow (ierr == 0, ExcPETScError(ierr)); + + last_action == LastAction::none; } diff --git a/tests/mpi/petsc_02.cc b/tests/mpi/petsc_02.cc index 401600b3ec..52ceeda37a 100644 --- a/tests/mpi/petsc_02.cc +++ b/tests/mpi/petsc_02.cc @@ -2,7 +2,7 @@ // $Id: simple_mpi_01.cc 23327 2011-02-11 03:19:07Z bangerth $ // Version: $Name$ // -// Copyright (C) 2009 by the deal.II authors +// Copyright (C) 2009, 2011 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -19,8 +19,8 @@ #include "../tests.h" #include -#include -#include +#include +#include #include #include @@ -39,18 +39,17 @@ void test() CompressedSimpleSparsityPattern csp(2); csp.add(0,0); csp.add(1,1); - + PETScWrappers::MPI::SparseMatrix mat; std::vector< unsigned int > local_rows(numprocs, 0); local_rows[0]=2; - + mat.reinit(MPI_COMM_WORLD, csp, local_rows, local_rows, myid); -#warning this here hangs, so it is uncommented until bug is fixed - //if (myid == 0 ) - mat.add(0, 0, 1.0); + if (myid == 0 ) + mat.add(0, 0, 1.0); + - mat.compress(); if (myid==0) diff --git a/tests/mpi/petsc_02/ncpu_1/cmp/generic b/tests/mpi/petsc_02/ncpu_1/cmp/generic index 76a3add384..721906ea79 100644 --- a/tests/mpi/petsc_02/ncpu_1/cmp/generic +++ b/tests/mpi/petsc_02/ncpu_1/cmp/generic @@ -1,4 +1,3 @@ DEAL::Running on 1 CPU(s). -DEAL::2nd try DEAL::done diff --git a/tests/mpi/petsc_02/ncpu_10/cmp/generic b/tests/mpi/petsc_02/ncpu_10/cmp/generic index 729abc0fba..7430da6964 100644 --- a/tests/mpi/petsc_02/ncpu_10/cmp/generic +++ b/tests/mpi/petsc_02/ncpu_10/cmp/generic @@ -1,4 +1,3 @@ DEAL::Running on 10 CPU(s). -DEAL::2nd try DEAL::done diff --git a/tests/mpi/petsc_02/ncpu_4/cmp/generic b/tests/mpi/petsc_02/ncpu_4/cmp/generic index 25e315071a..befd397045 100644 --- a/tests/mpi/petsc_02/ncpu_4/cmp/generic +++ b/tests/mpi/petsc_02/ncpu_4/cmp/generic @@ -1,4 +1,3 @@ DEAL::Running on 4 CPU(s). -DEAL::2nd try DEAL::done -- 2.39.5