]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Overhaul the documentation of Concepts.
authorDavid Wells <wellsd2@rpi.edu>
Wed, 9 Dec 2015 14:04:43 +0000 (09:04 -0500)
committerDavid Wells <wellsd2@rpi.edu>
Sat, 12 Dec 2015 23:47:07 +0000 (18:47 -0500)
doc/doxygen/headers/concepts.h [new file with mode: 0644]
doc/doxygen/headers/glossary.h
doc/doxygen/headers/laoperators.h
doc/doxygen/headers/matrices.h
doc/doxygen/headers/preconditioners.h
doc/doxygen/headers/vectors.h
doc/news/changes.h
include/deal.II/multigrid/mg_smoother.h

diff --git a/doc/doxygen/headers/concepts.h b/doc/doxygen/headers/concepts.h
new file mode 100644 (file)
index 0000000..1cc7cb0
--- /dev/null
@@ -0,0 +1,210 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+/**
+ * @defgroup Concepts Concepts, or expectations on template parameters
+ *
+ * Sometimes imposing constraints on the type of an object without requiring
+ * it to belong to a specific inheritance hierarchy is useful. These are
+ * usually referred to as <em>concepts</em> in the C++ community. This module
+ * lists the concepts commonly used in deal.II with brief descriptions of
+ * their intent. The convention in deal.II for listing constraints on a type
+ * is to provide the name of the concept as a <code>typename</code> in a
+ * template: for example, the type of a Vector depends on the type of the
+ * underlying field, but Vector is not a field element, so it is defined as a
+ * template:
+ * @code
+ * template <typename Number>
+ * class Vector;
+ * @endcode
+ * Where @ref ConceptNumber "Number" is understood to be an appropriate field
+ * element type.
+ *
+ * <dl>
+ *
+ * <dt class="concepts">@anchor ConceptContainerType <b>ContainerType</b></dt>
+ *
+ * <dd>
+
+ * There are several algorithms (e.g.,
+ * GridTools::find_active_cell_around_point) in deal.II that can operate on
+ * either a Triangulation or a DoFHandler, as both classes may be considered
+ * to be collections of cells: see the @ref GlossMeshAsAContainer
+ * "glossary entry" for a further discussion of this idea. %Functions that may
+ * be called with either class indicate this by accepting a template parameter
+ * like
+ * @code
+ * template <template <int, int> class Container>
+ * @endcode
+ * or
+ * @code
+ * template <typename Container>
+ * @endcode
+ * which is usually required to have a <code>typedef</code> named
+ * <code>active_cell_iterator</code>.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptDoFHandlerType <b>DoFHandlerType</b></dt>
+ *
+ * <dd>
+ * deal.II includes both DoFHandler and hp::DoFHandler as objects which manage
+ * degrees of freedom on a mesh. Though the two do not share any sort of
+ * inheritance relationship, they are similar enough that many functions just
+ * need something which resembles a DoFHandler to work correctly.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptMatrixType <b>MatrixType</b></dt>
+ *
+ * <dd>
+ * Many functions and classes in deal.II require an object which knows how to
+ * calculate matrix-vector products (the member function <code>vmult</code>),
+ * transposed matrix-vector products (the member function
+ * <code>Tvmult</code>), as well as the `multiply and add' equivalents
+ * <code>vmult_add</code> and <code>Tvmult_add</code>. Some functions only
+ * require <code>vmult</code> and <code>Tvmult</code>, but an object should
+ * implement all four member functions if the template requires a MatrixType
+ * argument. Writing classes that satisfy these conditions is a sufficiently
+ * common occurrence that the LinearOperator class was written to make things
+ * easier; see @ref LAOperators for more information.
+ *
+ * One way to think of <code>MatrixType</code> is to pretend it is a base
+ * class with the following signature (this is nearly the interface provided
+ * by SparseMatrix):
+ *
+ * @code
+ * class MatrixType
+ * {
+ * public:
+ *   template <typename VectorType>
+ *   virtual void vmult(VectorType &u, const VectorType &v) const =0;
+ *
+ *   template <typename VectorType>
+ *   virtual void Tvmult(VectorType &u, const VectorType &v) const =0;
+ *
+ *   template <typename VectorType>
+ *   virtual void vmult_add(VectorType &u, const VectorType &v) const =0;
+ *
+ *   template <typename VectorType>
+ *   virtual void Tvmult_add(VectorType &u, const VectorType &v) const =0;
+ * };
+ * @endcode
+ *
+ * Template functions in C++ cannot be virtual (which is the main reason why
+ * this approach is not used in deal.II), so implementing this interface with
+ * inheritance will not work, but it is still a good way to think about this
+ * template concept. One can use the PointerMatrixAux class to implement
+ * <code>vmult_add</code> and <code>Tvmult_add</code> instead of implementing
+ * them manually.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptNumber <b>Number</b></dt>
+ *
+ * <dd>
+ * This concept describes scalars which make sense as vector or matrix
+ * entries, which is usually some finite precision approximation of a field
+ * element. The canonical examples are <code>double</code> and
+ * <code>float</code>, but deal.II supports <code>std::complex&lt;T&gt;</code>
+ * for floating point type <code>T</code> in many places as well.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptPolynomialType <b>PolynomialType</b></dt>
+ *
+ * <dd>
+ * See the description in @ref Polynomials for more information. In some
+ * contexts, anything that satisfies the interface resembling
+ * @code
+ * template <int dim>
+ * class PolynomialType
+ * {
+ *   virtual void compute (const Point<dim>            &unit_point,
+ *                         std::vector<Tensor<1,dim> > &values,
+ *                         std::vector<Tensor<2,dim> > &grads,
+ *                         std::vector<Tensor<3,dim> > &grad_grads) const =0;
+ * }
+ * @endcode
+ *
+ * may be considered as a polynomial for the sake of implementing finite
+ * elements.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptPreconditionerType <b>PreconditionerType</b></dt>
+ *
+ * <dd>
+ * This is essentially a synonym for <code>MatrixType</code>, but usually only
+ * requires that <code>vmult()</code> and <code>Tvmult()</code> be
+ * defined. Most of the time defining <code>Tvmult()</code> is not
+ * necessary. One should think of <code>vmult()</code> as applying some
+ * approximation of the inverse of a linear operator to a vector, instead of
+ * the action of a linear operator to a vector, for the preconditioner
+ * classes.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptRelaxationType <b>RelaxationType</b></dt>
+ *
+ * <dd>
+ * This is an object capable of relaxation for multigrid methods. One can
+ * think of an object satisfying this constraint as having the following
+ * interface as well as the constraints required by
+ * @ref ConceptMatrixType "MatrixType":
+ * @code
+ * class RelaxationType
+ * {
+ * public:
+ *   template <typename VectorType>
+ *   virtual void step(VectorType &u, const VectorType &v) const =0;
+ *
+ *   template <typename VectorType>
+ *   virtual void Tstep(VectorType &u, const VectorType &v) const =0;
+ * };
+ * @endcode
+ * where these two member functions perform one step (or the transpose of such
+ * a step) of the smoothing scheme. In other words, the operation performed by
+ * these functions is
+ * $dst = dst - P^{-1} (A dst - rhs)$ and $dst = dst - P^{-T} (A dst - rhs)$.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptSparsityPatternType <b>SparsityPatternType</b></dt>
+ *
+ * <dd>
+ * Almost all functions (with the notable exception of
+ * SparsityTools::distribute_sparsity_pattern) which take a sparsity pattern
+ * as an argument can take either a regular SparsityPattern or a
+ * DynamicSparsityPattern, or even one of the block sparsity patterns. See
+ * @ref Sparsity for more information.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptStreamType <b>StreamType</b></dt>
+ *
+ * <dd>
+ * Deriving new stream classes in C++ is well-known to be difficult. To get
+ * around this, some functions accept a parameter which defines
+ * <code>operator&lt;&lt;</code>, which allows for easy output to any kind of
+ * output stream.
+ * </dd>
+ *
+ * <dt class="concepts">@anchor ConceptVectorType <b>VectorType</b></dt>
+ *
+ * <dd>
+ * deal.II supports many different vector classes, including bindings to
+ * vectors in other libraries. These are similar to standard library vectors
+ * (i.e., they define <code>begin()</code>, <code>end()</code>,
+ * <code>operator[]</code>, and <code>size()</code>) but also define numerical
+ * operations like <code>add()</code>. Some examples of VectorType include
+ * Vector, TrilinosWrappers::MPI::Vector, and BlockVector.
+ * </dd>
+ *
+ * </dl>
+ */
index d47ae7e707ca89f87dc650a8db829fdab25d70c8..6891129f085f808457151af55e96d177a2192f44 100644 (file)
  * @ref GlossGhostedVector "vectors with ghost elements".
  *
  *
+ * <dt class="glossary">@anchor GlossConcept <b>Concepts in deal.II</b></dt>
+ *
+ * <dd> There are several places in deal.II where we require that a type in a
+ * template match a certain interface or behave in a certain way: such
+ * constraints are called <em>concepts</em> in C++. See the discussion in
+ * @ref Concepts for more information and a list of concepts in deal.II.
+ * </dd>
+ *
  * <dt class="glossary">@anchor GlossDoF <b>Degree of freedom</b></dt>
  *
  * <dd> The term "degree of freedom" (often abbreviated as "DoF") is commonly
index 0aba9b36986e535c5e1a742c9b9dad3e33e13387..1af479d6cf0b0bf3716324f489578ab5a9cf5f1f 100644 (file)
  * linear operator is available. (For questions about C++11, see
  * @ref CPP11 .)
  *
- * This is done with a LinearOperator class that, similarly to the abstract
- * MATRIX interface, defines a minimal interface for <i>applying</i> a
- * linear operation on a vector.
+ * This is done with a LinearOperator class that, like
+ * @ref ConceptMatrixType "the MatrixType concept",
+ * defines a minimal interface for <i>applying</i> a linear operation on a
+ * vector.
+ *
  * @code
  *   std::function<void(Range &, const Domain &)> vmult;
  *   std::function<void(Range &, const Domain &)> vmult_add;
index 33ac66e49f0e9857a5efedf04714467e50feb744..ee54556d6f1b42e0a1cd9b127c3e820ff22cf8a5 100644 (file)
  * @ingroup LAC
  */
 
-/**
- * Template for matrix classes.
- *
- * @note This is a description of the expectations on the MATRIX
- * template argument. It is not a class by itself.
- *
- * Depending on where the MATRIX is used, its interface is expected to
- * conform to one of the following groups, which can also be found in
- * the overview.
- *
- * <h3>Solver interface</h3>
- *
- * Functions in this group are the minimal interface to use the MATRIX
- * in a linear Solver. Solvers use a matrix only as a linear operator,
- * that is, they map a vector to another. To this end, we either
- * multiply with the matrix itself or its transpose.
- *
- * The function vmult() is the bare necessity in this group. Some
- * solvers use Tvmult() as well, in which case it needs to be
- * implemented. Some derived matrices like PointerMatrix require its
- * existence, in which case it can be implemented empty with an
- * assertion <code>Assert(false, ExcNotImplemented())</code>.
- *
- * If vmult_add() and Tvmult_add() are missing, PointerMatrixAux can
- * be used to provide the missing functionality without implementing
- * it by hand.
- *
- * @ingroup LAC
- */
-template <class VECTOR>
-class MATRIX
-{
-  public:
-                                     /**
-                                      * @name Solver interface
-                                      */
-                                     /*@{*/
-                                     /**
-                                      * The matrix vector product $u = Av$.
-                                      */
-    void vmult(VECTOR& u, const VECTOR& v) const;
-                                     /**
-                                      * The matrix vector product $u = A^Tv$.
-                                      */
-    void Tvmult(VECTOR& u, const VECTOR& v) const;
-                                     /**
-                                      * The matrix vector product $u += Av$.
-                                      */
-    void vmult_add(VECTOR& u, const VECTOR& v) const;
-                                     /**
-                                      * The matrix vector product $u += A^Tv$.
-                                      */
-    void Tvmult_add(VECTOR& u, const VECTOR& v) const;
-                                     /*@}*/
-};
 
 /**
  * @defgroup Matrix1 Basic matrices
@@ -120,7 +65,8 @@ class MATRIX
  * @defgroup Matrix2 Derived matrices
  *
  * These matrices are built on top of the basic matrices. They perform special
- * operations using the interface defined in @ref Solvers.
+ * operations using the interface defined by
+ * @ref ConceptMatrixType "the MatrixType concept".
  *
  * @ingroup Matrices
  */
index c6816e9a73c171fcb21803089a96f74ad273851b..703514ab6cf5e0d93bd78fc04b75fb38f613e2e0 100644 (file)
@@ -29,7 +29,7 @@
  * Broadly speaking, preconditioners are operators, which are multiplied with
  * a matrix to improve conditioning. The idea is, that the preconditioned
  * system <i>P<sup>-1</sup>Ax = P<sup>-1</sup>b</i> is much easier to solve
- * than the original system <i>Ax = b</i>.What this means exactly depends on
+ * than the original system <i>Ax = b</i>. What this means exactly depends on
  * the structure of the matrix and cannot be discussed here in generality. For
  * symmetric, positive definite matrices <i>A</i> and <i>P</i>, it means that
  * the spectral condition number (the quotient of greatest and smallest
  * @f]
  * Accordingly, preconditioning amounts to applying a linear operator to the
  * residual, and consequently, the action of the preconditioner
- * <i>P<sup>-1</sup></i> is implemented as <tt>vmult()</tt>. The
- * generic interface is like for matrices
- * @code
- * class PRECONDITIONER
- * {
- *   template <class VECTOR>
- *   void vmult(VECTOR& dst, const VECTOR& src) const;
+ * <i>P<sup>-1</sup></i> is implemented as <tt>vmult()</tt>.
+ * Templates in deal.II that require a preconditioner indicate the
+ * requirement with
+ * @ref ConceptPreconditionerType "the PreconditionerType concept". In
+ * practice, one can usually treat any matrix-like object which defines
+ * <code>vmult()</code> and <code>Tvmult()</code> as a preconditioner. All
+ * preconditioner classes in this module implement this interface.
  *
- *   template <class VECTOR>
- *   void Tvmult(VECTOR& dst, const VECTOR& src) const;
- * }
- * @endcode
- * It is implemented in all the preconditioner classes in this module.
- * 
  * When used
  * in Krylov space methods, it is up to the method, whether it simply
  * replaces multiplications with <i>A</i> by those with
  * thus avoiding multiplication with <i>A</i> completely. We call
  * operators mapping the previous iterate <i>x<sup>k</sup></i> to the
  * next iterate in this way relaxation operators. Their generic
- * interface is
- * @code
- * class RELAXATION
- * {
- *   template <class VECTOR>
- *   void step(VECTOR& newstep, const VECTOR& oldstep, const VECTOR& rhs) const;
- *
- *   template <class VECTOR>
- *   void Tstep(VECTOR& newstep, const VECTOR& oldstep, const VECTOR& rhs) const;
- * }
- * @endcode
- * The classes with names starting with <tt>Relaxation</tt> in this
- *   module implement this interface, as well as the preconditioners
- *   PreconditionJacobi, PreconditionSOR, PreconditionSSORP,
- *   reconditionBlockJacobi, PreconditionBlockSOR, and
- *   PreconditionBlockSSOR.
+ * interface is given by @ref ConceptRelaxationType "the RelaxationType concept".
+ * The classes with names starting with <tt>Relaxation</tt> in this module
+ * implement this interface, as well as the preconditioners
+ * PreconditionJacobi, PreconditionSOR, PreconditionSSORP,
+ * reconditionBlockJacobi, PreconditionBlockSOR, and PreconditionBlockSSOR.
  *
  * <h3>The interface</h3>
  *
  * additional required parameters and sets up the internal structures
  * of the preconditioner.
  *
- * <h4>Application of the preconditioner</h4>
- *
- * Preconditioners in deal.II are just considered linear operators.
- * Therefore, in order to be used by deal.II solvers, preconditioners must
- * conform to the standard matrix interface, namely the functions
- * @code
- * void  vmult (VECTOR& dst, const VECTOR& src) const;
- * void  Tvmult (VECTOR& dst, const VECTOR& src) const;
- * @endcode
- * These functions apply the preconditioning operator to the source vector
- * $src$ and return the result in $dst$ as $dst=P^{-1}src$ or
- * $dst=P^{-T}src$. Preconditioned iterative dolvers use these
- * <tt>vmult()</tt> functions of the preconditioner.  Some solvers may also
- * use <tt>Tvmult()</tt>.
- *
  * <h4>Relaxation methods</h4>
  *
- * Additional to the interface described above, some preconditioners like SOR
- * and Jacobi have been known as iterative methods themselves. For them, an
- * additional interface exists, consisting of the functions
- * @code
- * void  step (VECTOR& dst, const VECTOR& rhs) const;
- * void  Tstep (VECTOR& dst, const VECTOR& rhs) const;
- * @endcode
- *
- * Here, $src$ is a residual vector and $dst$ is the iterate that is supposed
- * to be updated. In other words, the operation performed by these functions
- * is
- * $dst = dst - P^{-1} (A dst - rhs)$ and $dst = dst - P^{-T} (A dst - rhs)$.
- * The functions are called this way because they perform <i>one step</i> of a
- * fixed point (Richardson) iteration. Note that preconditioners store a
- * reference to the original matrix $A$ during initialization.
+ * Some preconditioners, like SOR and Jacobi, were used as iterative solvers
+ * long before they were used as preconditioners. Thus, they satisfy both
+ * @ref ConceptMatrixType "MatrixType" and
+ * @ref ConceptRelaxationType "RelaxationType" concepts.
  *
  * @ingroup LAC
  * @ingroup Matrices
index e429566ab0e07bd16f446d35464247700f415384..a0767f664fc61d16fcd3106330dbf1e895c7b3b3 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2003 - 2013 by the deal.II authors
+// Copyright (C) 2003 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -17,8 +17,9 @@
 /**
  * @defgroup Vectors Vector classes
  *
- * Here, we list all the classes that can be used as VECTOR in linear solvers
- * (see @ref Solvers) and for matrix-vector operations.
+ * Here, we list all the classes that satisfy the <code>VectorType</code>
+ * concept and may be used in linear solvers (see @ref Solvers) and for
+ * matrix-vector operations.
  *
  * @ingroup LAC
  */
index a33ee863d3cdf9cdd0375d1307a09fbcba54b8dc..f2417a9f1150f7325e422e9484f35f5ec30475ab 100644 (file)
@@ -192,6 +192,12 @@ inconvenience this causes.
   (Rene Gassmoeller, 2015/12/09)
   </li>
 
+  <li> New: There is a new module, @ref Concepts, which describes the meaning
+  behind template parameter type names.
+  <br>
+  (David Wells, 2015/12/09)
+  </li>
+
   <li> New: The WorkStream class's design and implementation are now much
   better documented in the form of a @ref workstream_paper "preprint".
   <br>
index 2906a57b3858c393d29da782b7645020f1b07696..ac7eb9b81c58ce0b00399630adc6dab412ae8414 100644 (file)
@@ -150,12 +150,8 @@ namespace mg
   /**
    * Smoother using relaxation classes.
    *
-   * A relaxation class is an object that has two member functions,
-   * @code
-   * void  step(VectorType& x, const VectorType& d) const;
-   * void Tstep(VectorType& x, const VectorType& d) const;
-   * @endcode
-   * performing one step of the smoothing scheme.
+   * A relaxation class is an object that satisfies the
+   * @ref ConceptRelaxationType "RelaxationType concept".
    *
    * This class performs smoothing on each level. The operation can be
    * controlled by several parameters. First, the relaxation parameter @p

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.