]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add typename and instantiation 13907/head
authorNiklas Wik <niiklaswiik@gmail.com>
Sat, 4 Jun 2022 15:10:36 +0000 (17:10 +0200)
committerNiklas Wik <niiklaswiik@gmail.com>
Sat, 4 Jun 2022 15:41:17 +0000 (17:41 +0200)
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 55d07a5ef1a707d6877cf7c69633d0ed3a4a061b..0d139e1b25ee38cb923888c29fd4e1454ab27955 100644 (file)
@@ -718,15 +718,15 @@ namespace VectorTools
    * @see
    * @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
    */
-  template <int dim>
+  template <int dim, typename number, typename number2 = number>
   void
   project_boundary_values_div_conforming(
-    const DoFHandler<dim, dim> & dof_handler,
-    const unsigned int           first_vector_component,
-    const Function<dim, double> &boundary_function,
-    const types::boundary_id     boundary_component,
-    AffineConstraints<double> &  constraints,
-    const Mapping<dim> &         mapping);
+    const DoFHandler<dim, dim> &  dof_handler,
+    const unsigned int            first_vector_component,
+    const Function<dim, number2> &boundary_function,
+    const types::boundary_id      boundary_component,
+    AffineConstraints<number> &   constraints,
+    const Mapping<dim> &          mapping);
 
   /**
    * Same as above for the hp-namespace.
@@ -736,14 +736,14 @@ namespace VectorTools
    * @see
    * @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
    */
-  template <int dim>
+  template <int dim, typename number, typename number2 = number>
   void
   project_boundary_values_div_conforming(
     const DoFHandler<dim, dim> &           dof_handler,
     const unsigned int                     first_vector_component,
-    const Function<dim, double> &          boundary_function,
+    const Function<dim, number2> &         boundary_function,
     const types::boundary_id               boundary_component,
-    AffineConstraints<double> &            constraints,
+    AffineConstraints<number> &            constraints,
     const hp::MappingCollection<dim, dim> &mapping_collection =
       hp::StaticMappingQ1<dim>::mapping_collection);
 
index ec71f01e8639a9bd1248bc0b8233ea0b8ab2ef48..b84f413ac0a75f7ca6078a656d89325bc3a23fe1 100644 (file)
@@ -2149,16 +2149,16 @@ namespace VectorTools
   {
     // This function computes the projection of the boundary function on the
     // boundary in 2d.
-    template <typename cell_iterator>
+    template <typename cell_iterator, typename number, typename number2>
     void
     compute_face_projection_div_conforming(
       const cell_iterator &                       cell,
       const unsigned int                          face,
       const FEFaceValues<2> &                     fe_values,
       const unsigned int                          first_vector_component,
-      const Function<2> &                         boundary_function,
+      const Function<2, number2> &                boundary_function,
       const std::vector<DerivativeForm<1, 2, 2>> &jacobians,
-      AffineConstraints<double> &                 constraints)
+      AffineConstraints<number> &                 constraints)
     {
       // Compute the integral over the product of the normal components of
       // the boundary function times the normal components of the shape
@@ -2171,9 +2171,9 @@ namespace VectorTools
                                                                       1,
                                                                       0,
                                                                       0};
-      std::vector<Vector<double>> values(fe_values.n_quadrature_points,
-                                         Vector<double>(2));
-      Vector<double>              dof_values(fe.n_dofs_per_face(face));
+      std::vector<Vector<number2>> values(fe_values.n_quadrature_points,
+                                          Vector<number2>(2));
+      Vector<number2>              dof_values(fe.n_dofs_per_face(face));
 
       // Get the values of the boundary function at the quadrature points.
       {
@@ -2186,7 +2186,7 @@ namespace VectorTools
       for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points;
            ++q_point)
         {
-          double tmp = 0.0;
+          number2 tmp = 0.0;
 
           for (unsigned int d = 0; d < 2; ++d)
             tmp += normals[q_point][d] * values[q_point](d);
@@ -2235,32 +2235,35 @@ namespace VectorTools
     }
 
     // dummy implementation of above function for all other dimensions
-    template <int dim, typename cell_iterator>
+    template <int dim,
+              typename cell_iterator,
+              typename number,
+              typename number2>
     void
     compute_face_projection_div_conforming(
       const cell_iterator &,
       const unsigned int,
       const FEFaceValues<dim> &,
       const unsigned int,
-      const Function<dim> &,
+      const Function<dim, number2> &,
       const std::vector<DerivativeForm<1, dim, dim>> &,
-      AffineConstraints<double> &)
+      AffineConstraints<number> &)
     {
       Assert(false, ExcNotImplemented());
     }
 
     // This function computes the projection of the boundary function on the
     // boundary in 3d.
-    template <typename cell_iterator>
+    template <typename cell_iterator, typename number, typename number2>
     void
     compute_face_projection_div_conforming(
       const cell_iterator &                       cell,
       const unsigned int                          face,
       const FEFaceValues<3> &                     fe_values,
       const unsigned int                          first_vector_component,
-      const Function<3> &                         boundary_function,
+      const Function<3, number2> &                boundary_function,
       const std::vector<DerivativeForm<1, 3, 3>> &jacobians,
-      std::vector<double> &                       dof_values,
+      std::vector<number> &                       dof_values,
       std::vector<types::global_dof_index> &      projected_dofs)
     {
       // Compute the intergral over the product of the normal components of
@@ -2272,9 +2275,9 @@ namespace VectorTools
       const unsigned int
         face_coordinate_directions[GeometryInfo<3>::faces_per_cell][2] = {
           {1, 2}, {1, 2}, {2, 0}, {2, 0}, {0, 1}, {0, 1}};
-      std::vector<Vector<double>> values(fe_values.n_quadrature_points,
-                                         Vector<double>(3));
-      Vector<double>              dof_values_local(fe.n_dofs_per_face(face));
+      std::vector<Vector<number2>> values(fe_values.n_quadrature_points,
+                                          Vector<number2>(3));
+      Vector<number2>              dof_values_local(fe.n_dofs_per_face(face));
 
       {
         const std::vector<Point<3>> &quadrature_points =
@@ -2286,7 +2289,7 @@ namespace VectorTools
       for (unsigned int q_point = 0; q_point < fe_values.n_quadrature_points;
            ++q_point)
         {
-          double tmp = 0.0;
+          number2 tmp = 0.0;
 
           for (unsigned int d = 0; d < 3; ++d)
             tmp += normals[q_point][d] * values[q_point](d);
@@ -2342,16 +2345,19 @@ namespace VectorTools
     // dummy implementation of above
     // function for all other
     // dimensions
-    template <int dim, typename cell_iterator>
+    template <int dim,
+              typename cell_iterator,
+              typename number,
+              typename number2>
     void
     compute_face_projection_div_conforming(
       const cell_iterator &,
       const unsigned int,
       const FEFaceValues<dim> &,
       const unsigned int,
-      const Function<dim> &,
+      const Function<dim, number2> &,
       const std::vector<DerivativeForm<1, dim, dim>> &,
-      std::vector<double> &,
+      std::vector<number> &,
       std::vector<types::global_dof_index> &)
     {
       Assert(false, ExcNotImplemented());
@@ -2359,15 +2365,15 @@ namespace VectorTools
   } // namespace internals
 
 
-  template <int dim>
+  template <int dim, typename number, typename number2>
   void
   project_boundary_values_div_conforming(
-    const DoFHandler<dim> &    dof_handler,
-    const unsigned int         first_vector_component,
-    const Function<dim> &      boundary_function,
-    const types::boundary_id   boundary_component,
-    AffineConstraints<double> &constraints,
-    const Mapping<dim> &       mapping)
+    const DoFHandler<dim> &       dof_handler,
+    const unsigned int            first_vector_component,
+    const Function<dim, number2> &boundary_function,
+    const types::boundary_id      boundary_component,
+    AffineConstraints<number> &   constraints,
+    const Mapping<dim> &          mapping)
   {
     const unsigned int spacedim = dim;
     // Interpolate the normal components
@@ -2430,7 +2436,9 @@ namespace VectorTools
                         {
                           AssertThrow(
                             dynamic_cast<const FE_RaviartThomas<dim> *>(&fe) !=
-                              nullptr,
+                                nullptr ||
+                              dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
+                                &fe) != nullptr,
                             typename FiniteElement<
                               dim>::ExcInterpolationNotImplemented());
                         }
@@ -2485,7 +2493,9 @@ namespace VectorTools
                         {
                           AssertThrow(
                             dynamic_cast<const FE_RaviartThomas<dim> *>(&fe) !=
-                              nullptr,
+                                nullptr ||
+                              dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
+                                &fe) != nullptr,
                             typename FiniteElement<
                               dim>::ExcInterpolationNotImplemented());
                         }
@@ -2529,14 +2539,14 @@ namespace VectorTools
   }
 
 
-  template <int dim>
+  template <int dim, typename number, typename number2>
   void
   project_boundary_values_div_conforming(
     const DoFHandler<dim> &                dof_handler,
     const unsigned int                     first_vector_component,
-    const Function<dim> &                  boundary_function,
+    const Function<dim, number2> &         boundary_function,
     const types::boundary_id               boundary_component,
-    AffineConstraints<double> &            constraints,
+    AffineConstraints<number> &            constraints,
     const hp::MappingCollection<dim, dim> &mapping_collection)
   {
     const unsigned int           spacedim = dim;
@@ -2618,7 +2628,7 @@ namespace VectorTools
         case 3:
           {
             const unsigned int                   n_dofs = dof_handler.n_dofs();
-            std::vector<double>                  dof_values(n_dofs);
+            std::vector<number2>                 dof_values(n_dofs);
             std::vector<types::global_dof_index> projected_dofs(n_dofs);
 
             for (unsigned int dof = 0; dof < n_dofs; ++dof)
index 18013351270a94afc69a4448994cf558853276b3..0604787d394f4318062ade6551223f04e22edf4e 100644 (file)
@@ -268,21 +268,37 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
         AffineConstraints<std::complex<double>> &,
         const hp::MappingCollection<deal_II_dimension> &);
 #    endif
+#  endif
+#endif
+    \}
+  }
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
+     number : REAL_SCALARS;
+     number2 : REAL_SCALARS)
+  {
+    namespace VectorTools
+    \{
+#if deal_II_dimension == deal_II_space_dimension
+
+#  if deal_II_dimension != 1
+
       template void
       project_boundary_values_div_conforming<deal_II_dimension>(
         const DoFHandler<deal_II_dimension> &,
         const unsigned int,
-        const Function<deal_II_dimension> &,
+        const Function<deal_II_dimension, number2> &,
         const types::boundary_id,
-        AffineConstraints<double> &,
+        AffineConstraints<number> &,
         const Mapping<deal_II_dimension> &);
+
       template void
       project_boundary_values_div_conforming<deal_II_dimension>(
         const DoFHandler<deal_II_dimension> &,
         const unsigned int,
-        const Function<deal_II_dimension> &,
+        const Function<deal_II_dimension, number2> &,
         const types::boundary_id,
-        AffineConstraints<double> &,
+        AffineConstraints<number> &,
         const hp::MappingCollection<deal_II_dimension> &);
 #  endif
 #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.