]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Restrict template arguments for FEValues::reinit(). 1420/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 24 Aug 2015 00:47:58 +0000 (19:47 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 24 Aug 2015 00:47:58 +0000 (19:47 -0500)
include/deal.II/fe/fe_values.h
source/fe/fe_values.cc

index c73190e91401f14f2ce49942b8a4654dfb649e51..e8c1b5e607abe0da820218bf21f4f6fc9b910e69 100644 (file)
@@ -2465,8 +2465,8 @@ 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 <class DH, bool level_dof_access>
-  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell);
+  template <template <int,int> class DH, bool level_dof_access>
+  void reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>,level_dof_access> > &cell);
 
   /**
    * Reinitialize the gradients, Jacobi determinants, etc for the given cell
@@ -2671,8 +2671,8 @@ public:
    * Reinitialize the gradients, Jacobi determinants, etc for the face with
    * number @p face_no of @p cell and the given finite element.
    */
-  template <class DH, bool level_dof_access>
-  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell,
+  template <template <int,int> class DH, bool level_dof_access>
+  void reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>,level_dof_access> > &cell,
                const unsigned int face_no);
 
   /**
@@ -2780,8 +2780,8 @@ 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 <class DH, bool level_dof_access>
-  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell,
+  template <template <int,int> class DH, bool level_dof_access>
+  void reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>,level_dof_access> > &cell,
                const unsigned int                    face_no,
                const unsigned int                    subface_no);
 
index 7003914a6d4b7b39bc02493396c56417c592815e..09b597b7516168e90a2966799365d896cee3f8aa 100644 (file)
@@ -3499,9 +3499,9 @@ void FEValues<dim,spacedim>::reinit (const typename Triangulation<dim,spacedim>:
 
 
 template <int dim, int spacedim>
-template <class DH, bool lda>
+template <template <int,int> class DH, bool lda>
 void
-FEValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell)
+FEValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > &cell)
 {
   // assert that the finite elements
   // passed to the constructor and
@@ -3515,7 +3515,7 @@ FEValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &c
   this->maybe_invalidate_previous_present_cell (cell);
   this->check_cell_similarity(cell);
 
-  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH, lda> > > >
+  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > > >
   (this->present_cell, cell);
 
   // this was the part of the work
@@ -3687,9 +3687,9 @@ FEFaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
 
 
 template <int dim, int spacedim>
-template <class DH, bool lda>
+template <template <int,int> class DH, bool lda>
 void
-FEFaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell,
+FEFaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > &cell,
                                     const unsigned int face_no)
 {
   // assert that the finite elements
@@ -3706,7 +3706,7 @@ FEFaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda>
           ExcIndexRange (face_no, 0, GeometryInfo<dim>::faces_per_cell));
 
   this->maybe_invalidate_previous_present_cell (cell);
-  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH, lda> > > >
+  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > > >
   (this->present_cell, cell);
 
   // this was the part of the work
@@ -3843,8 +3843,8 @@ FESubfaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
 
 
 template <int dim, int spacedim>
-template <class DH, bool lda>
-void FESubfaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell,
+template <template <int,int> class DH, bool lda>
+void FESubfaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > &cell,
                                             const unsigned int         face_no,
                                             const unsigned int         subface_no)
 {
@@ -3879,7 +3879,7 @@ void FESubfaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<D
                       "instead in these cases."));
 
   this->maybe_invalidate_previous_present_cell (cell);
-  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH, lda> > > >
+  reset_pointer_in_place_if_possible<typename FEValuesBase<dim,spacedim>::template CellIterator<TriaIterator<DoFCellAccessor<DH<dim,spacedim>, lda> > > >
   (this->present_cell, cell);
 
   // this was the part of the work

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.