]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Instantiate FETools::interpolate for general vector types.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 6 May 2013 22:18:42 +0000 (22:18 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 6 May 2013 22:18:42 +0000 (22:18 +0000)
git-svn-id: https://svn.dealii.org/trunk@29468 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/fe/fe_tools.h
deal.II/source/fe/fe_tools.inst.in

index 572d9d1cb93b6c845dd996b74a07c4e4cbd837f1..6c1458c3df989a22b16f999cd194b191cfcb2db1 100644 (file)
@@ -117,6 +117,12 @@ this function.
 
 <ol>
 
+<li> Changed: FETools::interpolate is now instantiated for all
+vector types, not just dealii::Vector and dealii::BlockVector.
+<br>
+(Wolfgang Bangerth, 2013/05/06)
+</li>
+
 <li> Fixed: Generate an error if the users tries to refine a cell
 that is already on the maximum level in a distributed triangulation.
 <br>
index 926eaa6e508399c4a2c4249cb0dd78a4f9a4da2e..38ac6ecf1ba933ecb1785160a1805af63040f4cf 100644 (file)
@@ -546,7 +546,7 @@ namespace FETools
    * course Q1 on each cell.
    *
    * For this case (continuous elements on grids with hanging nodes), please
-   * use the @p interpolate function with an additional @p ConstraintMatrix
+   * use the @p interpolate() function with an additional ConstraintMatrix
    * argument, see below, or make the field conforming yourself by calling the
    * @p distribute function of your hanging node constraints object.
    */
@@ -573,8 +573,8 @@ namespace FETools
    * interpolation. The same is true if @p fe1 is a continuous and @p fe2 is a
    * discontinuous finite element. For the case that @p fe1 is a discontinuous
    * and @p fe2 is a continuous finite element there is no point interpolation
-   * defined at the discontinuities.  Therefore the meanvalue is taken at the
-   * DoF values on the discontinuities.
+   * defined at the discontinuities.  Therefore the mean value is taken at the
+   * DoF values at the discontinuities.
    */
   template <int dim, int spacedim,
            template <int, int> class DH1,
index ab7eafa1838a97813c1f914d9ff1a347ccaf1368..22ee6025b016be0a5082706e7a2b830182acbf04 100644 (file)
 //---------------------------------------------------------------------------
 
 
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS; Vector : SERIAL_VECTORS)
+  {
+    namespace FETools
+      \{
+#if deal_II_dimension <= deal_II_space_dimension
+      template
+       void interpolate<deal_II_dimension,deal_II_space_dimension>
+       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector &,
+        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, Vector &);
+
+      template
+       void interpolate<deal_II_dimension,deal_II_space_dimension>
+       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector &,
+        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const ConstraintMatrix &,
+        Vector &);
+#endif
+      \}
+  }
+
 for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
   {
     namespace FETools
@@ -30,52 +49,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
       template
        void compute_embedding_matrices<deal_II_dimension, double, deal_II_space_dimension>
        (const FiniteElement<deal_II_dimension,deal_II_space_dimension> &,
-        std::vector<std::vector<FullMatrix<double> > >&,bool);
-      
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector<double> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, Vector<double> &);
-
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector<double> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const ConstraintMatrix &,
-        Vector<double> &);
-
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector<float> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, Vector<float> &);
-      template
-       
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const Vector<float> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const ConstraintMatrix &,
-        Vector<float> &);
-
-
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const BlockVector<double> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, BlockVector<double> &);
-      
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const BlockVector<double> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const ConstraintMatrix &,
-        BlockVector<double> &);
-
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const BlockVector<float> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, BlockVector<float> &);
-
-      template
-       void interpolate<deal_II_dimension,deal_II_space_dimension>
-       (const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const BlockVector<float> &,
-        const DoFHandler<deal_II_dimension,deal_II_space_dimension> &, const ConstraintMatrix &,
-        BlockVector<float> &);
+        std::vector<std::vector<FullMatrix<double> > >&,bool);      
 #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.