]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate MeshWorker::LocalIntegrator.
authorWolfgang Bangerth <bangerth@colostate.edu>
Tue, 21 May 2024 03:10:10 +0000 (21:10 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 23 May 2024 12:54:12 +0000 (06:54 -0600)
include/deal.II/meshworker/local_integrator.h
include/deal.II/meshworker/local_results.h
include/deal.II/meshworker/loop.h
include/deal.II/meshworker/vector_selector.h

index 785a72b197848872068d97a157d97dabb25ce069..ff0dace4bc858ebcab8ff9e34730bdc1f6d75f42 100644 (file)
@@ -49,10 +49,14 @@ namespace MeshWorker
    * If a function is not overloaded in a derived class, but its usage flag is
    * true, the function will cause an exception ExcPureFunction.
    *
+   * @deprecated This class is deprecated. It used to be the basis for
+   *   integration via the MeshWorker::integration_loop() function, but the
+   *   same functionality is available via MeshWorker::loop().
+   *
    * @ingroup MeshWorker
    */
   template <int dim, int spacedim = dim, typename number = double>
-  class LocalIntegrator : public Subscriptor
+  class DEAL_II_DEPRECATED_EARLY LocalIntegrator : public Subscriptor
   {
   public:
     /**
index 3a5a988bfff6af65d017aeac30fa165a592644d9..d23468e75eec962bd9c44b0c5a3f620972e09d0e 100644 (file)
@@ -39,14 +39,12 @@ class BlockIndices;
  * ubiquitous part of each finite element program.
  *
  * The workhorse of this namespace is the loop() function, which implements a
- * completely generic loop over all mesh cells. Since the calls to loop() are
- * error-prone due to its generality, for many applications it is advisable to
- * derive a class from MeshWorker::LocalIntegrator and use the less general
- * integration_loop() instead.
- *
- * The loop() depends on certain objects handed to it as arguments. These
+ * completely generic loop over all mesh cells.
+ * The loop() function depends on certain objects handed to it as arguments.
+ * These
  * objects are of two types, @p info objects like DoFInfo and IntegrationInfo and
- * worker objects like LocalWorker and IntegrationWorker.
+ * function objects ("workers") that perform the local operations on cells,
+ * faces, and boundaries.
  *
  * Worker objects usually do two different jobs: first, they compute the local
  * contribution of a cell or face to the global operation. Second, they
@@ -155,23 +153,6 @@ class BlockIndices;
  * };
  * @endcode
  *
- * <h3>Simplified interfaces</h3>
- *
- * Since the loop() is fairly general, a specialization integration_loop() is
- * available, which is a wrapper around loop() with a simplified interface.
- *
- * The integration_loop() function loop takes most of the information that it
- * needs to pass to loop() from an IntegrationInfoBox object. Its use is
- * explained in step-12, but in short it requires functions that do the local
- * integration on a cell, interior or boundary face, and it needs an object
- * (called "assembler") that copies these local contributions into the global
- * matrix and right hand side objects.
- *
- * Before we can run the integration loop, we have to initialize several data
- * structures in our IntegrationWorker and assembler objects. For instance, we
- * have to decide on the quadrature rule or we may need more than the default
- * update flags.
- *
  * @ingroup MeshWorker
  * @ingroup Integrators
  */
index 3fcea362f28d9abacbfbc3fbd7cb568d4670833f..c8a2c97c18ecbd3fefdbe3fec099afb31e5b22ed 100644 (file)
@@ -491,13 +491,17 @@ namespace MeshWorker
    * Simplified interface for loop() if specialized for integration, using the
    * virtual functions in LocalIntegrator.
    *
+   * @deprecated This function is deprecated, along with the LocalIntegrator
+   *   class. Use the MeshWorker::loop() function directly, with three function
+   *   objects that perform the cell, boundary, and interior face integration.
+   *
    * @ingroup MeshWorker
    */
   template <int dim,
             int spacedim,
             typename IteratorType,
             typename AssemblerType>
-  void
+  DEAL_II_DEPRECATED_EARLY void
   integration_loop(IteratorType                             begin,
                    std_cxx20::type_identity_t<IteratorType> end,
                    DoFInfo<dim, spacedim>                  &dof_info,
index a74f16d7d416427fe8ddce22ecb22de3d8355aaa..9b8e7e703ec1d42f2f3b48e6b86c7eab1e3a6a5f 100644 (file)
@@ -42,9 +42,6 @@ namespace MeshWorker
    * those, which are actually used in computing residuals etc. This class
    * organizes the selection.
    *
-   * It is used for instance in IntegrationWorker to determine which values,
-   * derivatives or second derivatives are actually computed.
-   *
    * @ingroup MeshWorker
    */
   class VectorSelector : public Subscriptor

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.