]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Pass TriaIterator<DoFCellAccessor> by reference ant not by value 270/head
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Mon, 24 Nov 2014 17:09:24 +0000 (18:09 +0100)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Mon, 24 Nov 2014 17:09:24 +0000 (18:09 +0100)
This makes the interface of FEValues::reinit consistent between Triangulation::cell_iterator and DoFHandler::cell_iterator.

include/deal.II/fe/fe_values.h
source/fe/fe_values.cc
source/fe/fe_values.inst.in

index 6018d3b0212bdd7b2f4ebef33f8a72926c77aa81..8493c301b965443f661821b28276f4f15a3e152a 100644 (file)
@@ -2582,7 +2582,7 @@ public:
    * 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);
+  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell);
 
   /**
    * Reinitialize the gradients, Jacobi determinants, etc for the given cell
@@ -2787,7 +2787,7 @@ public:
    * 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,
+  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell,
                const unsigned int face_no);
 
   /**
@@ -2896,7 +2896,7 @@ public:
    * 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,
+  void reinit (const TriaIterator<DoFCellAccessor<DH,level_dof_access> > &cell,
                const unsigned int                    face_no,
                const unsigned int                    subface_no);
 
index 6fac190ccab859827058329fe3f0ca2ea9960187..d0566b0d29bc0f7ca32e749bedad89e3bdcafb24 100644 (file)
@@ -3471,7 +3471,7 @@ void FEValues<dim,spacedim>::reinit (const typename Triangulation<dim,spacedim>:
 template <int dim, int spacedim>
 template <class DH, bool lda>
 void
-FEValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > cell)
+FEValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell)
 {
   // assert that the finite elements
   // passed to the constructor and
@@ -3646,7 +3646,7 @@ FEFaceValues<dim,spacedim>::initialize (const UpdateFlags update_flags)
 template <int dim, int spacedim>
 template <class DH, bool lda>
 void
-FEFaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > cell,
+FEFaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell,
                                     const unsigned int face_no)
 {
   // assert that the finite elements
@@ -3802,7 +3802,7 @@ 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,
+void FESubfaceValues<dim,spacedim>::reinit (const TriaIterator<DoFCellAccessor<DH, lda> > &cell,
                                             const unsigned int         face_no,
                                             const unsigned int         subface_no)
 {
index cdf87b8a99a1e23ad87bcf8e1a562f035846fea4..d57d4328647ee118866c479b3c45dcfedb704811 100644 (file)
@@ -52,11 +52,11 @@ for (dof_handler : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS; deal_II
 #if (deal_II_space_dimension == DIM_A) || (deal_II_space_dimension == DIM_B)
 
     template void FEValues<deal_II_dimension,deal_II_space_dimension>::reinit(
-    TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >);
+    const TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >&);
     template void FEFaceValues<deal_II_dimension,deal_II_space_dimension>::reinit(
-    TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >, unsigned int);
+    const TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >&, unsigned int);
     template void FESubfaceValues<deal_II_dimension,deal_II_space_dimension>::reinit(
-    TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >, unsigned int, unsigned int);
+    const TriaIterator<DoFCellAccessor<dof_handler<deal_II_dimension,deal_II_space_dimension>, lda> >&, unsigned int, unsigned int);
 #endif
 #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.