]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add std::complex make_periodicity_constraints()
authorDaniel Garcia-Sanchez <daniel.garcia-sanchez@insp.upmc.fr>
Fri, 2 Aug 2019 08:39:51 +0000 (10:39 +0200)
committerDaniel Garcia-Sanchez <daniel.garcia-sanchez@insp.upmc.fr>
Fri, 2 Aug 2019 08:39:51 +0000 (10:39 +0200)
include/deal.II/dofs/dof_tools.h
source/dofs/dof_tools_constraints.cc
source/dofs/dof_tools_constraints.inst.in

index 28bde4e4ab97ee56ff5ce1deba9e75868206edc7..104546738513c11ac8989ab9f0c46c2040f54d1d 100644 (file)
@@ -1153,12 +1153,12 @@ namespace DoFTools
    *
    * @author Matthias Maier, 2012 - 2015
    */
-  template <typename FaceIterator>
+  template <typename FaceIterator, typename number>
   void
   make_periodicity_constraints(
     const FaceIterator &                         face_1,
     const typename identity<FaceIterator>::type &face_2,
-    AffineConstraints<double> &                  constraints,
+    AffineConstraints<number> &                  constraints,
     const ComponentMask &            component_mask   = ComponentMask(),
     const bool                       face_orientation = true,
     const bool                       face_flip        = false,
@@ -1190,13 +1190,13 @@ namespace DoFTools
    *
    * @author Daniel Arndt, Matthias Maier, 2013 - 2015
    */
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(
     const std::vector<
       GridTools::PeriodicFacePair<typename DoFHandlerType::cell_iterator>>
       &                              periodic_faces,
-    AffineConstraints<double> &      constraints,
+    AffineConstraints<number> &      constraints,
     const ComponentMask &            component_mask = ComponentMask(),
     const std::vector<unsigned int> &first_vector_components =
       std::vector<unsigned int>());
@@ -1233,14 +1233,14 @@ namespace DoFTools
    *
    * @author Matthias Maier, 2012
    */
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(
     const DoFHandlerType &     dof_handler,
     const types::boundary_id   b_id1,
     const types::boundary_id   b_id2,
     const int                  direction,
-    AffineConstraints<double> &constraints,
+    AffineConstraints<number> &constraints,
     const ComponentMask &      component_mask = ComponentMask());
 
 
@@ -1269,13 +1269,13 @@ namespace DoFTools
    * @ref GlossPeriodicConstraints "Glossary entry on periodic boundary conditions"
    * for further information.
    */
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(
     const DoFHandlerType &     dof_handler,
     const types::boundary_id   b_id,
     const int                  direction,
-    AffineConstraints<double> &constraints,
+    AffineConstraints<number> &constraints,
     const ComponentMask &      component_mask = ComponentMask());
 
   /**
index 44a85b5545e9e8a780121136e791e9559e94c74e..3a9e3f65b4f85bef2d1fec1e69c4c2a454f39ac3 100644 (file)
@@ -1817,13 +1817,13 @@ namespace DoFTools
      * such "identity constraints" if the opposite constraint already
      * exists.
      */
-    template <typename FaceIterator>
+    template <typename FaceIterator, typename number>
     void
     set_periodicity_constraints(
       const FaceIterator &                         face_1,
       const typename identity<FaceIterator>::type &face_2,
       const FullMatrix<double> &                   transformation,
-      AffineConstraints<double> &                  affine_constraints,
+      AffineConstraints<number> &                  affine_constraints,
       const ComponentMask &                        component_mask,
       const bool                                   face_orientation,
       const bool                                   face_flip,
@@ -2073,7 +2073,7 @@ namespace DoFTools
           // a dependency cycle:
 
           bool   constraints_are_cyclic = true;
-          double cycle_factor           = factor;
+          number cycle_factor           = factor;
 
           for (auto test_dof = dof_right; test_dof != dof_left;)
             {
@@ -2220,12 +2220,12 @@ namespace DoFTools
   // Low level interface:
 
 
-  template <typename FaceIterator>
+  template <typename FaceIterator, typename number>
   void
   make_periodicity_constraints(
     const FaceIterator &                         face_1,
     const typename identity<FaceIterator>::type &face_2,
-    AffineConstraints<double> &                  affine_constraints,
+    AffineConstraints<number> &                  affine_constraints,
     const ComponentMask &                        component_mask,
     const bool                                   face_orientation,
     const bool                                   face_flip,
@@ -2449,13 +2449,13 @@ namespace DoFTools
 
 
 
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(
     const std::vector<
       GridTools::PeriodicFacePair<typename DoFHandlerType::cell_iterator>>
       &                              periodic_faces,
-    AffineConstraints<double> &      constraints,
+    AffineConstraints<number> &      constraints,
     const ComponentMask &            component_mask,
     const std::vector<unsigned int> &first_vector_components)
   {
@@ -2489,13 +2489,13 @@ namespace DoFTools
   // High level interface variants:
 
 
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(const DoFHandlerType &             dof_handler,
                                const types::boundary_id           b_id1,
                                const types::boundary_id           b_id2,
                                const int                          direction,
-                               dealii::AffineConstraints<double> &constraints,
+                               dealii::AffineConstraints<number> &constraints,
                                const ComponentMask &component_mask)
   {
     static const int space_dim = DoFHandlerType::space_dimension;
@@ -2522,12 +2522,12 @@ namespace DoFTools
 
 
 
-  template <typename DoFHandlerType>
+  template <typename DoFHandlerType, typename number>
   void
   make_periodicity_constraints(const DoFHandlerType &     dof_handler,
                                const types::boundary_id   b_id,
                                const int                  direction,
-                               AffineConstraints<double> &constraints,
+                               AffineConstraints<number> &constraints,
                                const ComponentMask &      component_mask)
   {
     static const int dim       = DoFHandlerType::dimension;
index 931ef21f2663074317e738c911b3bf8b94f0d904..82c18c6f166b18b6e282bd4ad7a8356901b8874c 100644 (file)
@@ -60,13 +60,36 @@ for (DH : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS)
       const FullMatrix<double> &,
       const std::vector<unsigned int> &);
 
-    template void DoFTools::make_periodicity_constraints<DH<deal_II_dimension>>(
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+    template void DoFTools::make_periodicity_constraints(
+      const DH<deal_II_dimension>::face_iterator &,
+      const DH<deal_II_dimension>::face_iterator &,
+      AffineConstraints<std::complex<double>> &,
+      const ComponentMask &,
+      bool,
+      bool,
+      bool,
+      const FullMatrix<double> &,
+      const std::vector<unsigned int> &);
+#endif
+
+    template void
+    DoFTools::make_periodicity_constraints<DH<deal_II_dimension>, double>(
       const std::vector<
         GridTools::PeriodicFacePair<DH<deal_II_dimension>::cell_iterator>> &,
       AffineConstraints<double> &,
       const ComponentMask &,
       const std::vector<unsigned int> &);
 
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+    template void DoFTools::make_periodicity_constraints<DH<deal_II_dimension>,
+                                                         std::complex<double>>(
+      const std::vector<
+        GridTools::PeriodicFacePair<DH<deal_II_dimension>::cell_iterator>> &,
+      AffineConstraints<std::complex<double>> &,
+      const ComponentMask &,
+      const std::vector<unsigned int> &);
+#endif
 
     template void DoFTools::make_periodicity_constraints(
       const DH<deal_II_dimension> &,
@@ -76,12 +99,31 @@ for (DH : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS)
       AffineConstraints<double> &,
       const ComponentMask &);
 
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+    template void DoFTools::make_periodicity_constraints(
+      const DH<deal_II_dimension> &,
+      const types::boundary_id,
+      const types::boundary_id,
+      const int,
+      AffineConstraints<std::complex<double>> &,
+      const ComponentMask &);
+#endif
+
     template void DoFTools::make_periodicity_constraints(
       const DH<deal_II_dimension> &,
       const types::boundary_id,
       const int,
       AffineConstraints<double> &,
       const ComponentMask &);
+
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+    template void DoFTools::make_periodicity_constraints(
+      const DH<deal_II_dimension> &,
+      const types::boundary_id,
+      const int,
+      AffineConstraints<std::complex<double>> &,
+      const ComponentMask &);
+#endif
   }
 
 for (DH : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS;

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.