]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove DoFHandlerType from FEValues templates. 10557/head
authorDavid Wells <drwells@email.unc.edu>
Fri, 19 Jun 2020 16:45:59 +0000 (12:45 -0400)
committerDavid Wells <drwells@email.unc.edu>
Fri, 19 Jun 2020 16:45:59 +0000 (12:45 -0400)
include/deal.II/fe/fe_values.h
source/fe/fe_values.cc
source/fe/fe_values.inst.in

index 351f1b1c4dd18adcefa77dd23c864d97124f4205..3c23680332ea9cceadb6bdf5db06f125894e47dc 100644 (file)
@@ -38,8 +38,6 @@
 #include <deal.II/grid/tria.h>
 #include <deal.II/grid/tria_iterator.h>
 
-#include <deal.II/hp/dof_handler.h>
-
 #include <algorithm>
 #include <memory>
 #include <type_traits>
@@ -3634,10 +3632,10 @@ public:
    * associated with this object. It is assumed that the finite element used
    * by the given cell is also the one used by this FEValues object.
    */
-  template <template <int, int> class DoFHandlerType, bool level_dof_access>
+  template <bool level_dof_access>
   void
-  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
-                                            level_dof_access>> &cell);
+  reinit(const TriaIterator<
+         DoFCellAccessor<DoFHandler<dim, spacedim>, level_dof_access>> &cell);
 
   /**
    * Reinitialize the gradients, Jacobi determinants, etc for the given cell
@@ -3854,11 +3852,11 @@ public:
    * Reinitialize the gradients, Jacobi determinants, etc for the face with
    * number @p face_no of @p cell and the given finite element.
    */
-  template <template <int, int> class DoFHandlerType, bool level_dof_access>
+  template <bool level_dof_access>
   void
-  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
-                                            level_dof_access>> &cell,
-         const unsigned int                                     face_no);
+  reinit(const TriaIterator<
+           DoFCellAccessor<DoFHandler<dim, spacedim>, level_dof_access>> &cell,
+         const unsigned int face_no);
 
   /**
    * Reinitialize the gradients, Jacobi determinants, etc for face @p face
@@ -3866,11 +3864,11 @@ public:
    *
    * @note @p face must be one of @p cell's face iterators.
    */
-  template <template <int, int> class DoFHandlerType, bool level_dof_access>
+  template <bool level_dof_access>
   void
-  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
-                                            level_dof_access>> &     cell,
-         const typename Triangulation<dim, spacedim>::face_iterator &face);
+  reinit(const TriaIterator<
+           DoFCellAccessor<DoFHandler<dim, spacedim>, level_dof_access>> &cell,
+         const typename Triangulation<dim, spacedim>::face_iterator &     face);
 
   /**
    * Reinitialize the gradients, Jacobi determinants, etc for the given face
@@ -4003,22 +4001,22 @@ public:
    * associated with this object. It is assumed that the finite element used
    * by the given cell is also the one used by this FESubfaceValues object.
    */
-  template <template <int, int> class DoFHandlerType, bool level_dof_access>
+  template <bool level_dof_access>
   void
-  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
-                                            level_dof_access>> &cell,
-         const unsigned int                                     face_no,
-         const unsigned int                                     subface_no);
+  reinit(const TriaIterator<
+           DoFCellAccessor<DoFHandler<dim, spacedim>, level_dof_access>> &cell,
+         const unsigned int face_no,
+         const unsigned int subface_no);
 
   /**
    * Alternative reinitialization function that takes, as arguments, iterators
    * to the face and subface instead of their numbers.
    */
-  template <template <int, int> class DoFHandlerType, bool level_dof_access>
+  template <bool level_dof_access>
   void
-  reinit(const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>,
-                                            level_dof_access>> &     cell,
-         const typename Triangulation<dim, spacedim>::face_iterator &face,
+  reinit(const TriaIterator<
+           DoFCellAccessor<DoFHandler<dim, spacedim>, level_dof_access>> &cell,
+         const typename Triangulation<dim, spacedim>::face_iterator &     face,
          const typename Triangulation<dim, spacedim>::face_iterator &subface);
 
   /**
index 583d3238b8706569d183701feadd7700ac60ac2e..e706b49197d175c8a2ff5d4fccad5c881bfca86f 100644 (file)
@@ -4524,10 +4524,10 @@ FEValues<dim, spacedim>::reinit(
 
 
 template <int dim, int spacedim>
-template <template <int, int> class DoFHandlerType, bool lda>
+template <bool lda>
 void
 FEValues<dim, spacedim>::reinit(
-  const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell)
+  const TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>> &cell)
 {
   // assert that the finite elements passed to the constructor and
   // used by the DoFHandler used by this cell, are the same
@@ -4540,7 +4540,7 @@ FEValues<dim, spacedim>::reinit(
 
   reset_pointer_in_place_if_possible<
     typename FEValuesBase<dim, spacedim>::template CellIterator<
-      TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>>>>(
+      TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>>>>(
     this->present_cell, cell);
 
   // this was the part of the work that is dependent on the actual
@@ -4731,11 +4731,11 @@ FEFaceValues<dim, spacedim>::initialize(const UpdateFlags update_flags)
 
 
 template <int dim, int spacedim>
-template <template <int, int> class DoFHandlerType, bool lda>
+template <bool lda>
 void
 FEFaceValues<dim, spacedim>::reinit(
-  const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
-  const unsigned int face_no)
+  const TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>> &cell,
+  const unsigned int                                                   face_no)
 {
   // assert that the finite elements passed to the constructor and
   // used by the DoFHandler used by this cell, are the same
@@ -4749,7 +4749,7 @@ FEFaceValues<dim, spacedim>::reinit(
   this->maybe_invalidate_previous_present_cell(cell);
   reset_pointer_in_place_if_possible<
     typename FEValuesBase<dim, spacedim>::template CellIterator<
-      TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>>>>(
+      TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>>>>(
     this->present_cell, cell);
 
   // this was the part of the work that is dependent on the actual
@@ -4761,11 +4761,11 @@ FEFaceValues<dim, spacedim>::reinit(
 
 
 template <int dim, int spacedim>
-template <template <int, int> class DoFHandlerType, bool lda>
+template <bool lda>
 void
 FEFaceValues<dim, spacedim>::reinit(
-  const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
-  const typename Triangulation<dim, spacedim>::face_iterator &             face)
+  const TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>> &cell,
+  const typename Triangulation<dim, spacedim>::face_iterator &         face)
 {
   const auto face_n = cell->face_iterator_to_index(face);
   reinit(cell, face_n);
@@ -4932,15 +4932,15 @@ FESubfaceValues<dim, spacedim>::initialize(const UpdateFlags update_flags)
 
 
 template <int dim, int spacedim>
-template <template <int, int> class DoFHandlerType, bool lda>
+template <bool lda>
 void
 FESubfaceValues<dim, spacedim>::reinit(
-  const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
-  const unsigned int face_no,
+  const TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>> &cell,
+  const unsigned int                                                   face_no,
   const unsigned int subface_no)
 {
   // assert that the finite elements passed to the constructor and
-  // used by the hp::DoFHandler used by this cell, are the same
+  // used by the DoFHandler used by this cell, are the same
   Assert(static_cast<const FiniteElementData<dim> &>(*this->fe) ==
            static_cast<const FiniteElementData<dim> &>(
              cell->get_dof_handler().get_fe(cell->active_fe_index())),
@@ -4968,7 +4968,7 @@ FESubfaceValues<dim, spacedim>::reinit(
   this->maybe_invalidate_previous_present_cell(cell);
   reset_pointer_in_place_if_possible<
     typename FEValuesBase<dim, spacedim>::template CellIterator<
-      TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>>>>(
+      TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>>>>(
     this->present_cell, cell);
 
   // this was the part of the work that is dependent on the actual
@@ -4980,12 +4980,12 @@ FESubfaceValues<dim, spacedim>::reinit(
 
 
 template <int dim, int spacedim>
-template <template <int, int> class DoFHandlerType, bool lda>
+template <bool lda>
 void
 FESubfaceValues<dim, spacedim>::reinit(
-  const TriaIterator<DoFCellAccessor<DoFHandlerType<dim, spacedim>, lda>> &cell,
-  const typename Triangulation<dim, spacedim>::face_iterator &             face,
-  const typename Triangulation<dim, spacedim>::face_iterator &subface)
+  const TriaIterator<DoFCellAccessor<DoFHandler<dim, spacedim>, lda>> &cell,
+  const typename Triangulation<dim, spacedim>::face_iterator &         face,
+  const typename Triangulation<dim, spacedim>::face_iterator &         subface)
 {
   reinit(cell,
          cell->face_iterator_to_index(face),
index efc80295f5c3efd9478e434ea98b8e179aab610a..317dedf94b8eb2a9cd484e58faf0c4b39faf5a69 100644 (file)
@@ -52,8 +52,7 @@ for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
   }
 
 
-for (DoFHandler : DOFHANDLER_TEMPLATE; deal_II_dimension : DIMENSIONS;
-     deal_II_space_dimension : SPACE_DIMENSIONS;
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS;
      lda : BOOL)
   {
 #if deal_II_dimension <= deal_II_space_dimension

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.