]> https://gitweb.dealii.org/ - dealii.git/commitdiff
address comments
authorMatthias Maier <tamiko@43-1.org>
Sun, 25 Oct 2015 19:14:22 +0000 (14:14 -0500)
committerMatthias Maier <tamiko@43-1.org>
Mon, 28 Dec 2015 19:21:38 +0000 (20:21 +0100)
doc/doxygen/headers/constraints.h
doc/news/changes.h
include/deal.II/lac/constrained_linear_operator.h [moved from include/deal.II/lac/constraint_linear_operator.h with 94% similarity]
tests/lac/constraint_linear_operator_01.cc

index 93f740e4ae71d816521e31b1e7d2cac31ff3560a..ab9c8179d65fa1a3f9bea23c26ab663666b2f9e3 100644 (file)
  * <h3>Applying constraints indirectly with a LinearOperator</h3>
  *
  * Sometimes it is either not desirable, or not possible to directly
- * condensate, or eliminate constraints from a system of linear equations.
- * In particular if there is no underlying matrix object that could be
+ * condense, or eliminate constraints from a system of linear equations. In
+ * particular if there is no underlying matrix object that could be
  * condensed (or taken care of constraints during assembly). This is
  * usually the case if the system is described by a LinearOperator.
  *
  * @f[
  *   (C^T A C + Id_c) \tilde x = C^T (b - A\,k)
  * @f]
- * instead [1]. Here, $A$ is a given (unconstrained) system matrix $A$ and
- * $b$ the corresponding right hand side of a system of linear equations
+ * instead [1] (M. S. Shephard: Linear multipoint constraints applied via
+ * transformation as part of a direct stiffness assembly process, 1985).
+ *
+ * Here, $A$ is a given (unconstrained) system matrix and $b$
+ * the corresponding right hand side of a system of linear equations
  * $A\,x=b$. The matrix $C$ describes the homogeneous part of the linear
  * constraints stored in a ConstraintMatrix and the vector $k$ is the
  * vector of corresponding inhomogeneities. More precisely, the
  * And finally, $Id_c$ denotes the identity on the subspace of constrained
  * degrees of freedom.
  *
- * The corresponding solution of $A\,x=b$ that obeys boundary conditions is
+ * The corresponding solution of $A\,x=b$ that obeys these constraints is
  * then recovered by distributing constraints: $x=C\tilde x+k$.
  *
  * The whole system can be set up and solved with the following snippet of
  * code:
  * @code
- * #include <deal.II/lac/constraint_linear_operator.h>
+ * #include <deal.II/lac/constrained_linear_operator.h>
  *
  * // ...
  *
  * solver.solve(op_amod, solution, rhs_mod, preconditioner);
  * constraint_matrix.distribute(solution);
  * @endcode
- *
- * [1] M. S. Shephard: Linear multipoint constraints applied via
- * transformation as part of a direct stiffness assembly process, 1985.
  */
index 41390297b76d13fc468a973a49db928c760b78da..7e03f5a70038caaa65b71bb88df0e57231984945 100644 (file)
@@ -184,9 +184,9 @@ inconvenience this causes.
 <h3>General</h3>
 <ol>
   <li> New: constrained_linear_operator() and constrained_right_hand_side()
-  implemented that provide a generic mechanism of applying constraints to a
-  LinearOperator. A detailed explanation with example code is given in the
-  @ref constraints module.
+  provide a generic mechanism of applying constraints to a LinearOperator.
+  A detailed explanation with example code is given in the @ref constraints
+  module.
   <br>
   (Mauro Bardelloni, Matthias Maier, 2015/10/25 - 2015/12/27)
   </li>
similarity index 94%
rename from include/deal.II/lac/constraint_linear_operator.h
rename to include/deal.II/lac/constrained_linear_operator.h
index 699d63a3b0e20c968d8456336bb91d2b5dc8054c..0597ae2baa40f81d7f27fda604ccf590937276eb 100644 (file)
@@ -13,8 +13,8 @@
 //
 // ---------------------------------------------------------------------
 
-#ifndef dealii__constraint_linear_operator_h
-#define dealii__constraint_linear_operator_h
+#ifndef dealii__constrained_linear_operator_h
+#define dealii__constrained_linear_operator_h
 
 #include <deal.II/lac/linear_operator.h>
 #include <deal.II/lac/packaged_operation.h>
@@ -32,15 +32,16 @@ DEAL_II_NAMESPACE_OPEN
 
 
 /**
- * This function takes a ConstraintMatrix @p cm and an operator exemplar @p
- * exemplar (this exemplar is usually a linear operator that describes the
- * system matrix) and returns a LinearOperator object associated with the
- * "homogeneous action" of the underlying ConstraintMatrix object:
+ * This function takes a ConstraintMatrix @p constraint_matrix and an
+ * operator exemplar @p exemplar (this exemplar is usually a linear
+ * operator that describes the system matrix) and returns a LinearOperator
+ * object associated with the "homogeneous action" of the underlying
+ * ConstraintMatrix object:
  *
  * Applying the LinearOperator object on a vector <code>u</code> results in
  * a vector <code>v</code> that stores the result of calling
  * ConstraintMatrix::distribute() on <code>u</code> - with one important
- * difference; inhomogeneities are note applied, but always treated as 0
+ * difference: inhomogeneities are note applied, but always treated as 0
  * instead.
  *
  * The LinearOperator object created by this function is primarily used
index 188f8baa4d712d60883188e722909d414676a4db..f0c6e4e5c71b76c5649941b2171e0d0f4a3ff45d 100644 (file)
@@ -45,7 +45,7 @@
 #include <deal.II/lac/constraint_matrix.h>
 #include <deal.II/lac/linear_operator.h>
 #include <deal.II/lac/packaged_operation.h>
-#include <deal.II/lac/constraint_linear_operator.h>
+#include <deal.II/lac/constrained_linear_operator.h>
 
 #include <fstream>
 #include <iostream>
@@ -192,8 +192,8 @@ void Step6<dim>::solve ()
 
   const auto A = linear_operator(system_matrix_lo);
   const auto M = constrained_linear_operator(constraints, A);
-  Vector<double> rhs =
-    constrained_right_hand_side(constraints, A, system_rhs_lo);
+  const auto rhs = constrained_right_hand_side(constraints, A,
+                                               system_rhs_lo);
 
   check_solver_within_range(
     solver.solve(M, solution_lo, rhs, preconditioner),

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.