]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make the internal data object of the mapping classes const.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 22 Jul 2015 21:49:36 +0000 (16:49 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 26 Jul 2015 14:16:07 +0000 (09:16 -0500)
Logically, they are input arguments to the fill_fe_values() functions, they are only
initialized in get_data() etc. Thus, make them const in the fill_fe_values functions.
There are places where mapping classes want to move data across from fill_fe_values
in terms of temporary arrays to the transform() function that may be called right
afterwards from FiniteElement::fill_fe_values(). In those cases, make these members
mutable.

include/deal.II/fe/mapping.h
include/deal.II/fe/mapping_cartesian.h
include/deal.II/fe/mapping_fe_field.h
include/deal.II/fe/mapping_q.h
include/deal.II/fe/mapping_q1.h
source/fe/mapping_cartesian.cc
source/fe/mapping_fe_field.cc
source/fe/mapping_q.cc
source/fe/mapping_q1.cc

index a70abd7b6d51b38ecc402a7c205393890bcda607..cd5d6cf75a5641538fe4704626949d469f9ee18e 100644 (file)
@@ -324,7 +324,7 @@ public:
      * The determinant of the Jacobian in each quadrature point. Filled if
      * #update_volume_elements.
      */
-    std::vector<double> volume_elements;
+    mutable std::vector<double> volume_elements;
 
     /**
      * The positions of the mapped (generalized) support points.
@@ -644,15 +644,14 @@ private:
   CellSimilarity::Similarity
   fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                   const Quadrature<dim>                         &quadrature,
-                  InternalDataBase                              &internal,
+                  const InternalDataBase                              &internal,
                   std::vector<Point<spacedim> >                 &quadrature_points,
                   std::vector<double>                           &JxW_values,
                   std::vector<DerivativeForm<1,dim,spacedim>  > &jacobians,
                   std::vector<DerivativeForm<2,dim,spacedim>  > &jacobian_grads,
                   std::vector<DerivativeForm<1,spacedim,dim>  > &inverse_jacobians,
                   std::vector<Point<spacedim> >                 &cell_normal_vectors,
-                  const CellSimilarity::Similarity                    cell_similarity
-                 ) const=0;
+                  const CellSimilarity::Similarity                    cell_similarity) const=0;
 
 
 
@@ -670,7 +669,7 @@ private:
   fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                        const unsigned int                            face_no,
                        const Quadrature<dim-1>                      &quadrature,
-                       InternalDataBase                             &internal,
+                       const InternalDataBase                             &internal,
                        std::vector<Point<spacedim> >                &quadrature_points,
                        std::vector<double>                          &JxW_values,
                        std::vector<Tensor<1,spacedim> >             &boundary_form,
@@ -686,7 +685,7 @@ private:
                           const unsigned int                face_no,
                           const unsigned int                sub_no,
                           const Quadrature<dim-1>          &quadrature,
-                          InternalDataBase                 &internal,
+                          const InternalDataBase                 &internal,
                           std::vector<Point<spacedim> >    &quadrature_points,
                           std::vector<double>              &JxW_values,
                           std::vector<Tensor<1,spacedim> > &boundary_form,
index cb3b7a8505c74708fb7702b7b5248ae8f95e7b00..2d151d82b2786beaa979e25163a9265cefc9b5bc 100644 (file)
@@ -73,7 +73,7 @@ public:
   CellSimilarity::Similarity
   fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                   const Quadrature<dim>                             &quadrature,
-                  typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
+                  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
                   std::vector<Point<spacedim> >                     &quadrature_points,
                   std::vector<double>                               &JxW_values,
                   std::vector<DerivativeForm<1,dim,spacedim> >      &jacobians,
@@ -87,7 +87,7 @@ public:
   fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                        const unsigned int               face_no,
                        const Quadrature<dim-1>&         quadrature,
-                       typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
+                       const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
                        std::vector<Point<dim> >        &quadrature_points,
                        std::vector<double>             &JxW_values,
                        std::vector<Tensor<1,dim> >     &boundary_form,
@@ -99,7 +99,7 @@ public:
                           const unsigned int              face_no,
                           const unsigned int              sub_no,
                           const Quadrature<dim-1>&        quadrature,
-                          typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
+                          const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
                           std::vector<Point<dim> >        &quadrature_points,
                           std::vector<double>             &JxW_values,
                           std::vector<Tensor<1,dim> >     &boundary_form,
@@ -160,6 +160,13 @@ public:
 protected:
   /**
    * Storage for internal data of the scaling.
+   *
+   * This includes data that is computed once when the object is created
+   * (in get_data()) as well as data the class wants to store from between
+   * the call to fill_fe_values(), fill_fe_face_values(), or
+   * fill_fe_subface_values() until possible later calls from the finite
+   * element to functions such as transform(). The latter class of
+   * member variables are marked as 'mutable'.
    */
   class InternalData : public Mapping<dim, spacedim>::InternalDataBase
   {
@@ -178,12 +185,12 @@ protected:
      * Length of the cell in different coordinate directions,
      * <i>h<sub>x</sub></i>, <i>h<sub>y</sub></i>, <i>h<sub>z</sub></i>.
      */
-    Tensor<1,dim> length;
+    mutable Tensor<1,dim> length;
 
     /**
      * The volume element
      */
-    double volume_element;
+    mutable double volume_element;
 
     /**
      * Vector of all quadrature points. Especially, all points on all faces.
@@ -198,7 +205,7 @@ protected:
                      const unsigned int face_no,
                      const unsigned int sub_no,
                      const CellSimilarity::Similarity cell_similarity,
-                     InternalData &data,
+                     const InternalData &data,
                      std::vector<Point<dim> > &quadrature_points,
                      std::vector<Point<dim> > &normal_vectors) const;
 
index 8324e956f4b474e12dd3b7e9470d3384376cf536..66ede966605c28e0bf213bcd4149abe4cc6d573c 100644 (file)
@@ -208,8 +208,8 @@ public:
      * Shape function at quadrature point. Shape functions are in tensor
      * product order, so vertices must be reordered to obtain transformation.
      */
-    double shape (const unsigned int qpoint,
-                  const unsigned int shape_nr) const;
+    const double &shape (const unsigned int qpoint,
+                         const unsigned int shape_nr) const;
 
     /**
      * Shape function at quadrature point. See above.
@@ -220,8 +220,8 @@ public:
     /**
      * Gradient of shape function in quadrature point. See above.
      */
-    Tensor<1,dim> derivative (const unsigned int qpoint,
-                              const unsigned int shape_nr) const;
+    const Tensor<1,dim> &derivative (const unsigned int qpoint,
+                                     const unsigned int shape_nr) const;
 
     /**
      * Gradient of shape function in quadrature point. See above.
@@ -277,7 +277,7 @@ public:
      *
      * Computed on each cell.
      */
-    std::vector<DerivativeForm<1,dim, spacedim > >  covariant;
+    mutable std::vector<DerivativeForm<1,dim, spacedim > >  covariant;
 
     /**
      * Tensors of contravariant transformation at each of the quadrature
@@ -286,7 +286,7 @@ public:
      *
      * Computed on each cell.
      */
-    std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
+    mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
 
     /**
      * Unit tangential vectors. Used for the computation of boundary forms and
@@ -306,7 +306,7 @@ public:
     /**
      * Auxiliary vectors for internal use.
      */
-    std::vector<std::vector<Tensor<1,spacedim> > > aux;
+    mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
 
     /**
      * Number of shape functions. If this is a Q1 mapping, then it is simply
@@ -390,7 +390,7 @@ public:
                      const unsigned int      npts,
                      const typename QProjector<dim>::DataSetDescriptor data_set,
                      const CellSimilarity::Similarity cell_similarity,
-                     InternalData           &data,
+                     const InternalData           &data,
                      std::vector<Point<spacedim> > &quadrature_points) const;
 
 
@@ -403,7 +403,7 @@ public:
                           const unsigned int      npts,
                           const typename QProjector<dim>::DataSetDescriptor data_set,
                           const std::vector<double>   &weights,
-                          InternalData           &mapping_data,
+                          const InternalData           &mapping_data,
                           std::vector<Point<spacedim> >    &quadrature_points,
                           std::vector<double>         &JxW_values,
                           std::vector<Tensor<1,spacedim> > &boundary_form,
@@ -428,7 +428,7 @@ protected:
   CellSimilarity::Similarity
   fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                   const Quadrature<dim>                                     &quadrature,
-                  typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
+                  const typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
                   typename std::vector<Point<spacedim> >                    &quadrature_points,
                   std::vector<double>                                       &JxW_values,
                   std::vector<DerivativeForm<1,dim,spacedim> >       &jacobians,
@@ -444,7 +444,7 @@ protected:
   fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                        const unsigned int face_no,
                        const Quadrature<dim-1>& quadrature,
-                       typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                       const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                        typename std::vector<Point<spacedim> >        &quadrature_points,
                        std::vector<double>             &JxW_values,
                        std::vector<Tensor<1,spacedim> >             &exterior_forms,
@@ -460,7 +460,7 @@ protected:
                           const unsigned int face_no,
                           const unsigned int sub_no,
                           const Quadrature<dim-1>& quadrature,
-                          typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                          const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                           typename std::vector<Point<spacedim> >        &quadrature_points,
                           std::vector<double>             &JxW_values,
                           std::vector<Tensor<1,spacedim> > &exterior_forms,
@@ -627,7 +627,7 @@ private:
 
 template<int dim, int spacedim, class DH, class VECTOR>
 inline
-double
+const double &
 MappingFEField<dim,spacedim,DH,VECTOR>::InternalData::shape (const unsigned int qpoint,
     const unsigned int shape_nr) const
 {
@@ -654,7 +654,7 @@ MappingFEField<dim,spacedim,DH,VECTOR>::InternalData::shape (const unsigned int
 
 template<int dim, int spacedim, class DH, class VECTOR>
 inline
-Tensor<1,dim>
+const Tensor<1,dim> &
 MappingFEField<dim,spacedim,DH,VECTOR>::InternalData::derivative (const unsigned int qpoint,
     const unsigned int shape_nr) const
 {
index 8bbe509fe52e61e7d50b42e75110bfb90d4ea2eb..3e44964e2360719a5377b4955395aad2bb959653 100644 (file)
@@ -173,7 +173,7 @@ public:
      * If this flag is @p true we are on an interior cell and the @p
      * mapping_q1_data is used.
      */
-    bool use_mapping_q1_on_current_cell;
+    mutable bool use_mapping_q1_on_current_cell;
 
     /**
      * On interior cells @p MappingQ1 is used.
@@ -189,7 +189,7 @@ protected:
   CellSimilarity::Similarity
   fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                   const Quadrature<dim>                            &quadrature,
-                  typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                   typename std::vector<Point<spacedim> >           &quadrature_points,
                   std::vector<double>                              &JxW_values,
                   std::vector<DerivativeForm<1,dim,spacedim> >     &jacobians,
@@ -205,7 +205,7 @@ protected:
   fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                        const unsigned int                           face_no,
                        const Quadrature<dim-1>&                     quadrature,
-                       typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                       const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                        typename std::vector<Point<spacedim> >       &quadrature_points,
                        std::vector<double>                          &JxW_values,
                        typename std::vector<Tensor<1,spacedim> >    &exterior_form,
@@ -221,7 +221,7 @@ protected:
                           const unsigned int                           face_no,
                           const unsigned int                           sub_no,
                           const Quadrature<dim-1>&                     quadrature,
-                          typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                          const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                           typename std::vector<Point<spacedim> >       &quadrature_points,
                           std::vector<double>                          &JxW_values,
                           typename std::vector<Tensor<1,spacedim> >    &exterior_form,
index df8276af0a0382f8487b02e4c94e638adc61c0e9..2b6e7a429d84f50fd17f9e595ce08d0fd1eb529a 100644 (file)
@@ -167,8 +167,8 @@ public:
      * Shape function at quadrature point. Shape functions are in tensor
      * product order, so vertices must be reordered to obtain transformation.
      */
-    double shape (const unsigned int qpoint,
-                  const unsigned int shape_nr) const;
+    const double &shape (const unsigned int qpoint,
+                         const unsigned int shape_nr) const;
 
     /**
      * Shape function at quadrature point. See above.
@@ -179,8 +179,8 @@ public:
     /**
      * Gradient of shape function in quadrature point. See above.
      */
-    Tensor<1,dim> derivative (const unsigned int qpoint,
-                              const unsigned int shape_nr) const;
+    const Tensor<1,dim> &derivative (const unsigned int qpoint,
+                                     const unsigned int shape_nr) const;
 
     /**
      * Gradient of shape function in quadrature point. See above.
@@ -191,8 +191,8 @@ public:
     /**
      * Second derivative of shape function in quadrature point. See above.
      */
-    Tensor<2,dim> second_derivative (const unsigned int qpoint,
-                                     const unsigned int shape_nr) const;
+    const Tensor<2,dim> &second_derivative (const unsigned int qpoint,
+                                            const unsigned int shape_nr) const;
 
     /**
      * Second derivative of shape function in quadrature point. See above.
@@ -236,7 +236,7 @@ public:
      *
      * Computed on each cell.
      */
-    std::vector<DerivativeForm<1,dim, spacedim > >  covariant;
+    mutable std::vector<DerivativeForm<1,dim, spacedim > >  covariant;
 
     /**
      * Tensors of contravariant transformation at each of the quadrature
@@ -245,7 +245,7 @@ public:
      *
      * Computed on each cell.
      */
-    std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
+    mutable std::vector< DerivativeForm<1,dim,spacedim> > contravariant;
 
     /**
      * Unit tangential vectors. Used for the computation of boundary forms and
@@ -265,18 +265,18 @@ public:
     /**
      * Auxiliary vectors for internal use.
      */
-    std::vector<std::vector<Tensor<1,spacedim> > > aux;
+    mutable std::vector<std::vector<Tensor<1,spacedim> > > aux;
 
     /**
      * Stores the support points of the mapping shape functions on the @p
      * cell_of_current_support_points.
      */
-    std::vector<Point<spacedim> > mapping_support_points;
+    mutable std::vector<Point<spacedim> > mapping_support_points;
 
     /**
      * Stores the cell of which the @p mapping_support_points are stored.
      */
-    typename Triangulation<dim,spacedim>::cell_iterator cell_of_current_support_points;
+    mutable typename Triangulation<dim,spacedim>::cell_iterator cell_of_current_support_points;
 
     /**
      * Default value of this flag is @p true. If <tt>*this</tt> is an object
@@ -308,7 +308,7 @@ public:
   CellSimilarity::Similarity
   fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                   const Quadrature<dim>                        &quadrature,
-                  typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                   typename std::vector<Point<spacedim> >       &quadrature_points,
                   std::vector<double>                          &JxW_values,
                   std::vector<DerivativeForm<1,dim,spacedim> > &jacobians,
@@ -324,7 +324,7 @@ public:
   fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                        const unsigned int                            face_no,
                        const Quadrature<dim-1>                      &quadrature,
-                       typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                       const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                        typename std::vector<Point<spacedim> >       &quadrature_points,
                        std::vector<double>                          &JxW_values,
                        typename std::vector<Tensor<1,spacedim> >    &boundary_form,
@@ -340,7 +340,7 @@ public:
                           const unsigned int                           face_no,
                           const unsigned int                           sub_no,
                           const Quadrature<dim-1>&                     quadrature,
-                          typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                          const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                           typename std::vector<Point<spacedim> >       &quadrature_points,
                           std::vector<double>                          &JxW_values,
                           typename std::vector<Tensor<1,spacedim> >    &boundary_form,
@@ -386,7 +386,7 @@ public:
                      const unsigned int      npts,
                      const DataSetDescriptor data_set,
                      const CellSimilarity::Similarity cell_similarity,
-                     InternalData           &data,
+                     const InternalData           &data,
                      std::vector<Point<spacedim> > &quadrature_points) const;
 
   /**
@@ -398,7 +398,7 @@ public:
                           const unsigned int                npts,
                           const DataSetDescriptor           data_set,
                           const std::vector<double>        &weights,
-                          InternalData                     &mapping_data,
+                          const InternalData                     &mapping_data,
                           std::vector<Point<spacedim> >    &quadrature_points,
                           std::vector<double>              &JxW_values,
                           std::vector<Tensor<1,spacedim> > &boundary_form,
@@ -633,7 +633,7 @@ struct StaticMappingQ1
 
 template<int dim, int spacedim>
 inline
-double
+const double &
 MappingQ1<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
                                               const unsigned int shape_nr) const
 {
@@ -660,7 +660,7 @@ MappingQ1<dim,spacedim>::InternalData::shape (const unsigned int qpoint,
 
 template<int dim, int spacedim>
 inline
-Tensor<1,dim>
+const Tensor<1,dim> &
 MappingQ1<dim,spacedim>::InternalData::derivative (const unsigned int qpoint,
                                                    const unsigned int shape_nr) const
 {
@@ -687,7 +687,7 @@ MappingQ1<dim,spacedim>::InternalData::derivative (const unsigned int qpoint,
 
 template <int dim, int spacedim>
 inline
-Tensor<2,dim>
+const Tensor<2,dim> &
 MappingQ1<dim,spacedim>::InternalData::second_derivative (const unsigned int qpoint,
                                                           const unsigned int shape_nr) const
 {
index 78600e6de3eb5cdcd0929ec3061b5121e2f932f1..a2ae39b48301873d52882be8277a897c6e55c2e8 100644 (file)
@@ -136,7 +136,7 @@ MappingCartesian<dim, spacedim>::compute_fill (const typename Triangulation<dim,
                                                const unsigned int        face_no,
                                                const unsigned int        sub_no,
                                                const CellSimilarity::Similarity cell_similarity,
-                                               InternalData             &data,
+                                               const InternalData             &data,
                                                std::vector<Point<dim> > &quadrature_points,
                                                std::vector<Point<dim> > &normal_vectors) const
 {
@@ -319,7 +319,7 @@ CellSimilarity::Similarity
 MappingCartesian<dim, spacedim>::
 fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                 const Quadrature<dim> &q,
-                typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                 std::vector<Point<spacedim> > &quadrature_points,
                 std::vector<double> &JxW_values,
                 std::vector< DerivativeForm<1,dim,spacedim> > &jacobians,
@@ -332,8 +332,8 @@ fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   // data for this class. fails with
   // an exception if that is not
   // possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0, ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0, ExcInternalError());
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   std::vector<Point<dim> > dummy;
 
@@ -398,7 +398,7 @@ MappingCartesian<dim, spacedim>::fill_fe_face_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const unsigned int             face_no,
   const Quadrature<dim-1>       &q,
-  typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
+  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
   std::vector<Point<dim> >      &quadrature_points,
   std::vector<double>           &JxW_values,
   std::vector<Tensor<1,dim> >   &boundary_forms,
@@ -410,9 +410,9 @@ MappingCartesian<dim, spacedim>::fill_fe_face_values (
   // data for this class. fails with
   // an exception if that is not
   // possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   compute_fill (cell, face_no, invalid_face_number,
                 CellSimilarity::none,
@@ -471,7 +471,7 @@ MappingCartesian<dim, spacedim>::fill_fe_subface_values (
   const unsigned int             face_no,
   const unsigned int             sub_no,
   const Quadrature<dim-1>       &q,
-  typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
+  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
   std::vector<Point<dim> >      &quadrature_points,
   std::vector<double>           &JxW_values,
   std::vector<Tensor<1,dim> >   &boundary_forms,
@@ -483,8 +483,8 @@ MappingCartesian<dim, spacedim>::fill_fe_subface_values (
   // data for this class. fails with
   // an exception if that is not
   // possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0, ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0, ExcInternalError());
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   compute_fill (cell, face_no, sub_no, CellSimilarity::none,
                 data,
index 703c5eb706ea5d7ed3c40114951bab37f61c5ce9..0acb705a55d4ce3d693074e38f82badc669e29d0 100644 (file)
@@ -376,7 +376,7 @@ CellSimilarity::Similarity
 MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const Quadrature<dim>                                     &q,
-  typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
+  const typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
   std::vector<Point<spacedim> >                             &quadrature_points,
   std::vector<double>                                       &JxW_values,
   std::vector<DerivativeForm<1,dim,spacedim>   >          &jacobians,
@@ -387,8 +387,8 @@ MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_values (
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0, ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0, ExcInternalError());
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   const unsigned int n_q_points=q.size();
   const  CellSimilarity::Similarity updated_cell_similarity
@@ -554,7 +554,7 @@ MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_face_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const unsigned int       face_no,
   const Quadrature<dim-1> &q,
-  typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
   std::vector<Point<spacedim> >     &quadrature_points,
   std::vector<double>          &JxW_values,
   std::vector<Tensor<1,spacedim> >             &exterior_forms,
@@ -564,9 +564,9 @@ MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_face_values (
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   const unsigned int n_q_points=q.size();
   this->compute_fill_face (cell, face_no, numbers::invalid_unsigned_int,
@@ -591,7 +591,7 @@ MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_subface_values (const typename T
     const unsigned int       face_no,
     const unsigned int       sub_no,
     const Quadrature<dim-1> &q,
-    typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+    const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
     std::vector<Point<spacedim> >     &quadrature_points,
     std::vector<double>          &JxW_values,
     std::vector<Tensor<1,spacedim> > &exterior_forms,
@@ -601,9 +601,9 @@ MappingFEField<dim,spacedim,VECTOR,DH>::fill_fe_subface_values (const typename T
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   const unsigned int n_q_points=q.size();
   this->compute_fill_face (cell, face_no, sub_no,
@@ -967,7 +967,7 @@ MappingFEField<dim,spacedim,VECTOR,DH>::compute_fill (
   const unsigned int  n_q_points,
   const typename QProjector<dim>::DataSetDescriptor  data_set,
   const CellSimilarity::Similarity cell_similarity,
-  InternalData  &data,
+  const InternalData  &data,
   std::vector<Point<spacedim> > &quadrature_points) const
 {
   const UpdateFlags update_flags(data.current_update_flags());
@@ -1064,7 +1064,7 @@ MappingFEField<dim,spacedim,VECTOR,DH>::compute_fill_face (
   const unsigned int      n_q_points,//npts
   const typename QProjector<dim>::DataSetDescriptor data_set,
   const std::vector<double>   &weights,
-  InternalData           &data,
+  const InternalData           &data,
   std::vector<Point<spacedim> >    &quadrature_points,
   std::vector<double>         &JxW_values,
   std::vector<Tensor<1,spacedim> > &boundary_forms,
index 5c60da03d7faa80e50804207fec727a6122889ff..048d6117159242b06a10cef1a19ab447f7a2bd7f 100644 (file)
@@ -250,7 +250,7 @@ CellSimilarity::Similarity
 MappingQ<dim,spacedim>::fill_fe_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const Quadrature<dim>                                     &q,
-  typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
+  const typename Mapping<dim,spacedim>::InternalDataBase          &mapping_data,
   std::vector<Point<spacedim> >                             &quadrature_points,
   std::vector<double>                                       &JxW_values,
   std::vector<DerivativeForm<1,dim,spacedim>   >    &jacobians,
@@ -261,8 +261,8 @@ MappingQ<dim,spacedim>::fill_fe_values (
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0, ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0, ExcInternalError());
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   // check whether this cell needs the full mapping or can be treated by a
   // reduced Q1 mapping, e.g. if the cell is in the interior of the domain
@@ -271,7 +271,7 @@ MappingQ<dim,spacedim>::fill_fe_values (
 
   // depending on this result, use this or the other data object for the
   // mapping.
-  typename MappingQ1<dim,spacedim>::InternalData *
+  const typename MappingQ1<dim,spacedim>::InternalData *
   p_data = (data.use_mapping_q1_on_current_cell
             ?
             &data.mapping_q1_data
@@ -309,7 +309,7 @@ MappingQ<dim,spacedim>::fill_fe_face_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const unsigned int                face_no,
   const Quadrature<dim-1>          &q,
-  typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
   std::vector<Point<spacedim> >    &quadrature_points,
   std::vector<double>              &JxW_values,
   std::vector<Tensor<1,spacedim> > &exterior_forms,
@@ -319,9 +319,9 @@ MappingQ<dim,spacedim>::fill_fe_face_values (
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   // check whether this cell needs the full mapping or can be treated by a
   // reduced Q1 mapping, e.g. if the cell is entirely in the interior of the
@@ -334,11 +334,12 @@ MappingQ<dim,spacedim>::fill_fe_face_values (
 
   // depending on this result, use this or the other data object for the
   // mapping
-  typename MappingQ1<dim,spacedim>::InternalData *p_data=0;
-  if (data.use_mapping_q1_on_current_cell)
-    p_data=&data.mapping_q1_data;
-  else
-    p_data=&data;
+  const typename MappingQ1<dim,spacedim>::InternalData *p_data
+    = (data.use_mapping_q1_on_current_cell
+       ?
+       &data.mapping_q1_data
+       :
+       &data);
 
   const unsigned int n_q_points=q.size();
   this->compute_fill_face (cell, face_no, numbers::invalid_unsigned_int,
@@ -363,7 +364,7 @@ MappingQ<dim,spacedim>::fill_fe_subface_values (const typename Triangulation<dim
                                                 const unsigned int       face_no,
                                                 const unsigned int       sub_no,
                                                 const Quadrature<dim-1> &q,
-                                                typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                                                const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                                                 std::vector<Point<spacedim> >     &quadrature_points,
                                                 std::vector<double>          &JxW_values,
                                                 std::vector<Tensor<1,spacedim> >  &exterior_forms,
@@ -373,9 +374,9 @@ MappingQ<dim,spacedim>::fill_fe_subface_values (const typename Triangulation<dim
 {
   // convert data object to internal data for this class. fails with an
   // exception if that is not possible
-  Assert (dynamic_cast<InternalData *> (&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *> (&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &> (mapping_data);
+  const InternalData &data = static_cast<const InternalData &> (mapping_data);
 
   // check whether this cell needs the full mapping or can be treated by a
   // reduced Q1 mapping, e.g. if the cell is entirely in the interior of the
@@ -388,11 +389,12 @@ MappingQ<dim,spacedim>::fill_fe_subface_values (const typename Triangulation<dim
 
   // depending on this result, use this or the other data object for the
   // mapping
-  typename MappingQ1<dim,spacedim>::InternalData *p_data=0;
-  if (data.use_mapping_q1_on_current_cell)
-    p_data=&data.mapping_q1_data;
-  else
-    p_data=&data;
+  const typename MappingQ1<dim,spacedim>::InternalData *p_data
+    = (data.use_mapping_q1_on_current_cell
+       ?
+       &data.mapping_q1_data
+       :
+       &data);
 
   const unsigned int n_q_points=q.size();
   this->compute_fill_face (cell, face_no, sub_no,
index f2b6ce4537a74cf2fbfd2497cf969ba1c3d9574f..a2d077e749584f4eae7a8b5ee41f9698a9ba51a2 100644 (file)
@@ -608,7 +608,7 @@ MappingQ1<dim,spacedim>::compute_fill (const typename Triangulation<dim,spacedim
                                        const unsigned int  n_q_points,
                                        const DataSetDescriptor  data_set,
                                        const CellSimilarity::Similarity cell_similarity,
-                                       InternalData  &data,
+                                       const InternalData  &data,
                                        std::vector<Point<spacedim> > &quadrature_points) const
 {
   const UpdateFlags update_flags(data.current_update_flags());
@@ -736,7 +736,7 @@ CellSimilarity::Similarity
 MappingQ1<dim,spacedim>::fill_fe_values (
   const typename Triangulation<dim,spacedim>::cell_iterator &cell,
   const Quadrature<dim>                        &q,
-  typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
   std::vector<Point<spacedim> >                &quadrature_points,
   std::vector<double>                          &JxW_values,
   std::vector<DerivativeForm<1,dim,spacedim> > &jacobians,
@@ -746,9 +746,9 @@ MappingQ1<dim,spacedim>::fill_fe_values (
   const CellSimilarity::Similarity                   cell_similarity) const
 {
   // ensure that the following static_cast is really correct:
-  Assert (dynamic_cast<InternalData *>(&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &>(mapping_data);
+  const InternalData &data = static_cast<const InternalData &>(mapping_data);
 
   const unsigned int n_q_points=q.size();
 
@@ -926,7 +926,7 @@ namespace internal
                        const unsigned int               subface_no,
                        const unsigned int               n_q_points,
                        const std::vector<double>        &weights,
-                       typename dealii::MappingQ1<dim,spacedim>::InternalData &data,
+                       const typename dealii::MappingQ1<dim,spacedim>::InternalData &data,
                        std::vector<double>              &JxW_values,
                        std::vector<Tensor<1,spacedim> > &boundary_forms,
                        std::vector<Point<spacedim> >    &normal_vectors,
@@ -1071,7 +1071,7 @@ MappingQ1<dim,spacedim>::compute_fill_face (
   const unsigned int               n_q_points,
   const DataSetDescriptor          data_set,
   const std::vector<double>        &weights,
-  InternalData                     &data,
+  const InternalData                     &data,
   std::vector<Point<spacedim> >    &quadrature_points,
   std::vector<double>              &JxW_values,
   std::vector<Tensor<1,spacedim> > &boundary_forms,
@@ -1095,7 +1095,7 @@ MappingQ1<dim,spacedim>::
 fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
                      const unsigned int                               face_no,
                      const Quadrature<dim-1>                          &q,
-                     typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                     const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                      std::vector<Point<spacedim> >                    &quadrature_points,
                      std::vector<double>                              &JxW_values,
                      std::vector<Tensor<1,spacedim> >                 &boundary_forms,
@@ -1105,9 +1105,9 @@ fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &
 {
   // ensure that the following cast
   // is really correct:
-  Assert (dynamic_cast<InternalData *>(&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &>(mapping_data);
+  const InternalData &data = static_cast<const InternalData &>(mapping_data);
 
   const unsigned int n_q_points = q.size();
 
@@ -1137,7 +1137,7 @@ fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterato
                         const unsigned int       face_no,
                         const unsigned int       sub_no,
                         const Quadrature<dim-1> &q,
-                        typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
+                        const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
                         std::vector<Point<spacedim> >     &quadrature_points,
                         std::vector<double>               &JxW_values,
                         std::vector<Tensor<1,spacedim> >  &boundary_forms,
@@ -1147,9 +1147,9 @@ fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterato
 {
   // ensure that the following cast
   // is really correct:
-  Assert (dynamic_cast<InternalData *>(&mapping_data) != 0,
+  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
           ExcInternalError());
-  InternalData &data = static_cast<InternalData &>(mapping_data);
+  const InternalData &data = static_cast<const InternalData &>(mapping_data);
 
   const unsigned int n_q_points = q.size();
 

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.