]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove DoFHandlerType from numerics/error_estimator.h 10574/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 22 Jun 2020 20:14:46 +0000 (22:14 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 22 Jun 2020 21:06:26 +0000 (23:06 +0200)
include/deal.II/numerics/error_estimator.h
include/deal.II/numerics/error_estimator.templates.h
source/numerics/error_estimator.inst.in
source/numerics/error_estimator_1d.cc
source/numerics/error_estimator_1d.inst.in

index e6312f5c74669fdff1b1eda7353fe95dde4fe4a9..a5a88a65858af64b96812b69095efdcc32bfd202 100644 (file)
@@ -33,6 +33,8 @@ DEAL_II_NAMESPACE_OPEN
 // Forward declarations
 #ifndef DOXYGEN
 template <int, int>
+class DoFHandler;
+template <int, int>
 class Mapping;
 template <int>
 class Quadrature;
@@ -335,12 +337,12 @@ public:
    * cell in the mesh as reported by
    * parallel::distributed::Triangulation::n_locally_owned_active_cells().
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<dim, spacedim> &mapping,
-    const DoFHandlerType &        dof,
-    const Quadrature<dim - 1> &   quadrature,
+    const Mapping<dim, spacedim> &   mapping,
+    const DoFHandler<dim, spacedim> &dof,
+    const Quadrature<dim - 1> &      quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -357,11 +359,11 @@ public:
    * Call the @p estimate function, see above, with
    * <tt>mapping=MappingQGeneric@<dim@>(1)</tt>.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &     dof,
-    const Quadrature<dim - 1> &quadrature,
+    const DoFHandler<dim, spacedim> &dof,
+    const Quadrature<dim - 1> &      quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -387,12 +389,12 @@ public:
    * construct of vector of references, so we had to use a vector of
    * pointers.)
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<dim, spacedim> &mapping,
-    const DoFHandlerType &        dof,
-    const Quadrature<dim - 1> &   quadrature,
+    const Mapping<dim, spacedim> &   mapping,
+    const DoFHandler<dim, spacedim> &dof,
+    const Quadrature<dim - 1> &      quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -409,11 +411,11 @@ public:
    * Call the @p estimate function, see above, with
    * <tt>mapping=MappingQGeneric@<dim@>(1)</tt>.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &     dof,
-    const Quadrature<dim - 1> &quadrature,
+    const DoFHandler<dim, spacedim> &dof,
+    const Quadrature<dim - 1> &      quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -431,12 +433,12 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<dim, spacedim> &  mapping,
-    const DoFHandlerType &          dof,
-    const hp::QCollection<dim - 1> &quadrature,
+    const Mapping<dim, spacedim> &   mapping,
+    const DoFHandler<dim, spacedim> &dof,
+    const hp::QCollection<dim - 1> & quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -454,11 +456,11 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &          dof,
-    const hp::QCollection<dim - 1> &quadrature,
+    const DoFHandler<dim, spacedim> &dof,
+    const hp::QCollection<dim - 1> & quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -476,12 +478,12 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<dim, spacedim> &  mapping,
-    const DoFHandlerType &          dof,
-    const hp::QCollection<dim - 1> &quadrature,
+    const Mapping<dim, spacedim> &   mapping,
+    const DoFHandler<dim, spacedim> &dof,
+    const hp::QCollection<dim - 1> & quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -499,11 +501,11 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &          dof,
-    const hp::QCollection<dim - 1> &quadrature,
+    const DoFHandler<dim, spacedim> &dof,
+    const hp::QCollection<dim - 1> & quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -621,12 +623,12 @@ public:
    * respective parameter for compatibility with the function signature in the
    * general case.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<1, spacedim> &mapping,
-    const DoFHandlerType &      dof,
-    const Quadrature<0> &       quadrature,
+    const Mapping<1, spacedim> &   mapping,
+    const DoFHandler<1, spacedim> &dof,
+    const Quadrature<0> &          quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -643,11 +645,11 @@ public:
    * Call the @p estimate function, see above, with
    * <tt>mapping=MappingQGeneric1<1>()</tt>.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &dof,
-    const Quadrature<0> & quadrature,
+    const DoFHandler<1, spacedim> &dof,
+    const Quadrature<0> &          quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -673,12 +675,12 @@ public:
    * construct of vector of references, so we had to use a vector of
    * pointers.)
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<1, spacedim> &mapping,
-    const DoFHandlerType &      dof,
-    const Quadrature<0> &       quadrature,
+    const Mapping<1, spacedim> &   mapping,
+    const DoFHandler<1, spacedim> &dof,
+    const Quadrature<0> &          quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -695,11 +697,11 @@ public:
    * Call the @p estimate function, see above, with
    * <tt>mapping=MappingQGeneric1<1>()</tt>.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &dof,
-    const Quadrature<0> & quadrature,
+    const DoFHandler<1, spacedim> &dof,
+    const Quadrature<0> &          quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -717,12 +719,12 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<1, spacedim> &mapping,
-    const DoFHandlerType &      dof,
-    const hp::QCollection<0> &  quadrature,
+    const Mapping<1, spacedim> &   mapping,
+    const DoFHandler<1, spacedim> &dof,
+    const hp::QCollection<0> &     quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -740,11 +742,11 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &    dof,
-    const hp::QCollection<0> &quadrature,
+    const DoFHandler<1, spacedim> &dof,
+    const hp::QCollection<0> &     quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                       neumann_bc,
@@ -762,12 +764,12 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const Mapping<1, spacedim> &mapping,
-    const DoFHandlerType &      dof,
-    const hp::QCollection<0> &  quadrature,
+    const Mapping<1, spacedim> &   mapping,
+    const DoFHandler<1, spacedim> &dof,
+    const hp::QCollection<0> &     quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
@@ -785,11 +787,11 @@ public:
    * Equivalent to the set of functions above, except that this one takes a
    * quadrature collection for hp finite element dof handlers.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector>
   static void
   estimate(
-    const DoFHandlerType &    dof,
-    const hp::QCollection<0> &quadrature,
+    const DoFHandler<1, spacedim> &dof,
+    const hp::QCollection<0> &     quadrature,
     const std::map<types::boundary_id,
                    const Function<spacedim, typename InputVector::value_type> *>
       &                                     neumann_bc,
index 6d0695630c5f91dec279b0444bdbf34b90d85060..921395dbd29e8535ca664bfffa8a3d7ea6661d6a 100644 (file)
@@ -89,12 +89,9 @@ namespace internal
    * resizing to a smaller size doesn't imply memory allocation, this is
    * fast.
    */
-  template <typename DoFHandlerType, typename number>
+  template <int dim, int spacedim, typename number>
   struct ParallelData
   {
-    static const unsigned int dim      = DoFHandlerType::dimension;
-    static const unsigned int spacedim = DoFHandlerType::space_dimension;
-
     /**
      * The finite element to be used.
      */
@@ -206,9 +203,9 @@ namespace internal
   };
 
 
-  template <typename DoFHandlerType, typename number>
+  template <int dim, int spacedim, typename number>
   template <class FE>
-  ParallelData<DoFHandlerType, number>::ParallelData(
+  ParallelData<dim, spacedim, number>::ParallelData(
     const FE &                                          fe,
     const dealii::hp::QCollection<dim - 1> &            face_quadratures,
     const dealii::hp::MappingCollection<dim, spacedim> &mapping,
@@ -265,9 +262,9 @@ namespace internal
 
 
 
-  template <typename DoFHandlerType, typename number>
+  template <int dim, int spacedim, typename number>
   void
-  ParallelData<DoFHandlerType, number>::resize(
+  ParallelData<dim, spacedim, number>::resize(
     const unsigned int active_fe_index)
   {
     const unsigned int n_q_points   = face_quadratures[active_fe_index].size();
@@ -304,16 +301,16 @@ namespace internal
    * object into a global such map. This is the copier stage of a WorkStream
    * pipeline.
    */
-  template <typename DoFHandlerType>
+  template <int dim, int spacedim>
   void
   copy_local_to_global(
-    const std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-      &local_face_integrals,
-    std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-      &face_integrals)
+    const std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+                   std::vector<double>> &local_face_integrals,
+    std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+             std::vector<double>> &      face_integrals)
   {
     // now copy locally computed elements into the global map
-    for (typename std::map<typename DoFHandlerType::face_iterator,
+    for (typename std::map<typename DoFHandler<dim, spacedim>::face_iterator,
                            std::vector<double>>::const_iterator p =
            local_face_integrals.begin();
          p != local_face_integrals.end();
@@ -339,13 +336,12 @@ namespace internal
    * Actually do the computation based on the evaluated gradients in
    * ParallelData.
    */
-  template <typename DoFHandlerType, typename number>
+  template <int dim, int spacedim, typename number>
   std::vector<double>
-  integrate_over_face(ParallelData<DoFHandlerType, number> &parallel_data,
-                      const typename DoFHandlerType::face_iterator &face,
-                      dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                                               DoFHandlerType::space_dimension>
-                        &fe_face_values_cell)
+  integrate_over_face(
+    ParallelData<dim, spacedim, number> &                    parallel_data,
+    const typename DoFHandler<dim, spacedim>::face_iterator &face,
+    dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_cell)
   {
     const unsigned int n_q_points = parallel_data.psi[0].size(),
                        n_components =
@@ -489,35 +485,26 @@ namespace internal
    * A factor to scale the integral for the face at the boundary. Used for
    * Neumann BC.
    */
-  template <typename DoFHandlerType>
+  template <int dim, int spacedim>
   double
   boundary_face_factor(
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    const unsigned int                                   face_no,
-    const dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                                   DoFHandlerType::space_dimension>
-      &fe_face_values_cell,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const unsigned int                                              face_no,
+    const dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_cell,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
     switch (strategy)
       {
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter_over_24:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter_over_24:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::face_diameter_over_twice_max_degree:
+        case KellyErrorEstimator<dim,
+                                 spacedim>::face_diameter_over_twice_max_degree:
           {
             const double cell_degree =
               fe_face_values_cell.get_fe_collection()[cell->active_fe_index()]
@@ -536,38 +523,27 @@ namespace internal
   /**
    * A factor to scale the integral for the regular face.
    */
-  template <typename DoFHandlerType>
+  template <int dim, int spacedim>
   double
   regular_face_factor(
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    const unsigned int                                   face_no,
-    const dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                                   DoFHandlerType::space_dimension>
-      &fe_face_values_cell,
-    const dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                                   DoFHandlerType::space_dimension>
-      &fe_face_values_neighbor,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const unsigned int                                              face_no,
+    const dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_cell,
+    const dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_neighbor,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
     switch (strategy)
       {
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter_over_24:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter_over_24:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::face_diameter_over_twice_max_degree:
+        case KellyErrorEstimator<dim,
+                                 spacedim>::face_diameter_over_twice_max_degree:
           {
             const double cell_degree =
               fe_face_values_cell.get_fe_collection()[cell->active_fe_index()]
@@ -590,40 +566,30 @@ namespace internal
   /**
    * A factor to scale the integral for the irregular face.
    */
-  template <typename DoFHandlerType>
+  template <int dim, int spacedim>
   double
   irregular_face_factor(
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    const typename DoFHandlerType::active_cell_iterator &neighbor_child,
-    const unsigned int                                   face_no,
-    const unsigned int                                   subface_no,
-    const dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                                   DoFHandlerType::space_dimension>
-      &fe_face_values,
-    dealii::hp::FESubfaceValues<DoFHandlerType::dimension,
-                                DoFHandlerType::space_dimension>
-      &fe_subface_values,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator
+      &                                            neighbor_child,
+    const unsigned int                             face_no,
+    const unsigned int                             subface_no,
+    const dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values,
+    dealii::hp::FESubfaceValues<dim, spacedim> &   fe_subface_values,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
     switch (strategy)
       {
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter_over_24:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter_over_24:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter:
           {
             return 1.0;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::face_diameter_over_twice_max_degree:
+        case KellyErrorEstimator<dim,
+                                 spacedim>::face_diameter_over_twice_max_degree:
           {
             const double cell_degree =
               fe_face_values.get_fe_collection()[cell->active_fe_index()]
@@ -647,32 +613,26 @@ namespace internal
    * A factor used when summing up all the contribution from different faces
    * of each cell.
    */
-  template <typename DoFHandlerType>
+  template <int dim, int spacedim>
   double
-  cell_factor(const typename DoFHandlerType::active_cell_iterator &cell,
-              const unsigned int /*face_no*/,
-              const DoFHandlerType & /*dof_handler*/,
-              const typename KellyErrorEstimator<
-                DoFHandlerType::dimension,
-                DoFHandlerType::space_dimension>::Strategy strategy)
+  cell_factor(
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const unsigned int /*face_no*/,
+    const DoFHandler<dim, spacedim> & /*dof_handler*/,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
     switch (strategy)
       {
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter_over_24:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter_over_24:
           {
             return cell->diameter() / 24;
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::cell_diameter:
+        case KellyErrorEstimator<dim, spacedim>::cell_diameter:
           {
             return cell->diameter();
           }
-        case KellyErrorEstimator<
-          DoFHandlerType::dimension,
-          DoFHandlerType::space_dimension>::face_diameter_over_twice_max_degree:
+        case KellyErrorEstimator<dim,
+                                 spacedim>::face_diameter_over_twice_max_degree:
           {
             return 1.0;
           }
@@ -692,31 +652,23 @@ namespace internal
    * boundary), or the other side's refinement level is the same as that of
    * this side, then handle the integration of these both cases together.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector, int dim, int spacedim>
   void
   integrate_over_regular_face(
     const std::vector<const InputVector *> &solutions,
-    ParallelData<DoFHandlerType, typename InputVector::value_type>
-      &parallel_data,
-    std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-      &                                                  local_face_integrals,
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    const unsigned int                                   face_no,
-    dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                             DoFHandlerType::space_dimension>
-      &fe_face_values_cell,
-    dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                             DoFHandlerType::space_dimension>
-      &fe_face_values_neighbor,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    ParallelData<dim, spacedim, typename InputVector::value_type>
+      &                            parallel_data,
+    std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+             std::vector<double>> &local_face_integrals,
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const unsigned int                                              face_no,
+    dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_cell,
+    dealii::hp::FEFaceValues<dim, spacedim> &fe_face_values_neighbor,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
-    const unsigned int dim = DoFHandlerType::dimension;
-    (void)dim;
-
-    const typename DoFHandlerType::face_iterator face = cell->face(face_no);
-    const unsigned int n_solution_vectors             = solutions.size();
+    const typename DoFHandler<dim, spacedim>::face_iterator face =
+      cell->face(face_no);
+    const unsigned int n_solution_vectors = solutions.size();
 
 
     // initialize data of the restriction
@@ -737,8 +689,8 @@ namespace internal
         Assert(cell->neighbor(face_no).state() == IteratorState::valid,
                ExcInternalError());
 
-        const typename DoFHandlerType::active_cell_iterator neighbor =
-          cell->neighbor(face_no);
+        const typename DoFHandler<dim, spacedim>::active_cell_iterator
+          neighbor = cell->neighbor(face_no);
 
         // find which number the current face has relative to the
         // neighboring cell
@@ -754,11 +706,11 @@ namespace internal
                                        neighbor_neighbor,
                                        cell->active_fe_index());
 
-        factor = regular_face_factor<DoFHandlerType>(cell,
-                                                     face_no,
-                                                     fe_face_values_cell,
-                                                     fe_face_values_neighbor,
-                                                     strategy);
+        factor = regular_face_factor<dim, spacedim>(cell,
+                                                    face_no,
+                                                    fe_face_values_cell,
+                                                    fe_face_values_neighbor,
+                                                    strategy);
 
         // get gradients on neighbor cell
         for (unsigned int n = 0; n < n_solution_vectors; ++n)
@@ -773,10 +725,10 @@ namespace internal
       }
     else
       {
-        factor = boundary_face_factor<DoFHandlerType>(cell,
-                                                      face_no,
-                                                      fe_face_values_cell,
-                                                      strategy);
+        factor = boundary_face_factor<dim, spacedim>(cell,
+                                                     face_no,
+                                                     fe_face_values_cell,
+                                                     strategy);
       }
 
     // now go to the generic function that does all the other things
@@ -794,33 +746,24 @@ namespace internal
    * over face @p face_no of @p cell, where the respective neighbor is
    * refined, so that the integration is a bit more complex.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector, int dim, int spacedim>
   void
   integrate_over_irregular_face(
     const std::vector<const InputVector *> &solutions,
-    ParallelData<DoFHandlerType, typename InputVector::value_type>
-      &parallel_data,
-    std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-      &                                                  local_face_integrals,
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    const unsigned int                                   face_no,
-    dealii::hp::FEFaceValues<DoFHandlerType::dimension,
-                             DoFHandlerType::space_dimension> &fe_face_values,
-    dealii::hp::FESubfaceValues<DoFHandlerType::dimension,
-                                DoFHandlerType::space_dimension>
-      &fe_subface_values,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    ParallelData<dim, spacedim, typename InputVector::value_type>
+      &                            parallel_data,
+    std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+             std::vector<double>> &local_face_integrals,
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    const unsigned int                                              face_no,
+    dealii::hp::FEFaceValues<dim, spacedim> &   fe_face_values,
+    dealii::hp::FESubfaceValues<dim, spacedim> &fe_subface_values,
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
-    const unsigned int dim = DoFHandlerType::dimension;
-    (void)dim;
-
-    const typename DoFHandlerType::cell_iterator neighbor =
-      cell->neighbor(face_no);
+    const auto neighbor = cell->neighbor(face_no);
     (void)neighbor;
-    const unsigned int n_solution_vectors             = solutions.size();
-    const typename DoFHandlerType::face_iterator face = cell->face(face_no);
+    const unsigned int n_solution_vectors = solutions.size();
+    const auto         face               = cell->face(face_no);
 
     Assert(neighbor.state() == IteratorState::valid, ExcInternalError());
     Assert(face->has_children(), ExcInternalError());
@@ -843,8 +786,8 @@ namespace internal
          ++subface_no)
       {
         // get an iterator pointing to the cell behind the present subface
-        const typename DoFHandlerType::active_cell_iterator neighbor_child =
-          cell->neighbor_child_on_subface(face_no, subface_no);
+        const typename DoFHandler<dim, spacedim>::active_cell_iterator
+          neighbor_child = cell->neighbor_child_on_subface(face_no, subface_no);
         Assert(neighbor_child->is_active(), ExcInternalError());
 
         // restrict the finite element on the present cell to the subface
@@ -860,13 +803,13 @@ namespace internal
                               cell->active_fe_index());
 
         const double factor =
-          irregular_face_factor<DoFHandlerType>(cell,
-                                                neighbor_child,
-                                                face_no,
-                                                subface_no,
-                                                fe_face_values,
-                                                fe_subface_values,
-                                                strategy);
+          irregular_face_factor<dim, spacedim>(cell,
+                                               neighbor_child,
+                                               face_no,
+                                               subface_no,
+                                               fe_face_values,
+                                               fe_subface_values,
+                                               strategy);
 
         // store the gradient of the solution in psi
         for (unsigned int n = 0; n < n_solution_vectors; ++n)
@@ -918,20 +861,17 @@ namespace internal
    * This function is only needed in two or three dimensions.  The error
    * estimator in one dimension is implemented separately.
    */
-  template <typename InputVector, typename DoFHandlerType>
+  template <typename InputVector, int dim, int spacedim>
   void
   estimate_one_cell(
-    const typename DoFHandlerType::active_cell_iterator &cell,
-    ParallelData<DoFHandlerType, typename InputVector::value_type>
-      &parallel_data,
-    std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-      &                                     local_face_integrals,
+    const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+    ParallelData<dim, spacedim, typename InputVector::value_type>
+      &                                     parallel_data,
+    std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+             std::vector<double>> &         local_face_integrals,
     const std::vector<const InputVector *> &solutions,
-    const typename KellyErrorEstimator<
-      DoFHandlerType::dimension,
-      DoFHandlerType::space_dimension>::Strategy strategy)
+    const typename KellyErrorEstimator<dim, spacedim>::Strategy strategy)
   {
-    const unsigned int dim                = DoFHandlerType::dimension;
     const unsigned int n_solution_vectors = solutions.size();
 
     const types::subdomain_id subdomain_id = parallel_data.subdomain_id;
@@ -943,7 +883,8 @@ namespace internal
     // loop over all faces of this cell
     for (const unsigned int face_no : GeometryInfo<dim>::face_indices())
       {
-        const typename DoFHandlerType::face_iterator face = cell->face(face_no);
+        const typename DoFHandler<dim, spacedim>::face_iterator face =
+          cell->face(face_no);
 
         // make sure we do work only once: this face may either be regular
         // or irregular. if it is regular and has a neighbor, then we visit
@@ -1060,12 +1001,12 @@ namespace internal
 // the following function is still independent of dimension, but it
 // calls dimension dependent functions
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const Mapping<dim, spacedim> &mapping,
-  const DoFHandlerType &        dof_handler,
-  const Quadrature<dim - 1> &   quadrature,
+  const Mapping<dim, spacedim> &   mapping,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const Quadrature<dim - 1> &      quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -1097,11 +1038,11 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const DoFHandlerType &     dof_handler,
-  const Quadrature<dim - 1> &quadrature,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const Quadrature<dim - 1> &      quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -1130,12 +1071,12 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const Mapping<dim, spacedim> &  mapping,
-  const DoFHandlerType &          dof_handler,
-  const hp::QCollection<dim - 1> &quadrature,
+  const Mapping<dim, spacedim> &   mapping,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const hp::QCollection<dim - 1> & quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -1167,11 +1108,11 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const DoFHandlerType &          dof_handler,
-  const hp::QCollection<dim - 1> &quadrature,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const hp::QCollection<dim - 1> & quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -1201,12 +1142,12 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const Mapping<dim, spacedim> &  mapping,
-  const DoFHandlerType &          dof_handler,
-  const hp::QCollection<dim - 1> &face_quadratures,
+  const Mapping<dim, spacedim> &   mapping,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const hp::QCollection<dim - 1> & face_quadratures,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
@@ -1281,13 +1222,14 @@ KellyErrorEstimator<dim, spacedim>::estimate(
   // the integrated jump of the gradient for each face.  At the end of the
   // function, we again loop over the cells and collect the contributions of
   // the different faces of the cell.
-  std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
+  std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+           std::vector<double>>
     face_integrals;
 
   // all the data needed in the error estimator by each of the threads is
   // gathered in the following structures
   const hp::MappingCollection<dim, spacedim> mapping_collection(mapping);
-  const internal::ParallelData<DoFHandlerType, typename InputVector::value_type>
+  const internal::ParallelData<dim, spacedim, typename InputVector::value_type>
     parallel_data(dof_handler.get_fe_collection(),
                   face_quadratures,
                   mapping_collection,
@@ -1298,28 +1240,29 @@ KellyErrorEstimator<dim, spacedim>::estimate(
                   &neumann_bc,
                   component_mask,
                   coefficients);
-  std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
+  std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+           std::vector<double>>
     sample_local_face_integrals;
 
   // now let's work on all those cells:
   WorkStream::run(
     dof_handler.begin_active(),
-    static_cast<typename DoFHandlerType::active_cell_iterator>(
+    static_cast<typename DoFHandler<dim, spacedim>::active_cell_iterator>(
       dof_handler.end()),
     [&solutions, strategy](
-      const typename DoFHandlerType::active_cell_iterator &cell,
-      internal::ParallelData<DoFHandlerType, typename InputVector::value_type>
-        &parallel_data,
-      std::map<typename DoFHandlerType::face_iterator, std::vector<double>>
-        &local_face_integrals) {
+      const typename DoFHandler<dim, spacedim>::active_cell_iterator &cell,
+      internal::ParallelData<dim, spacedim, typename InputVector::value_type>
+        &                            parallel_data,
+      std::map<typename DoFHandler<dim, spacedim>::face_iterator,
+               std::vector<double>> &local_face_integrals) {
       internal::estimate_one_cell(
         cell, parallel_data, local_face_integrals, solutions, strategy);
     },
     [&face_integrals](
-      const std::map<typename DoFHandlerType::face_iterator,
+      const std::map<typename DoFHandler<dim, spacedim>::face_iterator,
                      std::vector<double>> &local_face_integrals) {
-      internal::copy_local_to_global<DoFHandlerType>(local_face_integrals,
-                                                     face_integrals);
+      internal::copy_local_to_global<dim, spacedim>(local_face_integrals,
+                                                    face_integrals);
     },
     parallel_data,
     sample_local_face_integrals);
@@ -1352,7 +1295,7 @@ KellyErrorEstimator<dim, spacedim>::estimate(
             Assert(face_integrals.find(cell->face(face_no)) !=
                      face_integrals.end(),
                    ExcInternalError());
-            const double factor = internal::cell_factor<DoFHandlerType>(
+            const double factor = internal::cell_factor<dim, spacedim>(
               cell, face_no, dof_handler, strategy);
 
             for (unsigned int n = 0; n < n_solution_vectors; ++n)
@@ -1375,12 +1318,12 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const Mapping<dim, spacedim> &mapping,
-  const DoFHandlerType &        dof_handler,
-  const Quadrature<dim - 1> &   quadrature,
+  const Mapping<dim, spacedim> &   mapping,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const Quadrature<dim - 1> &      quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
@@ -1410,11 +1353,11 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const DoFHandlerType &     dof_handler,
-  const Quadrature<dim - 1> &quadrature,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const Quadrature<dim - 1> &      quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
@@ -1444,11 +1387,11 @@ KellyErrorEstimator<dim, spacedim>::estimate(
 
 
 template <int dim, int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<dim, spacedim>::estimate(
-  const DoFHandlerType &          dof_handler,
-  const hp::QCollection<dim - 1> &quadrature,
+  const DoFHandler<dim, spacedim> &dof_handler,
+  const hp::QCollection<dim - 1> & quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
index 31d0daa456044930de97a8d8a8f1dfc63e766a0f..ebef49ddb95580ddbbb2c6c485f297c0d650b2b3 100644 (file)
@@ -23,17 +23,15 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
   }
 
 for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
-     deal_II_space_dimension : SPACE_DIMENSIONS;
-     DH : DOFHANDLER_TEMPLATES)
+     deal_II_space_dimension : SPACE_DIMENSIONS)
   {
 #if deal_II_dimension != 1 && deal_II_dimension <= deal_II_space_dimension
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -50,9 +48,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -69,10 +66,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -89,9 +85,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -108,10 +103,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -128,9 +122,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -147,10 +140,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -167,9 +159,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
index d812572042bfced41140f0451e49494f83ef984d..911cdd2af2e171e61151b69a88ebe5e0145d3a54 100644 (file)
@@ -57,12 +57,12 @@ DEAL_II_NAMESPACE_OPEN
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const Mapping<1, spacedim> &mapping,
-  const DoFHandlerType &      dof_handler,
-  const Quadrature<0> &       quadrature,
+  const Mapping<1, spacedim> &   mapping,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const Quadrature<0> &          quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -95,11 +95,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const DoFHandlerType &dof_handler,
-  const Quadrature<0> & quadrature,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const Quadrature<0> &          quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -129,11 +129,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const DoFHandlerType &dof_handler,
-  const Quadrature<0> & quadrature,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const Quadrature<0> &          quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
@@ -163,12 +163,12 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const Mapping<1, spacedim> &mapping,
-  const DoFHandlerType &      dof_handler,
-  const hp::QCollection<0> &  quadrature,
+  const Mapping<1, spacedim> &   mapping,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const hp::QCollection<0> &     quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -200,11 +200,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const DoFHandlerType &    dof_handler,
-  const hp::QCollection<0> &quadrature,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const hp::QCollection<0> &     quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                       neumann_bc,
@@ -234,11 +234,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const DoFHandlerType &    dof_handler,
-  const hp::QCollection<0> &quadrature,
+  const DoFHandler<1, spacedim> &dof_handler,
+  const hp::QCollection<0> &     quadrature,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
     &                                     neumann_bc,
@@ -268,11 +268,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
   const Mapping<1, spacedim> & /*mapping*/,
-  const DoFHandlerType & /*dof_handler*/,
+  const DoFHandler<1, spacedim> & /*dof_handler*/,
   const hp::QCollection<0> &,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
@@ -292,11 +292,11 @@ KellyErrorEstimator<1, spacedim>::estimate(
 
 
 template <int spacedim>
-template <typename InputVector, typename DoFHandlerType>
+template <typename InputVector>
 void
 KellyErrorEstimator<1, spacedim>::estimate(
-  const Mapping<1, spacedim> &mapping,
-  const DoFHandlerType &      dof_handler,
+  const Mapping<1, spacedim> &   mapping,
+  const DoFHandler<1, spacedim> &dof_handler,
   const Quadrature<0> &,
   const std::map<types::boundary_id,
                  const Function<spacedim, typename InputVector::value_type> *>
@@ -454,7 +454,7 @@ KellyErrorEstimator<1, spacedim>::estimate(
         for (unsigned int n = 0; n < 2; ++n)
           {
             // find left or right active neighbor
-            typename DoFHandlerType::cell_iterator neighbor = cell->neighbor(n);
+            auto neighbor = cell->neighbor(n);
             if (neighbor.state() == IteratorState::valid)
               while (neighbor->has_children())
                 neighbor = neighbor->child(n == 0 ? 1 : 0);
index 1c04d51fbda1679c8390070442dd9c4a2213f005..fe504279952086ed3432dc20dd99d90de98b4e8e 100644 (file)
 
 
 for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
-     deal_II_space_dimension : SPACE_DIMENSIONS;
-     DH : DOFHANDLER_TEMPLATES)
+     deal_II_space_dimension : SPACE_DIMENSIONS)
   {
 #if deal_II_dimension == 1 && deal_II_dimension <= deal_II_space_dimension
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -41,9 +39,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -59,10 +56,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -78,9 +74,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const Quadrature<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -96,10 +91,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -115,9 +109,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -133,10 +126,9 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
+      VEC>(
       const Mapping<deal_II_dimension, deal_II_space_dimension> &,
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>
@@ -152,9 +144,8 @@ for (VEC : VECTOR_TYPES; deal_II_dimension : DIMENSIONS;
 
     template void
     KellyErrorEstimator<deal_II_dimension, deal_II_space_dimension>::estimate<
-      VEC,
-      DH<deal_II_dimension, deal_II_space_dimension>>(
-      const DH<deal_II_dimension, deal_II_space_dimension> &,
+      VEC>(
+      const DoFHandler<deal_II_dimension, deal_II_space_dimension> &,
       const hp::QCollection<deal_II_dimension - 1> &,
       const std::map<types::boundary_id,
                      const Function<deal_II_space_dimension, VEC::value_type> *>

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.