From 7ef6020c627c0c68f42036c1443e798c47a8d1b4 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Sun, 26 Jul 2015 18:22:44 -0500 Subject: [PATCH] Revise the Mapping::fill_fe_*_values() interfaces. These functions previously got a whole lot of arguments for the output arrays. Consolidate them by simply passing the structure in which they actually lie. While there, also do: - Sort these arguments in such a way that they make a bit more sense. - Provide a thorough documentation of these functions. --- include/deal.II/fe/mapping.h | 190 ++++++++++++++++------- include/deal.II/fe/mapping_cartesian.h | 58 +++---- include/deal.II/fe/mapping_fe_field.h | 52 +++---- include/deal.II/fe/mapping_q.h | 53 +++---- include/deal.II/fe/mapping_q1.h | 53 +++---- include/deal.II/fe/mapping_q1_eulerian.h | 19 ++- include/deal.II/fe/mapping_q_eulerian.h | 19 ++- source/fe/fe_values.cc | 29 +--- source/fe/mapping_cartesian.cc | 134 +++++++--------- source/fe/mapping_fe_field.cc | 137 ++++++++-------- source/fe/mapping_q.cc | 105 ++++++------- source/fe/mapping_q1.cc | 136 ++++++++-------- source/fe/mapping_q1_eulerian.cc | 27 ++-- source/fe/mapping_q_eulerian.cc | 27 ++-- 14 files changed, 507 insertions(+), 532 deletions(-) diff --git a/include/deal.II/fe/mapping.h b/include/deal.II/fe/mapping.h index 2cbe97b02a..d4b01051b3 100644 --- a/include/deal.II/fe/mapping.h +++ b/include/deal.II/fe/mapping.h @@ -694,83 +694,153 @@ private: get_subface_data (const UpdateFlags flags, const Quadrature& quadrature) const = 0; - /** - * Fill the transformation fields of @p FEValues. Given a grid cell and the - * quadrature points on the unit cell, it computes all values specified by - * @p flags. The arrays to be filled have to have the correct size. + * Compute information about the mapping from the reference cell + * to the real cell indicated by the first argument to this function. + * Derived classes will have to implement this function based on the + * kind of mapping they represent. It is called by FEValues::reinit(). * - * Values are split into two groups: first, @p quadrature_points and @p - * JxW_values are filled with the quadrature rule transformed to the cell in - * physical space. + * Conceptually, this function's represents the application of the + * mapping $\mathbf x=\mathbf F_K(\hat {\mathbf x})$ from reference + * coordinates $\mathbf\in [0,1]^d$ to real space coordinates + * $\mathbf x$ for a given cell $K$. Its purpose is to compute the following + * kinds of data: * - * The second group contains the matrices needed to transform vector-valued - * functions, namely @p jacobians, the derivatives @p jacobian_grads, and - * the inverse operations in @p inverse_jacobians. - */ - /* virtual void */ - /* fill_fe_values (const typename Triangulation::cell_iterator &cell, */ - /* const Quadrature &quadrature, */ - /* InternalDataBase &internal, */ - /* std::vector > &quadrature_points, */ - /* std::vector &JxW_values) const = 0; */ - - /** - * The function above adjusted with the variable cell_normal_vectors for the - * case of codimension 1 + * - Data that results from the application of the mapping itself, e.g., + * computing the location $\mathbf x_q = \mathbf F_K(\hat{\mathbf x}_q)$ + * of quadrature points on the real cell, and that is directly useful + * to users of FEValues, for example during assembly. + * - Data that is necessary for finite element implementations to compute + * their shape functions on the real cell. To this end, the + * FEValues::reinit() function calls FiniteElement::fill_fe_values() + * after the current function, and the output of this function serves + * as input to FiniteElement::fill_fe_values(). Examples of + * information that needs to be computed here for use by the + * finite element classes is the Jacobian of the mapping, + * $\hat\nabla \mathbf F_K(\hat{\mathbf x})$ or its inverse, + * for example to transform the gradients of shape functions on + * the reference cell to the gradients of shape functions on + * the real cell. + * + * The information computed by this function is used to fill the various + * member variables of the output argument of this function. Which of + * the member variables of that structure should be filled is determined + * by the update flags stored in the Mapping::InternalDataBase object + * passed to this function. + * + * @param[in] cell The cell of the triangulation for which this function + * is to compute a mapping from the reference cell to. + * @param[in] cell_similarity Whether or not the cell given as first + * argument is simply a translation, rotation, etc of the cell for + * which this function was called the most recent time. This + * information is computed simply by matching the vertices (as stored + * by the Triangulation) between the previous and the current cell. + * The value passed here may be modified by implementations of + * this function and should then be returned (see the discussion of the + * return value of this function). + * @param[in] quadrature A reference to the quadrature formula in use + * for the current evaluation. This quadrature object is the same + * as the one used when creating the @p internal_data object + * @param[in] internal_data A reference to an object previously + * created by get_data() and that may be used to store information + * the mapping can compute once on the reference cell. See the + * documentation of the Mapping::InternalDataBase class for an + * extensive description of the purpose of these objects. + * @param[out] output_data A reference to an object whose member + * variables should be computed. Not all of the members of this + * argument need to be filled; which ones need to be filled is + * determined by the update flags stored inside the + * @p internal_data object. + * @return An updated value of the @p cell_similarity argument to + * this function. The returned value will be used for the corresponding + * argument when FEValues::reinit() calls + * FiniteElement::fill_fe_values(). In most cases, derived classes will + * simply want to return the value passed for @p cell_similarity. + * However, implementations of this function may downgrade the + * level of cell similarity. This is, for example, the case for + * classes that take not only into account the locations of the + * vertices of a cell (as reported by the Triangulation), but also + * other information specific to the mapping. The purpose is that + * FEValues::reinit() can compute whether a cell is similar to the + * previous one only based on the cell's vertices, whereas the + * mapping may also consider displacement fields (e.g., in the + * MappingQ1Eulerian and MappingFEField classes). In such cases, + * the mapping may conclude that the previously computed + * cell similarity is too optimistic, and invalidate it for + * subsequent use in FiniteElement::fill_fe_values() by + * returning a less optimistic cell similarity value. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, - const Quadrature &quadrature, - const InternalDataBase &internal, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const=0; - - + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const InternalDataBase &internal_data, + FEValuesData &output_data) const = 0; /** - * Performs the same as @p fill_fe_values on a face. Additionally, @p - * boundary_form (see - * @ref GlossBoundaryForm - * ) and @p normal_vectors can be computed on surfaces. Since the boundary - * form already contains the determinant of the Jacobian of the - * transformation, it is sometimes more economic to use the boundary form - * instead of the product of the unit normal and the transformed quadrature - * weight. + * This function is the equivalent to Mapping::fill_fe_values(), + * but for faces of cells. See there for an extensive discussion + * of its purpose. It is called by FEFaceValues::reinit(). + * + * @param[in] cell The cell of the triangulation for which this function + * is to compute a mapping from the reference cell to. + * @param[in] face_no The number of the face of the given cell for which + * information is requested. + * @param[in] quadrature A reference to the quadrature formula in use + * for the current evaluation. This quadrature object is the same + * as the one used when creating the @p internal_data object + * @param[in] internal_data A reference to an object previously + * created by get_data() and that may be used to store information + * the mapping can compute once on the reference cell. See the + * documentation of the Mapping::InternalDataBase class for an + * extensive description of the purpose of these objects. + * @param[out] output_data A reference to an object whose member + * variables should be computed. Not all of the members of this + * argument need to be filled; which ones need to be filled is + * determined by the update flags stored inside the + * @p internal_data object. */ virtual void fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &quadrature, - const InternalDataBase &internal, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_form, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const = 0; + const unsigned int face_no, + const Quadrature &quadrature, + const InternalDataBase &internal_data, + FEValuesData &output_data) const = 0; /** - * See above. + * This function is the equivalent to Mapping::fill_fe_values(), + * but for subfaces (i.e., children of faces) of cells. + * See there for an extensive discussion + * of its purpose. It is called by FESubfaceValues::reinit(). + * + * @param[in] cell The cell of the triangulation for which this function + * is to compute a mapping from the reference cell to. + * @param[in] face_no The number of the face of the given cell for which + * information is requested. + * @param[in] subface_no The number of the child of a face of the + * given cell for which information is requested. + * @param[in] quadrature A reference to the quadrature formula in use + * for the current evaluation. This quadrature object is the same + * as the one used when creating the @p internal_data object + * @param[in] internal_data A reference to an object previously + * created by get_data() and that may be used to store information + * the mapping can compute once on the reference cell. See the + * documentation of the Mapping::InternalDataBase class for an + * extensive description of the purpose of these objects. + * @param[out] output_data A reference to an object whose member + * variables should be computed. Not all of the members of this + * argument need to be filled; which ones need to be filled is + * determined by the update flags stored inside the + * @p internal_data object. */ virtual void fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature &quadrature, - const InternalDataBase &internal, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_form, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const = 0; + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const InternalDataBase &internal_data, + FEValuesData &output_data) const = 0; /** * Give class @p FEValues access to the private get_...data and diff --git a/include/deal.II/fe/mapping_cartesian.h b/include/deal.II/fe/mapping_cartesian.h index ba9bf2d63c..f0ae837557 100644 --- a/include/deal.II/fe/mapping_cartesian.h +++ b/include/deal.II/fe/mapping_cartesian.h @@ -69,43 +69,43 @@ public: get_subface_data (const UpdateFlags flags, const Quadrature& quadrature) const; + /** + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. + */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, - const Quadrature &quadrature, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &, - const CellSimilarity::Similarity cell_similarity) const; - + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; + /** + * Compute mapping-related information for a face of a cell. + * See the documentation of Mapping::fill_fe_face_values() for + * a discussion of purpose and arguments of this function. + */ virtual void fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_form, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; + + /** + * Compute mapping-related information for a child of a face of a cell. + * See the documentation of Mapping::fill_fe_subface_values() for + * a discussion of purpose and arguments of this function. + */ virtual void fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_form, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; virtual void transform (const VectorSlice > > input, diff --git a/include/deal.II/fe/mapping_fe_field.h b/include/deal.II/fe/mapping_fe_field.h index d982a9029d..070a8fcae1 100644 --- a/include/deal.II/fe/mapping_fe_field.h +++ b/include/deal.II/fe/mapping_fe_field.h @@ -434,52 +434,42 @@ public: protected: /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, const Quadrature &quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const; + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a face of a cell. + * See the documentation of Mapping::fill_fe_face_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a child of a face of a cell. + * See the documentation of Mapping::fill_fe_subface_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; - + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** * This function and the next allow to generate the transform require by the diff --git a/include/deal.II/fe/mapping_q.h b/include/deal.II/fe/mapping_q.h index 71e515ede5..334dc97028 100644 --- a/include/deal.II/fe/mapping_q.h +++ b/include/deal.II/fe/mapping_q.h @@ -196,51 +196,42 @@ public: protected: /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, - const Quadrature &quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const; + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a face of a cell. + * See the documentation of Mapping::fill_fe_face_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - typename std::vector > &exterior_form, - typename std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a child of a face of a cell. + * See the documentation of Mapping::fill_fe_subface_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - typename std::vector > &exterior_form, - typename std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** * For dim=2,3. Append the support points of all shape functions diff --git a/include/deal.II/fe/mapping_q1.h b/include/deal.II/fe/mapping_q1.h index 6e4395bed6..d3eb9006a9 100644 --- a/include/deal.II/fe/mapping_q1.h +++ b/include/deal.II/fe/mapping_q1.h @@ -318,51 +318,42 @@ public: DataSetDescriptor; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, - const Quadrature &quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const; + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a face of a cell. + * See the documentation of Mapping::fill_fe_face_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - typename std::vector > &boundary_form, - typename std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** - * Implementation of the interface in Mapping. + * Compute mapping-related information for a child of a face of a cell. + * See the documentation of Mapping::fill_fe_subface_values() for + * a discussion of purpose and arguments of this function. */ virtual void fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature& quadrature, - const typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - typename std::vector > &boundary_form, - typename std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const; + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** * Compute shape values and/or derivatives. diff --git a/include/deal.II/fe/mapping_q1_eulerian.h b/include/deal.II/fe/mapping_q1_eulerian.h index f3b2da3097..8b31afdc44 100644 --- a/include/deal.II/fe/mapping_q1_eulerian.h +++ b/include/deal.II/fe/mapping_q1_eulerian.h @@ -122,21 +122,20 @@ public: protected: /** - * Implementation of the interface in MappingQ1. Overrides the function in - * the base class, since we cannot use any cell similarity for this class. + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. + * + * This function overrides the function in + * the base class since we cannot use any cell similarity for this class. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, const Quadrature &quadrature, - typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const; + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** * Reference to the vector of shifts. diff --git a/include/deal.II/fe/mapping_q_eulerian.h b/include/deal.II/fe/mapping_q_eulerian.h index e03f1f5035..34aecd29a3 100644 --- a/include/deal.II/fe/mapping_q_eulerian.h +++ b/include/deal.II/fe/mapping_q_eulerian.h @@ -137,21 +137,20 @@ public: protected: /** - * Implementation of the interface in MappingQ. Overrides the function in - * the base class, since we cannot use any cell similarity for this class. + * Compute mapping-related information for a cell. + * See the documentation of Mapping::fill_fe_values() for + * a discussion of purpose, arguments, and return value of this function. + * + * This function overrides the function in + * the base class since we cannot use any cell similarity for this class. */ virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, const Quadrature &quadrature, - typename Mapping::InternalDataBase &mapping_data, - typename std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &cell_normal_vectors, - const CellSimilarity::Similarity cell_similarity) const; + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const; /** * Reference to the vector of shifts. diff --git a/source/fe/fe_values.cc b/source/fe/fe_values.cc index 2cbedfbbba..909f6d11d3 100644 --- a/source/fe/fe_values.cc +++ b/source/fe/fe_values.cc @@ -3518,15 +3518,10 @@ void FEValues::do_reinit () // it this->cell_similarity = this->get_mapping().fill_fe_values(*this->present_cell, + this->cell_similarity, quadrature, *this->mapping_data, - this->quadrature_points, - this->JxW_values, - this->jacobians, - this->jacobian_grads, - this->inverse_jacobians, - this->normal_vectors, - this->cell_similarity); + *this); // then call the finite element and, with the data // already filled by the mapping, let it compute the @@ -3722,15 +3717,11 @@ void FEFaceValues::do_reinit (const unsigned int face_no) Assert(!(this->update_flags & update_jacobian_grads), ExcNotImplemented()); - this->get_mapping().fill_fe_face_values(*this->present_cell, face_no, + this->get_mapping().fill_fe_face_values(*this->present_cell, + face_no, this->quadrature, *this->mapping_data, - this->quadrature_points, - this->JxW_values, - this->boundary_forms, - this->normal_vectors, - this->jacobians, - this->inverse_jacobians); + *this); this->get_fe().fill_fe_face_values(this->get_mapping(), *this->present_cell, face_no, @@ -3958,15 +3949,11 @@ void FESubfaceValues::do_reinit (const unsigned int face_no, // now ask the mapping and the finite element to do the actual work this->get_mapping().fill_fe_subface_values(*this->present_cell, - face_no, subface_no, + face_no, + subface_no, this->quadrature, *this->mapping_data, - this->quadrature_points, - this->JxW_values, - this->boundary_forms, - this->normal_vectors, - this->jacobians, - this->inverse_jacobians); + *this); this->get_fe().fill_fe_subface_values(this->get_mapping(), *this->present_cell, diff --git a/source/fe/mapping_cartesian.cc b/source/fe/mapping_cartesian.cc index a2ae39b483..405d18b49e 100644 --- a/source/fe/mapping_cartesian.cc +++ b/source/fe/mapping_cartesian.cc @@ -318,28 +318,23 @@ template CellSimilarity::Similarity MappingCartesian:: fill_fe_values (const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector< DerivativeForm<1,dim,spacedim> > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &, - const CellSimilarity::Similarity cell_similarity) const + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal // data for this class. fails with // an exception if that is not // possible - Assert (dynamic_cast (&mapping_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); + const InternalData &data = static_cast (internal_data); std::vector > dummy; compute_fill (cell, invalid_face_number, invalid_face_number, cell_similarity, data, - quadrature_points, + output_data.quadrature_points, dummy); // compute Jacobian @@ -355,36 +350,37 @@ fill_fe_values (const typename Triangulation::cell_iterator &cell, J *= data.length[d]; data.volume_element = J; if (data.current_update_flags() & update_JxW_values) - for (unsigned int i=0; i(); + output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>(); for (unsigned int j=0; j(); + for (unsigned int i=0; i(); + // "compute" inverse Jacobian at the // quadrature points, which are all // the same if (data.current_update_flags() & update_inverse_jacobians) if (cell_similarity != CellSimilarity::translation) - for (unsigned int i=0; i(); + output_data.inverse_jacobians[i]=Tensor<2,dim>(); for (unsigned int j=0; j::cell_iterator &cell, template void -MappingCartesian::fill_fe_face_values ( - const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingCartesian:: +fill_fe_face_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal // data for this class. fails with // an exception if that is not // possible - Assert (dynamic_cast (&mapping_data) != 0, + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + const InternalData &data = static_cast (internal_data); compute_fill (cell, face_no, invalid_face_number, CellSimilarity::none, data, - quadrature_points, - normal_vectors); + output_data.quadrature_points, + output_data.normal_vectors); // first compute Jacobian // determinant, which is simply the @@ -430,12 +421,12 @@ MappingCartesian::fill_fe_face_values ( J *= data.length[d]; if (data.current_update_flags() & update_JxW_values) - for (unsigned int i=0; i::fill_fe_face_values ( } if (data.current_update_flags() & update_jacobians) - for (unsigned int i=0; i(); + output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>(); for (unsigned int d=0; d(); + output_data.inverse_jacobians[i] = DerivativeForm<1,dim,spacedim>(); for (unsigned int d=0; d::fill_fe_face_values ( template void -MappingCartesian::fill_fe_subface_values ( - const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingCartesian:: +fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal // data for this class. fails with // an exception if that is not // possible - Assert (dynamic_cast (&mapping_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); + const InternalData &data = static_cast (internal_data); - compute_fill (cell, face_no, sub_no, CellSimilarity::none, + compute_fill (cell, face_no, subface_no, CellSimilarity::none, data, - quadrature_points, - normal_vectors); + output_data.quadrature_points, + output_data.normal_vectors); // first compute Jacobian // determinant, which is simply the @@ -516,13 +502,13 @@ MappingCartesian::fill_fe_subface_values ( cell->face(face_no)->has_children() ? cell->face(face_no)->n_children() : GeometryInfo::max_children_per_face; - for (unsigned int i=0; i::fill_fe_subface_values ( } if (data.current_update_flags() & update_jacobians) - for (unsigned int i=0; i(); + output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>(); for (unsigned int d=0; d(); + output_data.inverse_jacobians[i] = DerivativeForm<1,spacedim,dim>(); for (unsigned int d=0; d::get_subface_data (const UpdateFlags upda // recalculate data even when cells are similar. template CellSimilarity::Similarity -MappingFEField::fill_fe_values ( - const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &normal_vectors, - const CellSimilarity::Similarity cell_similarity) const +MappingFEField:: +fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); + const InternalData &data = static_cast (internal_data); - const unsigned int n_q_points=q.size(); + const unsigned int n_q_points=quadrature.size(); const CellSimilarity::Similarity updated_cell_similarity = (get_degree() == 1 ? @@ -401,21 +396,22 @@ MappingFEField::fill_fe_values ( compute_fill (cell, n_q_points, QProjector::DataSetDescriptor::cell (), updated_cell_similarity, - data, quadrature_points); + data, + output_data.quadrature_points); const UpdateFlags update_flags(data.current_update_flags()); - const std::vector &weights=q.get_weights(); + const std::vector &weights=quadrature.get_weights(); // Multiply quadrature weights by absolute value of Jacobian determinants or // the area element g=sqrt(DX^t DX) in case of codim > 0 if (update_flags & (update_normal_vectors | update_JxW_values)) { - AssertDimension (JxW_values.size(), n_q_points); + AssertDimension (output_data.JxW_values.size(), n_q_points); Assert( !(update_flags & update_normal_vectors ) || - (normal_vectors.size() == n_q_points), - ExcDimensionMismatch(normal_vectors.size(), n_q_points)); + (output_data.normal_vectors.size() == n_q_points), + ExcDimensionMismatch(output_data.normal_vectors.size(), n_q_points)); if (cell_similarity != CellSimilarity::translation) @@ -433,7 +429,7 @@ MappingFEField::fill_fe_values ( Assert (det > 1e-12*Utilities::fixed_power(cell->diameter()/ std::sqrt(double(dim))), (typename Mapping::ExcDistortedMappedCell(cell->center(), det, point))); - JxW_values[point] = weights[point] * det; + output_data.JxW_values[point] = weights[point] * det; } // if dim==spacedim, then there is no cell normal to // compute. since this is for FEValues (and not FEFaceValues), @@ -450,13 +446,13 @@ MappingFEField::fill_fe_values ( for (unsigned int j=0; j::fill_fe_values ( ExcMessage("There is no cell normal in codim 2.")); if (dim==1) - cross_product(normal_vectors[point], -DX_t[0]); + cross_product(output_data.normal_vectors[point], -DX_t[0]); else //dim == 2 - cross_product(normal_vectors[point],DX_t[0],DX_t[1]); + cross_product(output_data.normal_vectors[point],DX_t[0],DX_t[1]); - normal_vectors[point] /= normal_vectors[point].norm(); + output_data.normal_vectors[point] /= output_data.normal_vectors[point].norm(); if (cell->direction_flag() == false) - normal_vectors[point] *= -1.; + output_data.normal_vectors[point] *= -1.; } } @@ -486,10 +482,10 @@ MappingFEField::fill_fe_values ( // copy values from InternalData to vector given by reference if (update_flags & update_jacobians) { - AssertDimension (jacobians.size(), n_q_points); + AssertDimension (output_data.jacobians.size(), n_q_points); if (cell_similarity != CellSimilarity::translation) for (unsigned int point=0; point::fill_fe_values ( // we only do it for cells, not faces. if (update_flags & update_jacobian_grads) { - AssertDimension (jacobian_grads.size(), n_q_points); + AssertDimension (output_data.jacobian_grads.size(), n_q_points); if (cell_similarity != CellSimilarity::translation) { - std::fill(jacobian_grads.begin(), - jacobian_grads.end(), + std::fill(output_data.jacobian_grads.begin(), + output_data.jacobian_grads.end(), DerivativeForm<2,dim,spacedim>()); const unsigned int data_set = QProjector::DataSetDescriptor::cell(); @@ -529,7 +525,7 @@ MappingFEField::fill_fe_values ( for (unsigned int i=0; i::fill_fe_values ( // copy values from InternalData to vector given by reference if (update_flags & update_inverse_jacobians) { - AssertDimension (inverse_jacobians.size(), n_q_points); + AssertDimension (output_data.inverse_jacobians.size(), n_q_points); if (cell_similarity != CellSimilarity::translation) for (unsigned int point=0; point::fill_fe_values ( template void -MappingFEField::fill_fe_face_values ( - const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingFEField:: +fill_fe_face_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + const InternalData &data = static_cast (internal_data); - const unsigned int n_q_points=q.size(); + const unsigned int n_q_points=quadrature.size(); this->compute_fill_face (cell, face_no, numbers::invalid_unsigned_int, n_q_points, QProjector::DataSetDescriptor:: @@ -578,49 +569,51 @@ MappingFEField::fill_fe_face_values ( cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points), - q.get_weights(), + quadrature.get_weights(), data, - quadrature_points, JxW_values, - exterior_forms, normal_vectors, jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } template void -MappingFEField::fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingFEField:: +fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + const InternalData &data = static_cast (internal_data); - const unsigned int n_q_points=q.size(); - this->compute_fill_face (cell, face_no, sub_no, + const unsigned int n_q_points=quadrature.size(); + this->compute_fill_face (cell, face_no, subface_no, n_q_points, QProjector::DataSetDescriptor:: - subface (face_no, sub_no, + subface (face_no, subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points, cell->subface_case(face_no)), - q.get_weights(), + quadrature.get_weights(), data, - quadrature_points, JxW_values, - exterior_forms, normal_vectors, jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } diff --git a/source/fe/mapping_q.cc b/source/fe/mapping_q.cc index 048d611715..f3bcb5394f 100644 --- a/source/fe/mapping_q.cc +++ b/source/fe/mapping_q.cc @@ -27,6 +27,7 @@ #include #include #include +#include #include @@ -247,22 +248,17 @@ MappingQ::get_subface_data (const UpdateFlags update_flags, // recalculate data even when cells are similar. template CellSimilarity::Similarity -MappingQ::fill_fe_values ( - const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &normal_vectors, - const CellSimilarity::Similarity cell_similarity) const +MappingQ:: +fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); + const InternalData &data = static_cast (internal_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 @@ -292,11 +288,11 @@ MappingQ::fill_fe_values ( CellSimilarity::invalid_next_cell : cell_similarity); - MappingQ1::fill_fe_values(cell, q, *p_data, - quadrature_points, JxW_values, - jacobians, jacobian_grads, inverse_jacobians, - normal_vectors, - updated_cell_similarity); + MappingQ1::fill_fe_values(cell, + updated_cell_similarity, + quadrature, + *p_data, + output_data); return updated_cell_similarity; } @@ -305,23 +301,18 @@ MappingQ::fill_fe_values ( template void -MappingQ::fill_fe_face_values ( - const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingQ:: +fill_fe_face_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + const InternalData &data = static_cast (internal_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 @@ -341,7 +332,7 @@ MappingQ::fill_fe_face_values ( : &data); - const unsigned int n_q_points=q.size(); + const unsigned int n_q_points = quadrature.size(); this->compute_fill_face (cell, face_no, numbers::invalid_unsigned_int, n_q_points, QProjector::DataSetDescriptor:: @@ -350,33 +341,32 @@ MappingQ::fill_fe_face_values ( cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points), - q.get_weights(), + quadrature.get_weights(), *p_data, - quadrature_points, JxW_values, - exterior_forms, normal_vectors, jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } template void -MappingQ::fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &exterior_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const +MappingQ:: +fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // convert data object to internal data for this class. fails with an // exception if that is not possible - Assert (dynamic_cast (&mapping_data) != 0, + Assert (dynamic_cast (&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast (mapping_data); + const InternalData &data = static_cast (internal_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 @@ -396,21 +386,24 @@ MappingQ::fill_fe_subface_values (const typename Triangulationcompute_fill_face (cell, face_no, sub_no, + const unsigned int n_q_points = quadrature.size(); + this->compute_fill_face (cell, face_no, subface_no, n_q_points, QProjector::DataSetDescriptor:: - subface (face_no, sub_no, + subface (face_no, subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points, cell->subface_case(face_no)), - q.get_weights(), + quadrature.get_weights(), *p_data, - quadrature_points, JxW_values, - exterior_forms, normal_vectors, jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } diff --git a/source/fe/mapping_q1.cc b/source/fe/mapping_q1.cc index 9aca842d58..c4cb6be7fb 100644 --- a/source/fe/mapping_q1.cc +++ b/source/fe/mapping_q1.cc @@ -734,31 +734,27 @@ MappingQ1::compute_mapping_support_points( template CellSimilarity::Similarity -MappingQ1::fill_fe_values ( - const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &normal_vectors, - const CellSimilarity::Similarity cell_similarity) const +MappingQ1:: +fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity cell_similarity, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // ensure that the following static_cast is really correct: - Assert (dynamic_cast(&mapping_data) != 0, + Assert (dynamic_cast(&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast(mapping_data); + const InternalData &data = static_cast(internal_data); - const unsigned int n_q_points=q.size(); + const unsigned int n_q_points=quadrature.size(); compute_fill (cell, n_q_points, DataSetDescriptor::cell (), cell_similarity, - data, quadrature_points); + data, + output_data.quadrature_points); const UpdateFlags update_flags(data.current_update_flags()); - const std::vector &weights=q.get_weights(); + const std::vector &weights=quadrature.get_weights(); // Multiply quadrature weights by absolute value of Jacobian determinants or // the area element g=sqrt(DX^t DX) in case of codim > 0 @@ -766,11 +762,11 @@ MappingQ1::fill_fe_values ( if (update_flags & (update_normal_vectors | update_JxW_values)) { - AssertDimension (JxW_values.size(), n_q_points); + AssertDimension (output_data.JxW_values.size(), n_q_points); Assert( !(update_flags & update_normal_vectors ) || - (normal_vectors.size() == n_q_points), - ExcDimensionMismatch(normal_vectors.size(), n_q_points)); + (output_data.normal_vectors.size() == n_q_points), + ExcDimensionMismatch(output_data.normal_vectors.size(), n_q_points)); if (cell_similarity != CellSimilarity::translation) @@ -790,7 +786,7 @@ MappingQ1::fill_fe_values ( std::sqrt(double(dim))), (typename Mapping::ExcDistortedMappedCell(cell->center(), det, point))); - JxW_values[point] = weights[point] * det; + output_data.JxW_values[point] = weights[point] * det; } // if dim==spacedim, then there is no cell normal to // compute. since this is for FEValues (and not FEFaceValues), @@ -807,14 +803,14 @@ MappingQ1::fill_fe_values ( for (unsigned int j=0; j::fill_fe_values ( Assert( codim==1 , ExcMessage("There is no cell normal in codim 2.")); if (dim==1) - cross_product(normal_vectors[point], + cross_product(output_data.normal_vectors[point], -DX_t[0]); else //dim == 2 - cross_product(normal_vectors[point],DX_t[0],DX_t[1]); + cross_product(output_data.normal_vectors[point],DX_t[0],DX_t[1]); - normal_vectors[point] /= normal_vectors[point].norm(); + output_data.normal_vectors[point] /= output_data.normal_vectors[point].norm(); if (cell->direction_flag() == false) - normal_vectors[point] *= -1.; + output_data.normal_vectors[point] *= -1.; } } @@ -848,23 +844,23 @@ MappingQ1::fill_fe_values ( // copy values from InternalData to vector given by reference if (update_flags & update_jacobians) { - AssertDimension (jacobians.size(), n_q_points); + AssertDimension (output_data.jacobians.size(), n_q_points); if (cell_similarity != CellSimilarity::translation) for (unsigned int point=0; point()); @@ -894,17 +890,17 @@ MappingQ1::fill_fe_values ( for (unsigned int i=0; i void MappingQ1:: fill_fe_face_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const + const unsigned int face_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // ensure that the following cast // is really correct: - Assert (dynamic_cast(&mapping_data) != 0, + Assert (dynamic_cast(&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast(mapping_data); + const InternalData &data = static_cast(internal_data); - const unsigned int n_q_points = q.size(); + const unsigned int n_q_points = quadrature.size(); compute_fill_face (cell, face_no, numbers::invalid_unsigned_int, n_q_points, @@ -1119,14 +1110,14 @@ fill_fe_face_values (const typename Triangulation::cell_iterator & cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points), - q.get_weights(), + quadrature.get_weights(), data, - quadrature_points, - JxW_values, - boundary_forms, - normal_vectors, - jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } @@ -1135,41 +1126,36 @@ template void MappingQ1:: fill_fe_subface_values (const typename Triangulation::cell_iterator &cell, - const unsigned int face_no, - const unsigned int sub_no, - const Quadrature &q, - const typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &boundary_forms, - std::vector > &normal_vectors, - std::vector > &jacobians, - std::vector > &inverse_jacobians) const + const unsigned int face_no, + const unsigned int subface_no, + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // ensure that the following cast // is really correct: - Assert (dynamic_cast(&mapping_data) != 0, + Assert (dynamic_cast(&internal_data) != 0, ExcInternalError()); - const InternalData &data = static_cast(mapping_data); + const InternalData &data = static_cast(internal_data); - const unsigned int n_q_points = q.size(); + const unsigned int n_q_points = quadrature.size(); - compute_fill_face (cell, face_no, sub_no, + compute_fill_face (cell, face_no, subface_no, n_q_points, - DataSetDescriptor::subface (face_no, sub_no, + DataSetDescriptor::subface (face_no, subface_no, cell->face_orientation(face_no), cell->face_flip(face_no), cell->face_rotation(face_no), n_q_points, cell->subface_case(face_no)), - q.get_weights(), + quadrature.get_weights(), data, - quadrature_points, - JxW_values, - boundary_forms, - normal_vectors, - jacobians, - inverse_jacobians); + output_data.quadrature_points, + output_data.JxW_values, + output_data.boundary_forms, + output_data.normal_vectors, + output_data.jacobians, + output_data.inverse_jacobians); } diff --git a/source/fe/mapping_q1_eulerian.cc b/source/fe/mapping_q1_eulerian.cc index dc7a768dd8..0903bd44da 100644 --- a/source/fe/mapping_q1_eulerian.cc +++ b/source/fe/mapping_q1_eulerian.cc @@ -112,26 +112,21 @@ MappingQ1Eulerian::clone () const template CellSimilarity::Similarity -MappingQ1Eulerian::fill_fe_values ( - const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector< DerivativeForm<1,dim,spacedim> > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &normal_vectors, - const CellSimilarity::Similarity) const +MappingQ1Eulerian:: +fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity , + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // call the function of the base class, but ignoring // any potentially detected cell similarity between // the current and the previous cell - MappingQ1::fill_fe_values (cell, q, mapping_data, - quadrature_points, JxW_values, jacobians, - jacobian_grads, inverse_jacobians, - normal_vectors, - CellSimilarity::invalid_next_cell); + MappingQ1::fill_fe_values (cell, + CellSimilarity::invalid_next_cell, + quadrature, + internal_data, + output_data); // also return the updated flag since any detected // similarity wasn't based on the mapped field, but // the original vertices which are meaningless diff --git a/source/fe/mapping_q_eulerian.cc b/source/fe/mapping_q_eulerian.cc index 69444b6994..47f74528c6 100644 --- a/source/fe/mapping_q_eulerian.cc +++ b/source/fe/mapping_q_eulerian.cc @@ -177,26 +177,21 @@ compute_mapping_support_points template CellSimilarity::Similarity -MappingQEulerian::fill_fe_values ( - const typename Triangulation::cell_iterator &cell, - const Quadrature &q, - typename Mapping::InternalDataBase &mapping_data, - std::vector > &quadrature_points, - std::vector &JxW_values, - std::vector > &jacobians, - std::vector > &jacobian_grads, - std::vector > &inverse_jacobians, - std::vector > &normal_vectors, - const CellSimilarity::Similarity ) const +MappingQEulerian:: +fill_fe_values (const typename Triangulation::cell_iterator &cell, + const CellSimilarity::Similarity , + const Quadrature &quadrature, + const typename Mapping::InternalDataBase &internal_data, + FEValuesData &output_data) const { // call the function of the base class, but ignoring // any potentially detected cell similarity between // the current and the previous cell - MappingQ::fill_fe_values (cell, q, mapping_data, - quadrature_points, JxW_values, jacobians, - jacobian_grads, inverse_jacobians, - normal_vectors, - CellSimilarity::invalid_next_cell); + MappingQ::fill_fe_values (cell, + CellSimilarity::invalid_next_cell, + quadrature, + internal_data, + output_data); // also return the updated flag since any detected // similarity wasn't based on the mapped field, but // the original vertices which are meaningless -- 2.39.5