]> https://gitweb.dealii.org/ - dealii.git/commitdiff
RCM: Update multigrid
authorMatthias Maier <tamiko@43-1.org>
Sat, 26 May 2018 00:13:38 +0000 (19:13 -0500)
committerMatthias Maier <tamiko@43-1.org>
Wed, 6 Jun 2018 15:19:38 +0000 (10:19 -0500)
include/deal.II/multigrid/mg_constrained_dofs.h
include/deal.II/multigrid/mg_transfer.h
include/deal.II/multigrid/mg_transfer_block.h
include/deal.II/multigrid/mg_transfer_component.h
source/multigrid/mg_transfer_prebuilt.cc
source/multigrid/multigrid.cc

index 22c643563916cc247e2ea99b8305a72165e13f06..6aefdebd38e5b13af865f03a8008054157aa5244 100644 (file)
@@ -20,6 +20,8 @@
 
 #include <deal.II/base/subscriptor.h>
 
+#include <deal.II/lac/affine_constraints.h>
+
 #include <deal.II/multigrid/mg_tools.h>
 
 #include <set>
@@ -53,14 +55,14 @@ public:
    * freedom at the boundary of the domain untouched assuming that no
    * Dirichlet boundary conditions for them exist.
    *
-   * Furthermore, this call sets up a ConstraintMatrix on each level that
-   * contains possible periodicity constraints in case those have been added to
-   * the underlying triangulation. The constraint matrix can be queried by
-   * get_level_constraint_matrix(level). Note that the current implementation of
-   * periodicity constraints in this class does not support rotation matrices in
-   * the periodicity definition, i.e., the respective argument in the
-   * GridTools::collect_periodic_faces() may not be different from the identity
-   * matrix.
+   * Furthermore, this call sets up an AffineConstraints object on each
+   * level that contains possible periodicity constraints in case those
+   * have been added to the underlying triangulation. The constraint matrix
+   * can be queried by get_level_constraint_matrix(level). Note that the
+   * current implementation of periodicity constraints in this class does
+   * not support rotation matrices in the periodicity definition, i.e., the
+   * respective argument in the GridTools::collect_periodic_faces() may not
+   * be different from the identity matrix.
    */
   template <int dim, int spacedim>
   void
@@ -160,7 +162,7 @@ public:
    * Return the level constraint matrix for a given level, containing
    * periodicity constraints (if enabled on the triangulation).
    */
-  const ConstraintMatrix &
+  const AffineConstraints<double> &
   get_level_constraint_matrix(const unsigned int level) const;
 
 private:
@@ -179,7 +181,7 @@ private:
    * Constraint matrices containing information regarding potential
    * periodic boundary conditions for each level .
    */
-  std::vector<ConstraintMatrix> level_constraints;
+  std::vector<AffineConstraints<double>> level_constraints;
 };
 
 
@@ -376,7 +378,7 @@ MGConstrainedDoFs::have_boundary_indices() const
 
 
 
-inline const ConstraintMatrix &
+inline const AffineConstraints<double> &
 MGConstrainedDoFs::get_level_constraint_matrix(const unsigned int level) const
 {
   AssertIndexRange(level, level_constraints.size());
index 2dbcfee7a20936461c74904f76b89b8c0dc07ff8..dcd6a2e621c84f66cf47b9965a85e8e0bfbd4234 100644 (file)
@@ -22,9 +22,9 @@
 
 #include <deal.II/dofs/dof_handler.h>
 
+#include <deal.II/lac/affine_constraints.h>
 #include <deal.II/lac/block_sparsity_pattern.h>
 #include <deal.II/lac/block_vector.h>
-#include <deal.II/lac/constraint_matrix.h>
 #include <deal.II/lac/la_parallel_vector.h>
 #include <deal.II/lac/petsc_parallel_sparse_matrix.h>
 #include <deal.II/lac/petsc_parallel_vector.h>
@@ -646,8 +646,8 @@ public:
    * @deprecated @p constraints is unused.
    */
   DEAL_II_DEPRECATED
-  MGTransferPrebuilt(const ConstraintMatrix & constraints,
-                     const MGConstrainedDoFs &mg_constrained_dofs);
+  MGTransferPrebuilt(const AffineConstraints<double> &constraints,
+                     const MGConstrainedDoFs &        mg_constrained_dofs);
 
   /**
    * Destructor.
@@ -667,8 +667,8 @@ public:
    */
   DEAL_II_DEPRECATED
   void
-  initialize_constraints(const ConstraintMatrix & constraints,
-                         const MGConstrainedDoFs &mg_constrained_dofs);
+  initialize_constraints(const AffineConstraints<double> &constraints,
+                         const MGConstrainedDoFs &        mg_constrained_dofs);
 
   /**
    * Reset the object to the state it had right after the default constructor.
index 3ec899a49b1fb57087de0b2dff0e24f769a2d02c..7854cd721daa37f0ffed6d60d9a990b38b8e6a3d 100644 (file)
@@ -77,8 +77,8 @@ public:
    * @deprecated @p constraints is unused.
    */
   DEAL_II_DEPRECATED
-  MGTransferBlockBase(const ConstraintMatrix & constraints,
-                      const MGConstrainedDoFs &mg_constrained_dofs);
+  MGTransferBlockBase(const AffineConstraints<double> &constraints,
+                      const MGConstrainedDoFs &        mg_constrained_dofs);
 
   /**
    * Memory used by this object.
@@ -329,8 +329,8 @@ public:
    * @deprecated @p constraints is unused.
    */
   DEAL_II_DEPRECATED
-  MGTransferBlockSelect(const ConstraintMatrix & constraints,
-                        const MGConstrainedDoFs &mg_constrained_dofs);
+  MGTransferBlockSelect(const AffineConstraints<double> &constraints,
+                        const MGConstrainedDoFs &        mg_constrained_dofs);
 
   /**
    * Destructor.
index 8ce0e025fcf619ea81300f595e950bd97b2ee835..72904d736cc11fb0ebde64ba72e11c221e55472a 100644 (file)
@@ -186,7 +186,7 @@ public:
   /**
    * Constructor with constraint matrices.
    */
-  MGTransferSelect(const ConstraintMatrix &constraints);
+  MGTransferSelect(const AffineConstraints<double> &constraints);
 
   /**
    * Destructor.
@@ -373,7 +373,7 @@ private:
    * The constraints of the global system.
    */
 public:
-  SmartPointer<const ConstraintMatrix> constraints;
+  SmartPointer<const AffineConstraints<double>> constraints;
 };
 
 /*@}*/
index 1fdc0ed80315316f084abbfae58b84c38691dc4a..45e69bb56e0da94f1b95a117a3f8c823669d011d 100644 (file)
@@ -56,7 +56,7 @@ MGTransferPrebuilt<VectorType>::MGTransferPrebuilt(
 
 template <typename VectorType>
 MGTransferPrebuilt<VectorType>::MGTransferPrebuilt(
-  const ConstraintMatrix & /*c*/,
+  const AffineConstraints<double> & /*c*/,
   const MGConstrainedDoFs &mg_c)
 {
   this->mg_constrained_dofs = &mg_c;
@@ -77,7 +77,7 @@ MGTransferPrebuilt<VectorType>::initialize_constraints(
 template <typename VectorType>
 void
 MGTransferPrebuilt<VectorType>::initialize_constraints(
-  const ConstraintMatrix & /*c*/,
+  const AffineConstraints<double> & /*c*/,
   const MGConstrainedDoFs &mg_c)
 {
   initialize_constraints(mg_c);
index e19da1f61dcc301ba8bceb94ebaf642a824082e0..c0c8105e17d3cd61d6e50002ad678f188db23636 100644 (file)
@@ -48,8 +48,9 @@ MGTransferBlockBase::MGTransferBlockBase(const MGConstrainedDoFs &mg_c) :
 
 
 
-MGTransferBlockBase::MGTransferBlockBase(const ConstraintMatrix & /*c*/,
-                                         const MGConstrainedDoFs &mg_c) :
+MGTransferBlockBase::MGTransferBlockBase(
+  const AffineConstraints<double> & /*c*/,
+  const MGConstrainedDoFs &mg_c) :
   n_mg_blocks(0),
   mg_constrained_dofs(&mg_c)
 {}
@@ -200,7 +201,7 @@ MGTransferSelect<number>::MGTransferSelect() :
 
 
 template <typename number>
-MGTransferSelect<number>::MGTransferSelect(const ConstraintMatrix &c) :
+MGTransferSelect<number>::MGTransferSelect(const AffineConstraints<double> &c) :
   selected_component(0),
   mg_selected_component(0),
   constraints(&c)
@@ -260,7 +261,7 @@ MGTransferBlockSelect<number>::MGTransferBlockSelect(
 
 template <typename number>
 MGTransferBlockSelect<number>::MGTransferBlockSelect(
-  const ConstraintMatrix & /*c*/,
+  const AffineConstraints<double> & /*c*/,
   const MGConstrainedDoFs &mg_c) :
   MGTransferBlockBase(mg_c),
   selected_block(0)

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.