From 8456be3581aca4b976f214ce69da31f77fdd0475 Mon Sep 17 00:00:00 2001 From: bangerth Date: Tue, 6 Aug 2013 21:44:57 +0000 Subject: [PATCH] Fix bug #78: BlockMatrixBase::set/add had race conditions when used from multiple threads because these functions used scratch arrays. One just needed to add a mutex to guard these data structures. git-svn-id: https://svn.dealii.org/trunk@30238 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/doc/news/changes.h | 8 + .../include/deal.II/lac/block_matrix_base.h | 158 ++++++++++++------ tests/lac/block_matrices_02.cc | 131 +++++++++++++++ tests/lac/block_matrices_02/cmp/generic | 25 +++ tests/lac/block_matrices_03.cc | 134 +++++++++++++++ tests/lac/block_matrices_03/cmp/generic | 25 +++ 6 files changed, 427 insertions(+), 54 deletions(-) create mode 100644 tests/lac/block_matrices_02.cc create mode 100644 tests/lac/block_matrices_02/cmp/generic create mode 100644 tests/lac/block_matrices_03.cc create mode 100644 tests/lac/block_matrices_03/cmp/generic diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index be39ececbd..e87472eac8 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -44,6 +44,14 @@ inconvenience this causes.

Specific improvements

    +
  1. + Fixed: The various block matrix classes are all derived from BlockMatrixBase + which had race conditions when the set() or add() functions were called from + different threads. This is now fixed. +
    + (Wolfgang Bangerth, 2013/08/05) +
  2. +
  3. Fixed: various fixes with assignment and reinit of PETScWrappers::MPI::Vector.
    diff --git a/deal.II/include/deal.II/lac/block_matrix_base.h b/deal.II/include/deal.II/lac/block_matrix_base.h index 273ad266be..6ed20d407e 100644 --- a/deal.II/include/deal.II/lac/block_matrix_base.h +++ b/deal.II/include/deal.II/lac/block_matrix_base.h @@ -20,6 +20,7 @@ #include #include +#include #include #include #include @@ -1246,30 +1247,72 @@ protected: private: - /** - * Temporary vector for counting the - * elements written into the - * individual blocks when doing a - * collective add or set. - */ - std::vector counter_within_block; /** - * Temporary vector for column - * indices on each block when writing - * local to global data on each - * sparse matrix. + * A structure containing some fields used by the + * set() and add() functions that is used to pre-sort + * the input fields. Since one can reasonably expect + * to call set() and add() from multiple threads at once + * as long as the matrix indices that are touched are + * disjoint, these temporary data fields need to be + * guarded by a mutex; the structure therefore contains such + * a mutex as a member variable. */ - std::vector > column_indices; + struct TemporaryData + { + /** + * Temporary vector for counting the + * elements written into the + * individual blocks when doing a + * collective add or set. + */ + std::vector counter_within_block; + + /** + * Temporary vector for column + * indices on each block when writing + * local to global data on each + * sparse matrix. + */ + std::vector > column_indices; + + /** + * Temporary vector for storing the + * local values (they need to be + * reordered when writing local to + * global). + */ + std::vector > column_values; + + /** + * A mutex variable used to guard access to the member + * variables of this structure; + */ + Threads::Mutex mutex; + + /** + * Copy operator. This is needed because the default copy + * operator of this class is deleted (since Threads::Mutex is + * not copyable) and hence the default copy operator of the + * enclosing class is also deleted. + * + * The implementation here simply does nothing -- TemporaryData + * objects are just scratch objects that are resized at the + * beginning of their use, so there is no point actually copying + * anything. + */ + TemporaryData & operator = (const TemporaryData &) + {} + }; /** - * Temporary vector for storing the - * local values (they need to be - * reordered when writing local to - * global). + * A set of scratch arrays that can be used by the add() + * and set() functions that take pointers to data to + * pre-sort indices before use. Access from multiple threads + * is synchronized via the mutex variable that is part of the + * structure. */ - std::vector > column_values; - + TemporaryData temporary_data; /** * Make the iterator class a @@ -1772,9 +1815,10 @@ BlockMatrixBase::memory_consumption () const MemoryConsumption::memory_consumption(row_block_indices)+ MemoryConsumption::memory_consumption(column_block_indices)+ MemoryConsumption::memory_consumption(sub_objects)+ - MemoryConsumption::memory_consumption(counter_within_block)+ - MemoryConsumption::memory_consumption(column_indices)+ - MemoryConsumption::memory_consumption(column_values); + MemoryConsumption::memory_consumption(temporary_data.counter_within_block)+ + MemoryConsumption::memory_consumption(temporary_data.column_indices)+ + MemoryConsumption::memory_consumption(temporary_data.column_values)+ + MemoryConsumption::memory_consumption(temporary_data.mutex); for (unsigned int r=0; r::set (const size_type row, { prepare_set_operation(); + // lock access to the temporary data structure to + // allow multiple threads to call this function concurrently + Threads::Mutex::ScopedLock lock (temporary_data.mutex); + // Resize scratch arrays - if (column_indices.size() < this->n_block_cols()) + if (temporary_data.column_indices.size() < this->n_block_cols()) { - column_indices.resize (this->n_block_cols()); - column_values.resize (this->n_block_cols()); - counter_within_block.resize (this->n_block_cols()); + temporary_data.column_indices.resize (this->n_block_cols()); + temporary_data.column_values.resize (this->n_block_cols()); + temporary_data.counter_within_block.resize (this->n_block_cols()); } // Resize sub-arrays to n_cols. This @@ -1994,19 +2042,19 @@ BlockMatrixBase::set (const size_type row, // whether the size of one is large // enough before actually going // through all of them. - if (column_indices[0].size() < n_cols) + if (temporary_data.column_indices[0].size() < n_cols) { for (unsigned int i=0; in_block_cols(); ++i) { - column_indices[i].resize(n_cols); - column_values[i].resize(n_cols); + temporary_data.column_indices[i].resize(n_cols); + temporary_data.column_values[i].resize(n_cols); } } // Reset the number of added elements // in each block to zero. for (unsigned int i=0; in_block_cols(); ++i) - counter_within_block[i] = 0; + temporary_data.counter_within_block[i] = 0; // Go through the column indices to // find out which portions of the @@ -2027,10 +2075,10 @@ BlockMatrixBase::set (const size_type row, const std::pair col_index = this->column_block_indices.global_to_local(col_indices[j]); - const size_type local_index = counter_within_block[col_index.first]++; + const size_type local_index = temporary_data.counter_within_block[col_index.first]++; - column_indices[col_index.first][local_index] = col_index.second; - column_values[col_index.first][local_index] = value; + temporary_data.column_indices[col_index.first][local_index] = col_index.second; + temporary_data.column_values[col_index.first][local_index] = value; } #ifdef DEBUG @@ -2038,7 +2086,7 @@ BlockMatrixBase::set (const size_type row, // the right length has been obtained. size_type length = 0; for (unsigned int i=0; in_block_cols(); ++i) - length += counter_within_block[i]; + length += temporary_data.counter_within_block[i]; Assert (length <= n_cols, ExcInternalError()); #endif @@ -2051,14 +2099,14 @@ BlockMatrixBase::set (const size_type row, row_index = this->row_block_indices.global_to_local (row); for (unsigned int block_col=0; block_col::add (const size_type row, return; } - // Resize scratch arrays - if (column_indices.size() < this->n_block_cols()) + // Lock scratch arrays, then resize them + Threads::Mutex::ScopedLock lock (temporary_data.mutex); + + if (temporary_data.column_indices.size() < this->n_block_cols()) { - column_indices.resize (this->n_block_cols()); - column_values.resize (this->n_block_cols()); - counter_within_block.resize (this->n_block_cols()); + temporary_data.column_indices.resize (this->n_block_cols()); + temporary_data.column_values.resize (this->n_block_cols()); + temporary_data.counter_within_block.resize (this->n_block_cols()); } // Resize sub-arrays to n_cols. This @@ -2239,19 +2289,19 @@ BlockMatrixBase::add (const size_type row, // whether the size of one is large // enough before actually going // through all of them. - if (column_indices[0].size() < n_cols) + if (temporary_data.column_indices[0].size() < n_cols) { for (unsigned int i=0; in_block_cols(); ++i) { - column_indices[i].resize(n_cols); - column_values[i].resize(n_cols); + temporary_data.column_indices[i].resize(n_cols); + temporary_data.column_values[i].resize(n_cols); } } // Reset the number of added elements // in each block to zero. for (unsigned int i=0; in_block_cols(); ++i) - counter_within_block[i] = 0; + temporary_data.counter_within_block[i] = 0; // Go through the column indices to // find out which portions of the @@ -2272,10 +2322,10 @@ BlockMatrixBase::add (const size_type row, const std::pair col_index = this->column_block_indices.global_to_local(col_indices[j]); - const size_type local_index = counter_within_block[col_index.first]++; + const size_type local_index = temporary_data.counter_within_block[col_index.first]++; - column_indices[col_index.first][local_index] = col_index.second; - column_values[col_index.first][local_index] = value; + temporary_data.column_indices[col_index.first][local_index] = col_index.second; + temporary_data.column_values[col_index.first][local_index] = value; } #ifdef DEBUG @@ -2283,7 +2333,7 @@ BlockMatrixBase::add (const size_type row, // the right length has been obtained. size_type length = 0; for (unsigned int i=0; in_block_cols(); ++i) - length += counter_within_block[i]; + length += temporary_data.counter_within_block[i]; Assert (length <= n_cols, ExcInternalError()); #endif @@ -2296,14 +2346,14 @@ BlockMatrixBase::add (const size_type row, row_index = this->row_block_indices.global_to_local (row); for (unsigned int block_col=0; block_col +#include +#include +#include +#include +#include +#include + + + +void do_set (const bool even_or_odd, + BlockSparseMatrix &bsm) +{ + BlockSparseMatrix::size_type col_indices[5]; + for (unsigned int i=0; i<5 ; ++i) + if (even_or_odd) + col_indices[i] = 2*i; + else + col_indices[i] = 2*i+1; + + BlockSparseMatrix::value_type values[5]; + for (unsigned int i=0; i<5 ; ++i) + if (even_or_odd) + values[i] = 1; + else + values[i] = 2; + + bsm.set (0, 5, col_indices, values, false); +} + + +void test () +{ + std::ofstream logfile("block_matrices_02/output"); + deallog << std::fixed; + deallog << std::setprecision(2); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + BlockSparsityPattern bsp(2,2); + // set sizes + for (unsigned int i=0; i<2; ++i) + for (unsigned int j=0; j<2; ++j) + bsp.block(i,j).reinit ( 5, 5, 5); + bsp.collect_sizes (); + + // make a full matrix + for (unsigned int row=0; row<10; ++row) + for (unsigned int i=0; i<10; ++i) + bsp.add (row, i); + bsp.compress (); + + BlockSparseMatrix bsm (bsp); + + Threads::ThreadGroup<> tg; + for (unsigned int i=0; i<100; ++i) + { + tg += Threads::new_thread (&do_set, true, bsm); + tg += Threads::new_thread (&do_set, false, bsm); + } + tg.join_all (); + + bsm.print_formatted (deallog.get_file_stream()); +} + + + + +int main () +{ + try + { + test (); + } + catch (std::exception &e) + { + std::cerr << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + std::cerr << "Exception on processing: " << e.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + // abort + return 2; + } + catch (...) + { + std::cerr << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + std::cerr << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + // abort + return 3; + }; + + + return 0; +} diff --git a/tests/lac/block_matrices_02/cmp/generic b/tests/lac/block_matrices_02/cmp/generic new file mode 100644 index 0000000000..3dca993a6d --- /dev/null +++ b/tests/lac/block_matrices_02/cmp/generic @@ -0,0 +1,25 @@ + +Component (0,0) +1.000e+00 2.000e+00 1.000e+00 2.000e+00 1.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (0,1) +2.000e+00 1.000e+00 2.000e+00 1.000e+00 2.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (1,0) +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (1,1) +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 diff --git a/tests/lac/block_matrices_03.cc b/tests/lac/block_matrices_03.cc new file mode 100644 index 0000000000..16d97674e6 --- /dev/null +++ b/tests/lac/block_matrices_03.cc @@ -0,0 +1,134 @@ +// --------------------------------------------------------------------- +// $Id$ +// +// Copyright (C) 2000 - 2013 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// the versions of BlockMatrixBase::set/add that take a pointer to a set of +// values were not thread safe, even if they were writing into disjoint sets +// of elements. this test verifies that the add() function is now indeed +// thread safe under these conditions +// +// we test this by calling add() from a multitude of threads at once. without +// the patch that fixed this, one can elicit a great deal of different memory +// corruption errors and assertions from this test, depending on luck and +// phase of the moon :-) + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include + + + +void do_add (const bool even_or_odd, + BlockSparseMatrix &bsm) +{ + BlockSparseMatrix::size_type col_indices[5]; + for (unsigned int i=0; i<5 ; ++i) + if (even_or_odd) + col_indices[i] = 2*i; + else + col_indices[i] = 2*i+1; + + BlockSparseMatrix::value_type values[5]; + for (unsigned int i=0; i<5 ; ++i) + if (even_or_odd) + values[i] = 1; + else + values[i] = 2; + + bsm.add (0, 5, col_indices, values, false); +} + + +void test () +{ + std::ofstream logfile("block_matrices_03/output"); + deallog << std::fixed; + deallog << std::setprecision(2); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + BlockSparsityPattern bsp(2,2); + // set sizes + for (unsigned int i=0; i<2; ++i) + for (unsigned int j=0; j<2; ++j) + bsp.block(i,j).reinit ( 5, 5, 5); + bsp.collect_sizes (); + + // make a full matrix + for (unsigned int row=0; row<10; ++row) + for (unsigned int i=0; i<10; ++i) + bsp.add (row, i); + bsp.compress (); + + BlockSparseMatrix bsm (bsp); + + Threads::ThreadGroup<> tg; + for (unsigned int i=0; i<100; ++i) + { + tg += Threads::new_thread (&do_add, true, bsm); + tg += Threads::new_thread (&do_add, false, bsm); + } + tg.join_all (); + + // divide whole matrix by 100 to get back to the effect of a single set + bsm /= 100; + + bsm.print_formatted (deallog.get_file_stream()); +} + + + + +int main () +{ + try + { + test (); + } + catch (std::exception &e) + { + std::cerr << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + std::cerr << "Exception on processing: " << e.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + // abort + return 2; + } + catch (...) + { + std::cerr << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + std::cerr << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + // abort + return 3; + }; + + + return 0; +} diff --git a/tests/lac/block_matrices_03/cmp/generic b/tests/lac/block_matrices_03/cmp/generic new file mode 100644 index 0000000000..3dca993a6d --- /dev/null +++ b/tests/lac/block_matrices_03/cmp/generic @@ -0,0 +1,25 @@ + +Component (0,0) +1.000e+00 2.000e+00 1.000e+00 2.000e+00 1.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (0,1) +2.000e+00 1.000e+00 2.000e+00 1.000e+00 2.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (1,0) +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +Component (1,1) +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 +0.000e+00 0.000e+00 0.000e+00 0.000e+00 0.000e+00 -- 2.39.5