]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Match the template arguments of matrix and vector for MatrixCreator::create_mass_matr...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 17 Jun 2015 17:26:41 +0000 (12:26 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 17 Jun 2015 19:34:02 +0000 (14:34 -0500)
doc/news/changes.h
include/deal.II/numerics/matrix_tools.h
source/numerics/matrix_tools.cc
source/numerics/matrix_tools.inst.in

index e2f743a95d0bda5f3884db07ee5f0b0b4b778fc2..3bd537bfe5014a0c50ecb3a5d08ec16b7e035113 100644 (file)
@@ -38,6 +38,14 @@ inconvenience this causes.
 </p>
 
 <ol>
+  <li>Changed: MatrixCreator::create_mass_matrix() took matrix and vector
+  objects where the scalar type of the matrix was a template argument
+  but the scalar type of the vector was always <code>double</code>. This
+  has been changed so that the two always need to match.
+  <br>
+  (David Davydov, Wolfgang Bangerth, 2015/06/17)
+  </li>
+
   <li>Changed: Functions such as FEValuesBase::get_function_values() that
   extracted the values of functions at the quadrature points of a cell
   implicitly always assumed that <code>double</code> is a reasonable
index b72cc11c8157b3f3919e5b0c2721d75e8b4269b0..9daa1ebd2c3ae1c28b389b6ca56be6f8fcbb7ac8 100644 (file)
@@ -273,7 +273,7 @@ namespace MatrixCreator
                            const Quadrature<dim>    &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim> &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const a = 0,
                            const ConstraintMatrix   &constraints = ConstraintMatrix());
 
@@ -286,7 +286,7 @@ namespace MatrixCreator
                            const Quadrature<dim>    &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim> &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const a = 0,
                            const ConstraintMatrix   &constraints = ConstraintMatrix());
 
@@ -320,7 +320,7 @@ namespace MatrixCreator
                            const hp::QCollection<dim> &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim> &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const a = 0,
                            const ConstraintMatrix   &constraints = ConstraintMatrix());
 
@@ -332,7 +332,7 @@ namespace MatrixCreator
                            const hp::QCollection<dim> &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim> &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const a = 0,
                            const ConstraintMatrix   &constraints = ConstraintMatrix());
 
@@ -360,7 +360,6 @@ namespace MatrixCreator
    * shape functions.
    */
   template <int dim, int spacedim>
-
   void create_boundary_mass_matrix (const Mapping<dim, spacedim>       &mapping,
                                     const DoFHandler<dim,spacedim>    &dof,
                                     const Quadrature<dim-1>  &q,
@@ -400,7 +399,6 @@ namespace MatrixCreator
    * <tt>mapping=MappingQ1@<dim@>()</tt>.
    */
   template <int dim, int spacedim>
-
   void create_boundary_mass_matrix (const DoFHandler<dim,spacedim>    &dof,
                                     const Quadrature<dim-1>  &q,
                                     SparseMatrix<double>     &matrix,
@@ -414,7 +412,6 @@ namespace MatrixCreator
    * Same function as above, but for hp objects.
    */
   template <int dim, int spacedim>
-
   void create_boundary_mass_matrix (const hp::MappingCollection<dim,spacedim>       &mapping,
                                     const hp::DoFHandler<dim,spacedim>    &dof,
                                     const hp::QCollection<dim-1>  &q,
@@ -442,7 +439,6 @@ namespace MatrixCreator
    * Same function as above, but for hp objects.
    */
   template <int dim, int spacedim>
-
   void create_boundary_mass_matrix (const hp::DoFHandler<dim,spacedim>    &dof,
                                     const hp::QCollection<dim-1>  &q,
                                     SparseMatrix<double>     &matrix,
index 8a385ed620e97f1b06e889d6fe8378d9662d930f..ef56f60fdbfb254bbb80f9dfe27242eddb18c85c 100644 (file)
@@ -201,11 +201,12 @@ namespace MatrixCreator
       };
 
 
+      template <typename number>
       struct CopyData
       {
         std::vector<types::global_dof_index> dof_indices;
-        FullMatrix<double>        cell_matrix;
-        dealii::Vector<double>    cell_rhs;
+        FullMatrix<number>        cell_matrix;
+        dealii::Vector<number>    cell_rhs;
         const ConstraintMatrix   *constraints;
       };
     }
@@ -213,10 +214,11 @@ namespace MatrixCreator
 
     template <int dim,
               int spacedim,
-              typename CellIterator>
+              typename CellIterator,
+              typename number>
     void mass_assembler (const CellIterator &cell,
                          MatrixCreator::internal::AssemblerData::Scratch<dim,spacedim> &data,
-                         MatrixCreator::internal::AssemblerData::CopyData              &copy_data)
+                         MatrixCreator::internal::AssemblerData::CopyData<number>      &copy_data)
     {
       data.x_fe_values.reinit (cell);
       const FEValues<dim,spacedim> &fe_values = data.x_fe_values.get_present_fe_values ();
@@ -397,7 +399,7 @@ namespace MatrixCreator
               typename CellIterator>
     void laplace_assembler (const CellIterator &cell,
                             MatrixCreator::internal::AssemblerData::Scratch<dim,spacedim> &data,
-                            MatrixCreator::internal::AssemblerData::CopyData              &copy_data)
+                            MatrixCreator::internal::AssemblerData::CopyData<double>      &copy_data)
     {
       data.x_fe_values.reinit (cell);
       const FEValues<dim,spacedim> &fe_values = data.x_fe_values.get_present_fe_values ();
@@ -573,9 +575,10 @@ namespace MatrixCreator
 
 
 
-    template <typename MatrixType,
+    template <typename number,
+              typename MatrixType,
               typename VectorType>
-    void copy_local_to_global (const AssemblerData::CopyData &data,
+    void copy_local_to_global (const AssemblerData::CopyData<number> &data,
                                MatrixType *matrix,
                                VectorType *right_hand_side)
     {
@@ -666,7 +669,7 @@ namespace MatrixCreator
                     coefficient, /*rhs_function=*/0,
                     q_collection, mapping_collection);
 
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<number> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -675,10 +678,10 @@ namespace MatrixCreator
 
     WorkStream::run (dof.begin_active(),
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
-                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
+                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator,number>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<number>, Vector<double> >,
-                                      std_cxx11::_1, &matrix, (Vector<double> *)0),
+                                      copy_local_to_global<number,SparseMatrix<number>, Vector<number> >,
+                                      std_cxx11::_1, &matrix, (Vector<number> *)0),
                      assembler_data,
                      copy_data);
   }
@@ -704,7 +707,7 @@ namespace MatrixCreator
                            const Quadrature<dim>    &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim>      &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const coefficient,
                            const ConstraintMatrix   &constraints)
   {
@@ -722,7 +725,7 @@ namespace MatrixCreator
                     update_JxW_values | update_quadrature_points,
                     coefficient, &rhs,
                     q_collection, mapping_collection);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<number> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -731,9 +734,9 @@ namespace MatrixCreator
 
     WorkStream::run (dof.begin_active(),
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
-                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
+                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator,number>,
                      std_cxx11::bind(&MatrixCreator::internal::
-                                     copy_local_to_global<SparseMatrix<number>, Vector<double> >,
+                                     copy_local_to_global<number,SparseMatrix<number>, Vector<number> >,
                                      std_cxx11::_1, &matrix, &rhs_vector),
                      assembler_data,
                      copy_data);
@@ -746,7 +749,7 @@ namespace MatrixCreator
                            const Quadrature<dim>    &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim>      &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const coefficient,
                            const ConstraintMatrix   &constraints)
   {
@@ -776,7 +779,7 @@ namespace MatrixCreator
                     (coefficient != 0 ? update_quadrature_points : UpdateFlags(0)),
                     coefficient, /*rhs_function=*/0,
                     q, mapping);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<number> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -785,10 +788,10 @@ namespace MatrixCreator
 
     WorkStream::run (dof.begin_active(),
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
-                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
+                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator,number>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<number>, Vector<double> >,
-                                      std_cxx11::_1, &matrix, (Vector<double> *)0),
+                                      copy_local_to_global<number,SparseMatrix<number>, Vector<number> >,
+                                      std_cxx11::_1, &matrix, (Vector<number> *)0),
                      assembler_data,
                      copy_data);
   }
@@ -814,7 +817,7 @@ namespace MatrixCreator
                            const hp::QCollection<dim>    &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim>      &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const coefficient,
                            const ConstraintMatrix   &constraints)
   {
@@ -829,7 +832,7 @@ namespace MatrixCreator
                     update_JxW_values | update_quadrature_points,
                     coefficient, &rhs,
                     q, mapping);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<number> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -838,9 +841,9 @@ namespace MatrixCreator
 
     WorkStream::run (dof.begin_active(),
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
-                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
+                     &MatrixCreator::internal::mass_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator,number>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<number>, Vector<double> >,
+                                      copy_local_to_global<number,SparseMatrix<number>, Vector<number> >,
                                       std_cxx11::_1, &matrix, &rhs_vector),
                      assembler_data,
                      copy_data);
@@ -853,7 +856,7 @@ namespace MatrixCreator
                            const hp::QCollection<dim> &q,
                            SparseMatrix<number>     &matrix,
                            const Function<spacedim> &rhs,
-                           Vector<double>           &rhs_vector,
+                           Vector<number>           &rhs_vector,
                            const Function<spacedim> *const coefficient,
                            const ConstraintMatrix   &constraints)
   {
@@ -1678,7 +1681,7 @@ namespace MatrixCreator
                     (coefficient != 0 ? update_quadrature_points : UpdateFlags(0)),
                     coefficient, /*rhs_function=*/0,
                     q_collection, mapping_collection);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<double> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -1689,7 +1692,7 @@ namespace MatrixCreator
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      copy_local_to_global<double,SparseMatrix<double>, Vector<double> >,
                                       std_cxx11::_1,
                                       &matrix,
                                       (Vector<double> *)NULL),
@@ -1735,7 +1738,7 @@ namespace MatrixCreator
                     update_JxW_values | update_quadrature_points,
                     coefficient, &rhs,
                     q_collection, mapping_collection);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<double> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -1746,7 +1749,7 @@ namespace MatrixCreator
                      static_cast<typename DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename DoFHandler<dim,spacedim>::active_cell_iterator>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      copy_local_to_global<double,SparseMatrix<double>, Vector<double> >,
                                       std_cxx11::_1,
                                       &matrix,
                                       &rhs_vector),
@@ -1790,7 +1793,7 @@ namespace MatrixCreator
                     (coefficient != 0 ? update_quadrature_points : UpdateFlags(0)),
                     coefficient, /*rhs_function=*/0,
                     q, mapping);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<double> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -1801,7 +1804,7 @@ namespace MatrixCreator
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      copy_local_to_global<double,SparseMatrix<double>, Vector<double> >,
                                       std_cxx11::_1,
                                       &matrix,
                                       (Vector<double> *)0),
@@ -1844,7 +1847,7 @@ namespace MatrixCreator
                     update_JxW_values | update_quadrature_points,
                     coefficient, &rhs,
                     q, mapping);
-    MatrixCreator::internal::AssemblerData::CopyData copy_data;
+    MatrixCreator::internal::AssemblerData::CopyData<double> copy_data;
     copy_data.cell_matrix.reinit (assembler_data.fe_collection.max_dofs_per_cell(),
                                   assembler_data.fe_collection.max_dofs_per_cell());
     copy_data.cell_rhs.reinit (assembler_data.fe_collection.max_dofs_per_cell());
@@ -1855,7 +1858,7 @@ namespace MatrixCreator
                      static_cast<typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>(dof.end()),
                      &MatrixCreator::internal::laplace_assembler<dim, spacedim, typename hp::DoFHandler<dim,spacedim>::active_cell_iterator>,
                      std_cxx11::bind (&MatrixCreator::internal::
-                                      copy_local_to_global<SparseMatrix<double>, Vector<double> >,
+                                      copy_local_to_global<double,SparseMatrix<double>, Vector<double> >,
                                       std_cxx11::_1,
                                       &matrix,
                                       &rhs_vector),
index 1a787b7336d9a41bad73c54dfa14e570507e70e9..a2e9346dc5bed08b1300eb1bfadd57afb3d9c9c2 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2010 - 2014 by the deal.II authors
+// Copyright (C) 2010 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -101,6 +101,42 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
        const Function<deal_II_space_dimension> * const,
        std::vector<unsigned int>);
 
+// same for float
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const Mapping<deal_II_dimension,deal_II_space_dimension>       &mapping,
+       const DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const Quadrature<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const Quadrature<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const Mapping<deal_II_dimension,deal_II_space_dimension>       &mapping,
+       const DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const Quadrature<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension>      &rhs,
+       Vector<float>           &rhs_vector,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const Quadrature<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension>      &rhs,
+       Vector<float>           &rhs_vector,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+
 #endif
   }
 
@@ -109,8 +145,8 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
 // where applicable and move to codimension cases above also when applicable
 for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
   {
-#if deal_II_dimension <= deal_II_space_dimension
 // hp versions of functions
+#if deal_II_dimension <= deal_II_space_dimension
     template
       void MatrixCreator::create_mass_matrix<deal_II_dimension,double,deal_II_space_dimension>
       (const hp::MappingCollection<deal_II_dimension,deal_II_space_dimension>       &mapping,
@@ -227,8 +263,56 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS
        const Function<deal_II_dimension> * const coefficient,
        const ConstraintMatrix   &constraints);
 
+#endif
+
+// same for float
+#if deal_II_dimension <= deal_II_space_dimension
+
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const hp::MappingCollection<deal_II_dimension,deal_II_space_dimension>       &mapping,
+       const hp::DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const hp::QCollection<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension,float,deal_II_space_dimension>
+      (const hp::MappingCollection<deal_II_dimension,deal_II_space_dimension>       &mapping,
+       const hp::DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const hp::QCollection<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension>      &rhs,
+       Vector<float>           &rhs_vector,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
 
 
 #endif
+
+#if deal_II_dimension == deal_II_space_dimension
+
+       
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension>
+      (const hp::DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const hp::QCollection<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+
+    template
+      void MatrixCreator::create_mass_matrix<deal_II_dimension>
+      (const hp::DoFHandler<deal_II_dimension,deal_II_space_dimension>    &dof,
+       const hp::QCollection<deal_II_dimension>    &q,
+       SparseMatrix<float>     &matrix,
+       const Function<deal_II_space_dimension>      &rhs,
+       Vector<float>           &rhs_vector,
+       const Function<deal_II_space_dimension> * const coefficient,
+       const ConstraintMatrix   &constraints);
+
+#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.