]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add missing versions of VectorTools::interpolate_boundary_values for p::MappingCollection 11367/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 13 Dec 2020 14:45:09 +0000 (15:45 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 13 Dec 2020 14:45:09 +0000 (15:45 +0100)
include/deal.II/numerics/vector_tools_boundary.h
include/deal.II/numerics/vector_tools_boundary.templates.h
source/numerics/vector_tools_boundary.inst.in

index 69cd87a0f7846e072e536071f376bec7a867e44e..6593be957bd2c94509873e62dd2e436e9b86c361 100644 (file)
@@ -139,6 +139,20 @@ namespace VectorTools
     std::map<types::global_dof_index, number> &boundary_values,
     const ComponentMask &component_mask = ComponentMask());
 
+  /**
+   * Like the previous function, but take a mapping collection to go with
+   * DoFHandler objects with hp-capabilities.
+   */
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const types::boundary_id                    boundary_component,
+    const Function<spacedim, number> &          boundary_function,
+    std::map<types::global_dof_index, number> & boundary_values,
+    const ComponentMask &component_mask = ComponentMask());
+
   /**
    * Call the other interpolate_boundary_values() function, see above, with
    * <tt>mapping=MappingQGeneric@<dim,spacedim@>(1)</tt>. The same comments
@@ -245,6 +259,20 @@ namespace VectorTools
     AffineConstraints<number> &constraints,
     const ComponentMask &      component_mask = ComponentMask());
 
+  /**
+   * Like the previous function, but take a mapping collection to go with
+   * DoFHandler objects with hp-capabilities.
+   */
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const std::map<types::boundary_id, const Function<spacedim, number> *>
+      &                        function_map,
+    AffineConstraints<number> &constraints,
+    const ComponentMask &      component_mask = ComponentMask());
+
   /**
    * Same function as above, but taking only one pair of boundary indicator
    * and corresponding boundary function. The same comments apply as for the
@@ -266,6 +294,20 @@ namespace VectorTools
     AffineConstraints<number> &       constraints,
     const ComponentMask &             component_mask = ComponentMask());
 
+  /**
+   * Like the previous function, but take a mapping collection to go with
+   * DoFHandler objects with hp-capabilities.
+   */
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const types::boundary_id                    boundary_component,
+    const Function<spacedim, number> &          boundary_function,
+    AffineConstraints<number> &                 constraints,
+    const ComponentMask &component_mask = ComponentMask());
+
   /**
    * Call the other interpolate_boundary_values() function, see above, with
    * <tt>mapping=MappingQGeneric@<dim,spacedim@>(1)</tt>. The same comments
index 7801f0d42e5be9f8ba9ff43ff2045f632485a3ce..cacf2894433870b4f05f23e3b30737d29f6ad497 100644 (file)
@@ -440,6 +440,25 @@ namespace VectorTools
 
 
 
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const types::boundary_id                    boundary_component,
+    const Function<spacedim, number> &          boundary_function,
+    std::map<types::global_dof_index, number> & boundary_values,
+    const ComponentMask &                       component_mask)
+  {
+    std::map<types::boundary_id, const Function<spacedim, number> *>
+      function_map;
+    function_map[boundary_component] = &boundary_function;
+    interpolate_boundary_values(
+      mapping, dof, function_map, boundary_values, component_mask);
+  }
+
+
+
   template <int dim, int spacedim, typename number>
   void
   interpolate_boundary_values(
@@ -530,6 +549,54 @@ namespace VectorTools
 
 
 
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const std::map<types::boundary_id, const Function<spacedim, number> *>
+      &                        function_map,
+    AffineConstraints<number> &constraints,
+    const ComponentMask &      component_mask_)
+  {
+    std::map<types::global_dof_index, number> boundary_values;
+    interpolate_boundary_values(
+      mapping, dof, function_map, boundary_values, component_mask_);
+    typename std::map<types::global_dof_index, number>::const_iterator
+      boundary_value = boundary_values.begin();
+    for (; boundary_value != boundary_values.end(); ++boundary_value)
+      {
+        if (constraints.can_store_line(boundary_value->first) &&
+            !constraints.is_constrained(boundary_value->first))
+          {
+            constraints.add_line(boundary_value->first);
+            constraints.set_inhomogeneity(boundary_value->first,
+                                          boundary_value->second);
+          }
+      }
+  }
+
+
+
+  template <int dim, int spacedim, typename number>
+  void
+  interpolate_boundary_values(
+    const hp::MappingCollection<dim, spacedim> &mapping,
+    const DoFHandler<dim, spacedim> &           dof,
+    const types::boundary_id                    boundary_component,
+    const Function<spacedim, number> &          boundary_function,
+    AffineConstraints<number> &                 constraints,
+    const ComponentMask &                       component_mask)
+  {
+    std::map<types::boundary_id, const Function<spacedim, number> *>
+      function_map;
+    function_map[boundary_component] = &boundary_function;
+    interpolate_boundary_values(
+      mapping, dof, function_map, constraints, component_mask);
+  }
+
+
+
   template <int dim, int spacedim, typename number>
   void
   interpolate_boundary_values(
index 468ceeb8055ff090bacdd9b54adab8a9d67484cc..f2e24d69f3a6fa4d81ffd49d9f6b33bb0938f393 100644 (file)
@@ -40,6 +40,18 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
 
       template void
       interpolate_boundary_values(
+        const hp::MappingCollection<deal_II_dimension, deal_II_space_dimension>
+          &,
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
+        const std::map<types::boundary_id,
+                       const Function<deal_II_space_dimension, number> *> &,
+        std::map<types::global_dof_index, number> &,
+        const ComponentMask &);
+
+      template void
+      interpolate_boundary_values(
+        const hp::MappingCollection<deal_II_dimension, deal_II_space_dimension>
+          &,
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
         const types::boundary_id,
         const Function<deal_II_space_dimension, number> &,
@@ -49,25 +61,13 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
       template void
       interpolate_boundary_values(
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
-        const std::map<types::boundary_id,
-                       const Function<deal_II_space_dimension, number> *> &,
+        const types::boundary_id,
+        const Function<deal_II_space_dimension, number> &,
         std::map<types::global_dof_index, number> &,
         const ComponentMask &);
-    \}
-#endif
-  }
-
 
-for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
-     number : REAL_AND_COMPLEX_SCALARS)
-  {
-#if deal_II_dimension <= deal_II_space_dimension
-    namespace VectorTools
-    \{
       template void
       interpolate_boundary_values(
-        const hp::MappingCollection<deal_II_dimension, deal_II_space_dimension>
-          &,
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
         const std::map<types::boundary_id,
                        const Function<deal_II_space_dimension, number> *> &,
@@ -104,6 +104,26 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
         AffineConstraints<number> &,
         const ComponentMask &);
 
+      template void
+      interpolate_boundary_values(
+        const hp::MappingCollection<deal_II_dimension, deal_II_space_dimension>
+          &,
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
+        const std::map<types::boundary_id,
+                       const Function<deal_II_space_dimension, number> *> &,
+        AffineConstraints<number> &,
+        const ComponentMask &);
+
+      template void
+      interpolate_boundary_values(
+        const hp::MappingCollection<deal_II_dimension, deal_II_space_dimension>
+          &,
+        const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
+        const types::boundary_id,
+        const Function<deal_II_space_dimension, number> &,
+        AffineConstraints<number> &,
+        const ComponentMask &);
+
       template void
       interpolate_boundary_values(
         const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,

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.