]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make VectorTools::compute_mean_value() use ReadVector.
authorDavid Wells <drwells@email.unc.edu>
Fri, 2 Jun 2023 20:08:10 +0000 (16:08 -0400)
committerDavid Wells <drwells@email.unc.edu>
Sat, 1 Jul 2023 14:22:08 +0000 (10:22 -0400)
include/deal.II/numerics/vector_tools_mean_value.h
include/deal.II/numerics/vector_tools_mean_value.templates.h
source/numerics/vector_tools_mean_value.inst.in

index c366a22ebc53b13a8ce117bd80274b7ae7dad826..4282935aaa8ecae9031eaf74d1b4b5f379d36190 100644 (file)
@@ -19,6 +19,8 @@
 
 #include <deal.II/base/config.h>
 
+#include <deal.II/lac/read_vector.h>
+
 #include <vector>
 
 DEAL_II_NAMESPACE_OPEN
@@ -161,47 +163,39 @@ namespace VectorTools
    * finite element function with mean value zero. In fact, it only works for
    * Lagrangian elements. For all other elements, you will need to compute the
    * mean value and subtract it right inside the evaluation routine.
-   *
-   * @dealiiConceptRequires{concepts::is_dealii_vector_type<VectorType>}
    */
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type compute_mean_value(
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(
     const hp::MappingCollection<dim, spacedim> &mapping_collection,
     const DoFHandler<dim, spacedim> &           dof,
     const hp::QCollection<dim> &                q_collection,
-    const VectorType &                          v,
+    const ReadVector<Number> &                  v,
     const unsigned int                          component);
 
   /**
    * Calls the other compute_mean_value() function, see above, for the non-hp
    * case. That means, it requires a single FiniteElement, a single Quadrature,
    * and a single Mapping object.
-   *
-   * @dealiiConceptRequires{concepts::is_dealii_vector_type<VectorType>}
    */
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type
-    compute_mean_value(const Mapping<dim, spacedim> &   mapping,
-                       const DoFHandler<dim, spacedim> &dof,
-                       const Quadrature<dim> &          quadrature,
-                       const VectorType &               v,
-                       const unsigned int               component);
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(const Mapping<dim, spacedim> &   mapping,
+                     const DoFHandler<dim, spacedim> &dof,
+                     const Quadrature<dim> &          quadrature,
+                     const ReadVector<Number> &       v,
+                     const unsigned int               component);
 
   /**
    * Call the other compute_mean_value() function, see above, with
    * <tt>mapping=MappingQ@<dim@>(1)</tt>.
-   *
-   * @dealiiConceptRequires{concepts::is_dealii_vector_type<VectorType>}
    */
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type
-    compute_mean_value(const DoFHandler<dim, spacedim> &dof,
-                       const Quadrature<dim> &          quadrature,
-                       const VectorType &               v,
-                       const unsigned int               component);
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(const DoFHandler<dim, spacedim> &dof,
+                     const Quadrature<dim> &          quadrature,
+                     const ReadVector<Number> &       v,
+                     const unsigned int               component);
   /** @} */
 } // namespace VectorTools
 
index 7b3be8aed4d27a869ca03194ae3c07c004ee6d1c..654b3727d10a27d4cdc31bfe5f91bbd8320eb145 100644 (file)
@@ -290,17 +290,15 @@ namespace VectorTools
 
 
 
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type compute_mean_value(
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(
     const hp::MappingCollection<dim, spacedim> &mapping_collection,
     const DoFHandler<dim, spacedim> &           dof,
     const hp::QCollection<dim> &                q_collection,
-    const VectorType &                          v,
+    const ReadVector<Number> &                  v,
     const unsigned int                          component)
   {
-    using Number = typename VectorType::value_type;
-
     const hp::FECollection<dim, spacedim> &fe_collection =
       dof.get_fe_collection();
     const unsigned int n_components = fe_collection.n_components();
@@ -371,14 +369,13 @@ namespace VectorTools
   }
 
 
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type
-    compute_mean_value(const Mapping<dim, spacedim> &   mapping,
-                       const DoFHandler<dim, spacedim> &dof,
-                       const Quadrature<dim> &          quadrature,
-                       const VectorType &               v,
-                       const unsigned int               component)
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(const Mapping<dim, spacedim> &   mapping,
+                     const DoFHandler<dim, spacedim> &dof,
+                     const Quadrature<dim> &          quadrature,
+                     const ReadVector<Number> &       v,
+                     const unsigned int               component)
   {
     return compute_mean_value(hp::MappingCollection<dim, spacedim>(mapping),
                               dof,
@@ -388,13 +385,12 @@ namespace VectorTools
   }
 
 
-  template <int dim, typename VectorType, int spacedim>
-  DEAL_II_CXX20_REQUIRES(concepts::is_dealii_vector_type<VectorType>)
-  typename VectorType::value_type
-    compute_mean_value(const DoFHandler<dim, spacedim> &dof,
-                       const Quadrature<dim> &          quadrature,
-                       const VectorType &               v,
-                       const unsigned int               component)
+  template <int dim, typename Number, int spacedim>
+  Number
+  compute_mean_value(const DoFHandler<dim, spacedim> &dof,
+                     const Quadrature<dim> &          quadrature,
+                     const ReadVector<Number> &       v,
+                     const unsigned int               component)
   {
     return compute_mean_value(get_default_linear_mapping(
                                 dof.get_triangulation()),
index bb9c726e9dc0f31651dcf21729653dea9ac6354e..4555c80c3f2bcdbd6308294ce8ea64b49d19002d 100644 (file)
 //
 // ---------------------------------------------------------------------
 
-
-for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
+for (S : REAL_AND_COMPLEX_SCALARS; deal_II_dimension : DIMENSIONS;
      deal_II_space_dimension : SPACE_DIMENSIONS)
   {
 #if deal_II_dimension <= deal_II_space_dimension
     namespace VectorTools
     \{
-
-      template VEC::value_type
-      compute_mean_value<deal_II_dimension>(
+      template S
+      compute_mean_value(
         const Mapping<deal_II_dimension, deal_II_space_dimension> &,
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
         const Quadrature<deal_II_dimension> &,
-        const VEC &,
+        const ReadVector<S> &,
         const unsigned int);
 
-      template VEC::value_type
-      compute_mean_value<deal_II_dimension>(
+      template S
+      compute_mean_value(
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
         const Quadrature<deal_II_dimension> &,
-        const VEC &,
+        const ReadVector<S> &,
         const unsigned int);
+    \}
+#endif
+  }
 
+for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
+     deal_II_space_dimension : SPACE_DIMENSIONS)
+  {
+#if deal_II_dimension <= deal_II_space_dimension
+    namespace VectorTools
+    \{
       template void
-      add_constant<>(VEC &                                      solution,
-                     const DoFHandler<deal_II_dimension,
-                                      deal_II_space_dimension> &dof_handler,
-                     const unsigned int                         component,
-                     const typename VEC::value_type constant_adjustment);
+      add_constant(VEC &solution,
+                   const DoFHandler<deal_II_dimension, deal_II_space_dimension>
+                     &                            dof_handler,
+                   const unsigned int             component,
+                   const typename VEC::value_type constant_adjustment);
     \}
 #endif
   }

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.