]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated function MGTools::apply_boundary_values.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 25 Jan 2015 18:20:39 +0000 (12:20 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 25 Jan 2015 18:39:50 +0000 (12:39 -0600)
doc/news/changes.h
include/deal.II/multigrid/mg_tools.h
source/multigrid/mg_tools.cc

index a38e73378176481f40647e62f8820abfea3dcc97..977ecd63dfe2ce0b7ae3e36ab4fd855494aff5b5 100644 (file)
@@ -204,6 +204,7 @@ inconvenience this causes.
   - GridTools::extract_boundary_mesh.
   - PreconditionBlock::size.
   - MGTools::count_dofs_per_component.
+  - MGTools::apply_boundary_values.
   <br>
   This release also removes the deprecated class MGDoFHandler. The
   functionality of this class had previously been incorporated into
index ce4e0a522d009c6843b0b97ff2502b3ff4795f7d..85277d5bf2f36d1dfd9237ebac4cbdd4f53aceea 100644 (file)
@@ -203,25 +203,6 @@ namespace MGTools
                       std::vector<IndexSet>                 &boundary_indices,
                       const ComponentMask               &component_mask = ComponentMask());
 
-  /**
-   * @deprecated
-   */
-  template <typename number>
-  void
-  apply_boundary_values (const std::set<types::global_dof_index> &boundary_dofs,
-                         SparseMatrix<number> &matrix,
-                         const bool preserve_symmetry,
-                         const bool ignore_zeros = false) DEAL_II_DEPRECATED;
-
-  /**
-   * @deprecated
-   */
-  template <typename number>
-  void
-  apply_boundary_values (const std::set<types::global_dof_index> &boundary_dofs,
-                         BlockSparseMatrix<number> &matrix,
-                         const bool preserve_symmetry) DEAL_II_DEPRECATED;
-
   /**
    * For each level in a multigrid hierarchy, produce an IndexSet that
    * indicates which of the degrees of freedom are along interfaces of this
index d23cf4a6d342e584d03f75aa55e99a8e6490bb6e..68b1a2963cefaf4b811587ec8537b8bc90e85b4d 100644 (file)
@@ -1576,89 +1576,10 @@ namespace MGTools
       }
 
   }
-
-
-  template <typename number>
-  void
-  apply_boundary_values (
-    const std::set<types::global_dof_index> &boundary_dofs,
-    SparseMatrix<number> &matrix,
-    const bool preserve_symmetry,
-    const bool /*ignore_zeros*/)
-  {
-    // this function is not documented and not tested in the testsuite
-    // so it isn't quite clear what it's supposed to do. it also isn't
-    // used anywhere else in the library. in avoiding the use of
-    // a deprecated function, I therefore threw away the original function
-    // and replaced it by the following, which I believe should work
-
-    std::map<types::global_dof_index, double> boundary_values;
-    for (std::set<types::global_dof_index>::const_iterator p=boundary_dofs.begin();
-         p != boundary_dofs.end(); ++p)
-      boundary_values[*p] = 0;
-
-    Vector<number> dummy(matrix.m());
-    MatrixTools::apply_boundary_values (boundary_values,
-                                        matrix, dummy, dummy,
-                                        preserve_symmetry);
-  }
-
-
-
-  template <typename number>
-  void
-  apply_boundary_values (
-    const std::set<types::global_dof_index> &boundary_dofs,
-    BlockSparseMatrix<number> &matrix,
-    const bool preserve_symmetry)
-  {
-    Assert (matrix.n_block_rows() == matrix.n_block_cols(),
-            ExcNotQuadratic());
-    Assert (matrix.get_sparsity_pattern().get_row_indices() ==
-            matrix.get_sparsity_pattern().get_column_indices(),
-            ExcNotQuadratic());
-
-    // this function is not documented and not tested in the testsuite
-    // so it isn't quite clear what it's supposed to do. it also isn't
-    // used anywhere else in the library. in avoiding the use of
-    // a deprecated function, I therefore threw away the original function
-    // and replaced it by the following, which I believe should work
-
-    std::map<types::global_dof_index, double> boundary_values;
-    for (std::set<types::global_dof_index>::const_iterator p=boundary_dofs.begin();
-         p != boundary_dofs.end(); ++p)
-      boundary_values[*p] = 0;
-
-    BlockVector<number> dummy(matrix.n_block_rows());
-    for (unsigned int i=0; i<matrix.n_block_rows(); ++i)
-      dummy.block(i).reinit (matrix.block(i,i).m());
-    dummy.collect_sizes();
-
-    MatrixTools::apply_boundary_values (boundary_values,
-                                        matrix, dummy, dummy,
-                                        preserve_symmetry);
-  }
 }
 
 
 // explicit instantiations
 #include "mg_tools.inst"
 
-namespace MGTools
-{
-  template void apply_boundary_values (
-    const std::set<types::global_dof_index> &,
-    SparseMatrix<float> &, const bool, const bool);
-  template void apply_boundary_values (
-    const std::set<types::global_dof_index> &,
-    SparseMatrix<double> &, const bool, const bool);
-  template void apply_boundary_values (
-    const std::set<types::global_dof_index> &,
-    BlockSparseMatrix<float> &, const bool);
-  template void apply_boundary_values (
-    const std::set<types::global_dof_index> &,
-    BlockSparseMatrix<double> &, const bool);
-}
-
-
 DEAL_II_NAMESPACE_CLOSE

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.