From 9e5dc470f4eae1ee4e888d64d65a9c87d138c499 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Fri, 27 Feb 2015 13:56:25 -0600 Subject: [PATCH] Remove deprecated function TrilinosWrappers::SparseMatrix::compress() without argument. --- include/deal.II/lac/trilinos_sparse_matrix.h | 14 ------------ source/lac/trilinos_block_sparse_matrix.cc | 2 +- source/lac/trilinos_sparse_matrix.cc | 24 ++++++++------------ source/numerics/matrix_tools.cc | 14 ++---------- 4 files changed, 13 insertions(+), 41 deletions(-) diff --git a/include/deal.II/lac/trilinos_sparse_matrix.h b/include/deal.II/lac/trilinos_sparse_matrix.h index 3427b09a15..6e4ebb7906 100644 --- a/include/deal.II/lac/trilinos_sparse_matrix.h +++ b/include/deal.II/lac/trilinos_sparse_matrix.h @@ -1083,11 +1083,6 @@ namespace TrilinosWrappers */ void compress (::dealii::VectorOperation::values operation); - /** - * @deprecated: use compress() with VectorOperation instead. - */ - void compress () DEAL_II_DEPRECATED; - /** * Set the element (i,j) to @p value. * @@ -2312,15 +2307,6 @@ namespace TrilinosWrappers - inline - void - SparseMatrix::compress () - { - compress(::dealii::VectorOperation::unknown); - } - - - inline SparseMatrix & SparseMatrix::operator = (const double d) diff --git a/source/lac/trilinos_block_sparse_matrix.cc b/source/lac/trilinos_block_sparse_matrix.cc index f06364549e..8366cee339 100644 --- a/source/lac/trilinos_block_sparse_matrix.cc +++ b/source/lac/trilinos_block_sparse_matrix.cc @@ -258,7 +258,7 @@ namespace TrilinosWrappers void BlockSparseMatrix::collect_sizes () { - compress(); + // simply forward to the (non-public) function of the base class BaseClass::collect_sizes (); } diff --git a/source/lac/trilinos_sparse_matrix.cc b/source/lac/trilinos_sparse_matrix.cc index 7a8ea4296e..abf80a79f4 100644 --- a/source/lac/trilinos_sparse_matrix.cc +++ b/source/lac/trilinos_sparse_matrix.cc @@ -126,9 +126,6 @@ namespace TrilinosWrappers return; } - // otherwise first flush Trilinos caches - matrix->compress (); - // get a representation of the present // row int ncols; @@ -365,7 +362,7 @@ namespace TrilinosWrappers { Assert(sparsity_pattern.trilinos_sparsity_pattern().Filled() == true, ExcMessage("The Trilinos sparsity pattern has not been compressed.")); - compress(); + compress(VectorOperation::insert); } @@ -605,7 +602,7 @@ namespace TrilinosWrappers // In the end, the matrix needs to be compressed in order to be really // ready. - compress(); + compress(VectorOperation::insert); } @@ -765,7 +762,7 @@ namespace TrilinosWrappers // In the end, the matrix needs to be compressed in order to be really // ready. - compress(); + compress(VectorOperation::insert); } @@ -788,7 +785,7 @@ namespace TrilinosWrappers else nonlocal_matrix.reset (); - compress(); + compress(VectorOperation::insert); last_action = Zero; } @@ -811,7 +808,7 @@ namespace TrilinosWrappers else nonlocal_matrix.reset(); - compress(); + compress(VectorOperation::unknown); } @@ -951,7 +948,7 @@ namespace TrilinosWrappers &values[0], false); } - compress(); + compress(VectorOperation::insert); } @@ -985,7 +982,7 @@ namespace TrilinosWrappers my_nonzeros*sizeof (TrilinosScalar)); } - compress(); + compress(VectorOperation::insert); } @@ -1009,7 +1006,7 @@ namespace TrilinosWrappers ((last_action == Add) && (operation!=::dealii::VectorOperation::insert)) || ((last_action == Insert) && (operation!=::dealii::VectorOperation::add)), - ExcMessage("operation and argument to compress() do not match")); + ExcMessage("Operation and argument to compress() do not match")); } // flush buffers @@ -1102,7 +1099,6 @@ namespace TrilinosWrappers SparseMatrix::clear_rows (const std::vector &rows, const TrilinosScalar new_diag_value) { - compress(); for (size_type row=0; row(&transposed_mat.trilinos_matrix()), A_,false); diff --git a/source/numerics/matrix_tools.cc b/source/numerics/matrix_tools.cc index 4b6e04253d..eeec3f0743 100644 --- a/source/numerics/matrix_tools.cc +++ b/source/numerics/matrix_tools.cc @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 1998 - 2014 by the deal.II authors +// Copyright (C) 1998 - 2015 by the deal.II authors // // This file is part of the deal.II library. // @@ -2462,8 +2462,6 @@ namespace MatrixTools const unsigned int n_blocks = matrix.n_block_rows(); - matrix.compress(); - // We need to find the subdivision // into blocks for the boundary values. // To this end, generate a vector of @@ -2558,12 +2556,6 @@ namespace MatrixTools Assert (local_range == solution.local_range(), ExcInternalError()); - // we have to read and write from this - // matrix (in this order). this will only - // work if we compress the matrix first, - // done here - matrix.compress (); - // determine the first nonzero diagonal // entry from within the part of the // matrix that we can see. if we can't @@ -2626,7 +2618,7 @@ namespace MatrixTools } // clean up - matrix.compress (); + matrix.compress (VectorOperation::insert); solution.compress (VectorOperation::insert); right_hand_side.compress (VectorOperation::insert); } @@ -2652,8 +2644,6 @@ namespace MatrixTools const unsigned int n_blocks = matrix.n_block_rows(); - matrix.compress(); - // We need to find the subdivision // into blocks for the boundary values. // To this end, generate a vector of -- 2.39.5