]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Write a glossary entry on how to use the new facility
authorMatthias Maier <tamiko@43-1.org>
Sat, 24 Oct 2015 17:27:25 +0000 (12:27 -0500)
committerMatthias Maier <tamiko@43-1.org>
Mon, 28 Dec 2015 11:36:24 +0000 (12:36 +0100)
doc/doxygen/headers/constraints.h

index a2d6d1ed90b6f921e88247918119bc4b8bfdb67c..0a5ed06cb7dfdcaaa26ab2bd0376f0b4ae199aad 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2013 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
  * Either behavior can also be achieved by building two separate
  * ConstraintMatrix objects and calling ConstraintMatrix::merge function with
  * a particular second argument.
+ *
+ *
+ * <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
+ * condensed (or taken care of constraints during assembly). This is
+ * usually the case if the system is described by a LinearOperator.
+ *
+ * In this case we can solve the modified system
+ * @f[
+ *   (C^T \cdot A \cdot 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
+ * $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
+ * ConstraintMatrix::distribute() operation applied on a vector $x$ is the
+ * operation
+ * @f[
+    $x$ \leftarrow C\,x+k.
+ * @f]
+ * 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
+ * then recovered by distributing constraints to $\tilde x$: $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>
+ *
+ * // ...
+ *
+ * // system_matrix     - unconstrained and assembled system matrix
+ * // right_hand_side   - unconstrained and assembled right hand side
+ * // constraint_matrix - a ConstraintMatrix object
+ * // solver            - an appropriate, iterative solver
+ * // preconditioner    - a preconditioner
+ *
+ * const auto op_a = linear_operator(system_matrix);
+ * const auto op_amod = constrained_linear_operator(constraint_matrix, op_a);
+ * Vector<double> rhs_mod = constrained_right_hand_side(constraint_matrix,
+ *                                                      op_a,
+ *                                                      right_hand_side);
+ *
+ * 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.
  */

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.