]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated function TrilinosWrappers::SparseMatrix::compress() without argument.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 27 Feb 2015 19:56:25 +0000 (13:56 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 6 Mar 2015 22:48:42 +0000 (16:48 -0600)
include/deal.II/lac/trilinos_sparse_matrix.h
source/lac/trilinos_block_sparse_matrix.cc
source/lac/trilinos_sparse_matrix.cc
source/numerics/matrix_tools.cc

index 3427b09a1535ed9a135ae92954234360a1c9db52..6e4ebb79065572bebe62fee42fa6c5022135a198 100644 (file)
@@ -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>i,j</i>) to @p value.
      *
@@ -2312,15 +2307,6 @@ namespace TrilinosWrappers
 
 
 
-  inline
-  void
-  SparseMatrix::compress ()
-  {
-    compress(::dealii::VectorOperation::unknown);
-  }
-
-
-
   inline
   SparseMatrix &
   SparseMatrix::operator = (const double d)
index f06364549e69c688cb88978c44990fab07deb9f0..8366cee339212f7b7f952d69a76ff9e152cd64c3 100644 (file)
@@ -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 ();
   }
 
index 7a8ea4296e9f3b646b1df9da6f92df154aa06357..abf80a79f4d92f70ffb36ecbd243d99e92f83003 100644 (file)
@@ -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<size_type> &rows,
                             const TrilinosScalar          new_diag_value)
   {
-    compress();
     for (size_type row=0; row<rows.size(); ++row)
       clear_row(rows[row], new_diag_value);
 
@@ -1111,7 +1107,7 @@ namespace TrilinosWrappers
     // data, so we need to flush the
     // buffers to make sure that the
     // right data is used.
-    compress();
+    compress(VectorOperation::insert);
   }
 
 
@@ -1441,7 +1437,7 @@ namespace TrilinosWrappers
                                     GID, values[j]);
 #endif
             }
-          transposed_mat.compress();
+          transposed_mat.compress(VectorOperation::insert);
           ML_Operator_WrapEpetraCrsMatrix
           (const_cast<Epetra_CrsMatrix *>(&transposed_mat.trilinos_matrix()),
            A_,false);
index 4b6e04253da1ac40271214181094c002311acc45..eeec3f07438a0b7048ee7818fa6a6cd23726e956 100644 (file)
@@ -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

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.