]> https://gitweb.dealii.org/ - dealii.git/commitdiff
step-37: Convert ConstraintMatrix to AffineConstraints.
authorDavid Wells <drwells@email.unc.edu>
Sun, 26 Aug 2018 17:28:39 +0000 (13:28 -0400)
committerDavid Wells <drwells@email.unc.edu>
Mon, 27 Aug 2018 17:45:24 +0000 (13:45 -0400)
examples/step-37/doc/intro.dox
examples/step-37/doc/results.dox
examples/step-37/step-37.cc

index 6541e07af586faba4d7fe08a89ea350bc882e3cb..ad5fb814547e63c2868aa13b7a583c55f9457bb7 100644 (file)
@@ -147,7 +147,7 @@ Matrixfree<dim>::vmult (Vector<double>       &dst,
 
 Here we neglected boundary conditions as well as any hanging nodes we may
 have, though neither would be very difficult to include using the
-ConstraintMatrix class. Note how we first generate the local matrix in the
+AffineConstraints class. Note how we first generate the local matrix in the
 usual way as a sum over all quadrature points for each local matrix entry.
 To form the actual product as expressed in the above formula, we
 extract the values of <code>src</code> of the cell-related degrees of freedom
index b1049dbaab28bdab93b2bd82a45072b0b08b8097..e35e05a9502d52d03aedaa1c5cc87ac88e149515 100644 (file)
@@ -483,7 +483,7 @@ constrained by Dirichlet conditions.
 In the implementation in deal.II, the integrals $(\nabla \varphi_i,\nabla \varphi_j)_\Omega$
 on the right hand side are already contained in the local matrix contributions
 we assemble on each cell. When using
-ConstraintMatrix::distributed_local_to_global() as first described in the
+AffineConstraints::distributed_local_to_global() as first described in the
 step-6 and step-7 tutorial programs, we can account for the contribution of
 inhomogeneous constraints <i>j</i> by multiplying the columns <i>j</i> and
 rows <i>i</i> of the local matrix according to the integrals $(\varphi_i,
@@ -507,10 +507,10 @@ unrelated to the real entries.
 In a matrix-free method, we need to take a different approach, since the @p
 LaplaceOperator class represents the matrix-vector product of a
 <b>homogeneous</b> operator (the left-hand side of the last formula).  It does
-not matter whether the ConstraintMatrix passed to the MatrixFree::reinit()
-contains inhomogeneous constraints or not, the MatrixFree::cell_loop() call
-will only resolve the homogeneous part of the constraints as long as it
-represents a <b>linear</b> operator.
+not matter whether the AffineConstraints object passed to the
+MatrixFree::reinit() contains inhomogeneous constraints or not, the
+MatrixFree::cell_loop() call will only resolve the homogeneous part of the
+constraints as long as it represents a <b>linear</b> operator.
 
 In our matrix-free code, the right hand side computation where the
 contribution of inhomogeneous conditions ends up is completely decoupled from
@@ -530,7 +530,7 @@ where we only set the Dirichlet values:
       solution(it->first) = it->second;
 @endcode
 or, equivalently, if we already had filled the inhomogeneous constraints into
-a constraint matrix object,
+an AffineConstraints object,
 @code
   solution = 0;
   constraints.distribute(solution);
@@ -613,13 +613,13 @@ assemble_residual() that computes the (weak) form of the residual, whereas the
 @p LaplaceOperator::apply_add() function would get the linearization of the
 residual with respect to the solution variable.
 
-<h5> Use LaplaceOperator with a second ConstraintMatrix without Dirichlet conditions </h5>
+<h5> Use LaplaceOperator with a second AffineConstraints object without Dirichlet conditions </h5>
 
 A second alternative to get the right hand side that re-uses the
 @p LaplaceOperator::apply_add() function is to instead add a second constraint
 matrix that skips Dirichlet constraints on the read operation. To do this, we
 initialize a MatrixFree object in a more extended way with two different
-DoFHandler / ConstraintMatrix combinations. The 0-th component includes
+DoFHandler-AffineConstraints combinations. The 0-th component includes
 Dirichlet conditions for solving the linear system, whereas 1-st component
 does read also from Dirichlet-constrained degrees of freedom for the right
 hand side assembly:
@@ -680,7 +680,7 @@ void LaplaceProblem<dim>::assemble_rhs ()
 }
 @endcode
 
-Instead of adding a second DoFHandler / ConstraintMatrix pair to the same
+Instead of adding a second DoFHandler-AffineConstraints pair to the same
 MatrixFree::reinit() call, one could of course also construct an independent
 MatrixFree object that feeds the second @p LaplaceOperator instance, see also
 the discussion in MatrixFreeOperators::Base.
index aeed9f0c40f63ba4bbcadd90f19d5f8a6adcc54a..a7ea2e6f79f89b6ec4c2ae5c5e9434630d7ea081 100644 (file)
@@ -415,7 +415,7 @@ namespace Step37
   // $v_\mathrm{cell}$ as mentioned in the introduction need to be added into
   // the result vector (and constraints are applied). This is done with a call
   // to @p distribute_local_to_global, the same name as the corresponding
-  // function in the ConstraintMatrix (only that we now store the local vector
+  // function in the AffineConstraints (only that we now store the local vector
   // in the FEEvaluation object, as are the indices between local and global
   // degrees of freedom).  </ol>
   template <int dim, int fe_degree, typename number>
@@ -476,7 +476,7 @@ namespace Step37
   // Note that after the cell loop, the constrained degrees of freedom need to
   // be touched once more for sensible vmult() operators: Since the assembly
   // loop automatically resolves constraints (just as the
-  // ConstraintMatrix::distribute_local_to_global call does), it does not
+  // AffineConstraints::distribute_local_to_global() call does), it does not
   // compute any contribution for constrained degrees of freedom, leaving the
   // respective entries zero. This would represent a matrix that had empty
   // rows and columns for constrained degrees of freedom. However, iterative
@@ -552,7 +552,7 @@ namespace Step37
   // <tt>unsigned int</tt> in place of the source vector to confirm with the
   // cell_loop interface. After the loop, we need to set the vector entries
   // subject to Dirichlet boundary conditions to one (either those on the
-  // boundary described by the ConstraintMatrix object inside MatrixFree or
+  // boundary described by the AffineConstraints object inside MatrixFree or
   // the indices at the interface between different grid levels in adaptive
   // multigrid). This is done through the function
   // MatrixFreeOperators::Base::set_constrained_entries_to_one() and matches
@@ -718,7 +718,7 @@ namespace Step37
     FE_Q<dim>       fe;
     DoFHandler<dim> dof_handler;
 
-    ConstraintMatrix constraints;
+    AffineConstraints<double> constraints;
     using SystemMatrixType =
       LaplaceOperator<dim, degree_finite_element, double>;
     SystemMatrixType system_matrix;
@@ -885,7 +885,7 @@ namespace Step37
         DoFTools::extract_locally_relevant_level_dofs(dof_handler,
                                                       level,
                                                       relevant_dofs);
-        ConstraintMatrix level_constraints;
+        AffineConstraints<double> level_constraints;
         level_constraints.reinit(relevant_dofs);
         level_constraints.add_lines(
           mg_constrained_dofs.get_boundary_indices(level));

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.