#include <deal.II/fe/fe_update_flags.h>
#include <deal.II/fe/fe_values_extractors.h>
#include <deal.II/fe/mapping.h>
+#include <deal.II/fe/mapping_related_data.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/vector.h>
#include <deal.II/base/config.h>
-#include <deal.II/base/derivative_form.h>
-#include <deal.II/base/point.h>
#include <deal.II/base/table.h>
#include <deal.II/base/tensor.h>
{
namespace FEValuesImplementation
{
- /**
- * A class that stores all of the mapping related data used in
- * dealii::FEValues, dealii::FEFaceValues, and dealii::FESubfaceValues
- * objects. Objects of this kind will be given as <i>output</i> argument
- * when dealii::FEValues::reinit() calls Mapping::fill_fe_values() for a
- * given cell, face, or subface.
- *
- * The data herein will then be provided as <i>input</i> argument in the
- * following call to FiniteElement::fill_fe_values().
- *
- * @ingroup feaccess
- */
- template <int dim, int spacedim = dim>
- class MappingRelatedData
- {
- public:
- /**
- * Initialize all vectors to correct size.
- */
- void
- initialize(const unsigned int n_quadrature_points,
- const UpdateFlags flags);
-
- /**
- * Compute and return an estimate for the memory consumption (in bytes)
- * of this object.
- */
- std::size_t
- memory_consumption() const;
-
- /**
- * Store an array of weights times the Jacobi determinant at the
- * quadrature points. This function is reset each time reinit() is
- * called. The Jacobi determinant is actually the reciprocal value of
- * the Jacobi matrices stored in this class, see the general
- * documentation of this class for more information.
- *
- * However, if this object refers to an FEFaceValues or FESubfaceValues
- * object, then the JxW_values correspond to the Jacobian of the
- * transformation of the face, not the cell, i.e. the dimensionality is
- * that of a surface measure, not of a volume measure. In this case, it
- * is computed from the boundary forms, rather than the Jacobian matrix.
- */
- std::vector<double> JxW_values;
-
- /**
- * Array of the Jacobian matrices at the quadrature points.
- */
- std::vector<DerivativeForm<1, dim, spacedim>> jacobians;
-
- /**
- * Array of the derivatives of the Jacobian matrices at the quadrature
- * points.
- */
- std::vector<DerivativeForm<2, dim, spacedim>> jacobian_grads;
-
- /**
- * Array of the inverse Jacobian matrices at the quadrature points.
- */
- std::vector<DerivativeForm<1, spacedim, dim>> inverse_jacobians;
-
- /**
- * Array of the derivatives of the Jacobian matrices at the quadrature
- * points, pushed forward to the real cell coordinates.
- */
- std::vector<Tensor<3, spacedim>> jacobian_pushed_forward_grads;
-
- /**
- * Array of the second derivatives of the Jacobian matrices at the
- * quadrature points.
- */
- std::vector<DerivativeForm<3, dim, spacedim>> jacobian_2nd_derivatives;
-
- /**
- * Array of the second derivatives of the Jacobian matrices at the
- * quadrature points, pushed forward to the real cell coordinates.
- */
- std::vector<Tensor<4, spacedim>> jacobian_pushed_forward_2nd_derivatives;
-
- /**
- * Array of the third derivatives of the Jacobian matrices at the
- * quadrature points.
- */
- std::vector<DerivativeForm<4, dim, spacedim>> jacobian_3rd_derivatives;
-
- /**
- * Array of the third derivatives of the Jacobian matrices at the
- * quadrature points, pushed forward to the real cell coordinates.
- */
- std::vector<Tensor<5, spacedim>> jacobian_pushed_forward_3rd_derivatives;
-
- /**
- * Array of quadrature points. This array is set up upon calling
- * reinit() and contains the quadrature points on the real element,
- * rather than on the reference element.
- */
- std::vector<Point<spacedim>> quadrature_points;
-
- /**
- * List of outward normal vectors at the quadrature points.
- */
- std::vector<Tensor<1, spacedim>> normal_vectors;
-
- /**
- * List of boundary forms at the quadrature points.
- */
- std::vector<Tensor<1, spacedim>> boundary_forms;
- };
-
-
/**
* A class that stores all of the shape function related data used in
* dealii::FEValues, dealii::FEFaceValues, and dealii::FESubfaceValues
#include <deal.II/fe/fe_update_flags.h>
#include <deal.II/fe/fe_values_extractors.h>
#include <deal.II/fe/mapping.h>
+#include <deal.II/fe/mapping_related_data.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/base/derivative_form.h>
#include <deal.II/fe/fe_update_flags.h>
+#include <deal.II/fe/mapping_related_data.h>
#include <deal.II/grid/tria.h>
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_fe_mapping_related_data_h
+#define dealii_fe_mapping_related_data_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/derivative_form.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/tensor.h>
+
+#include <deal.II/fe/fe_update_flags.h>
+
+#include <vector>
+
+
+DEAL_II_NAMESPACE_OPEN
+
+/** @addtogroup feaccess */
+/** @{ */
+
+namespace internal
+{
+ namespace FEValuesImplementation
+ {
+ /**
+ * A class that stores all of the mapping related data used in
+ * dealii::FEValues, dealii::FEFaceValues, and dealii::FESubfaceValues
+ * objects. Objects of this kind will be given as <i>output</i> argument
+ * when dealii::FEValues::reinit() calls Mapping::fill_fe_values() for a
+ * given cell, face, or subface.
+ *
+ * The data herein will then be provided as <i>input</i> argument in the
+ * following call to FiniteElement::fill_fe_values().
+ *
+ * @ingroup feaccess
+ */
+ template <int dim, int spacedim = dim>
+ class MappingRelatedData
+ {
+ public:
+ /**
+ * Initialize all vectors to correct size.
+ */
+ void
+ initialize(const unsigned int n_quadrature_points,
+ const UpdateFlags flags);
+
+ /**
+ * Compute and return an estimate for the memory consumption (in bytes)
+ * of this object.
+ */
+ std::size_t
+ memory_consumption() const;
+
+ /**
+ * Store an array of weights times the Jacobi determinant at the
+ * quadrature points. This function is reset each time reinit() is
+ * called. The Jacobi determinant is actually the reciprocal value of
+ * the Jacobi matrices stored in this class, see the general
+ * documentation of this class for more information.
+ *
+ * However, if this object refers to an FEFaceValues or FESubfaceValues
+ * object, then the JxW_values correspond to the Jacobian of the
+ * transformation of the face, not the cell, i.e. the dimensionality is
+ * that of a surface measure, not of a volume measure. In this case, it
+ * is computed from the boundary forms, rather than the Jacobian matrix.
+ */
+ std::vector<double> JxW_values;
+
+ /**
+ * Array of the Jacobian matrices at the quadrature points.
+ */
+ std::vector<DerivativeForm<1, dim, spacedim>> jacobians;
+
+ /**
+ * Array of the derivatives of the Jacobian matrices at the quadrature
+ * points.
+ */
+ std::vector<DerivativeForm<2, dim, spacedim>> jacobian_grads;
+
+ /**
+ * Array of the inverse Jacobian matrices at the quadrature points.
+ */
+ std::vector<DerivativeForm<1, spacedim, dim>> inverse_jacobians;
+
+ /**
+ * Array of the derivatives of the Jacobian matrices at the quadrature
+ * points, pushed forward to the real cell coordinates.
+ */
+ std::vector<Tensor<3, spacedim>> jacobian_pushed_forward_grads;
+
+ /**
+ * Array of the second derivatives of the Jacobian matrices at the
+ * quadrature points.
+ */
+ std::vector<DerivativeForm<3, dim, spacedim>> jacobian_2nd_derivatives;
+
+ /**
+ * Array of the second derivatives of the Jacobian matrices at the
+ * quadrature points, pushed forward to the real cell coordinates.
+ */
+ std::vector<Tensor<4, spacedim>> jacobian_pushed_forward_2nd_derivatives;
+
+ /**
+ * Array of the third derivatives of the Jacobian matrices at the
+ * quadrature points.
+ */
+ std::vector<DerivativeForm<4, dim, spacedim>> jacobian_3rd_derivatives;
+
+ /**
+ * Array of the third derivatives of the Jacobian matrices at the
+ * quadrature points, pushed forward to the real cell coordinates.
+ */
+ std::vector<Tensor<5, spacedim>> jacobian_pushed_forward_3rd_derivatives;
+
+ /**
+ * Array of quadrature points. This array is set up upon calling
+ * reinit() and contains the quadrature points on the real element,
+ * rather than on the reference element.
+ */
+ std::vector<Point<spacedim>> quadrature_points;
+
+ /**
+ * List of outward normal vectors at the quadrature points.
+ */
+ std::vector<Tensor<1, spacedim>> normal_vectors;
+
+ /**
+ * List of boundary forms at the quadrature points.
+ */
+ std::vector<Tensor<1, spacedim>> boundary_forms;
+ };
+ } // namespace FEValuesImplementation
+} // namespace internal
+
+
+/** @} */
+
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
#include <deal.II/fe/fe_update_flags.h>
#include <deal.II/fe/mapping_cartesian.h>
#include <deal.II/fe/mapping_q.h>
+#include <deal.II/fe/mapping_related_data.h>
#include <memory>
mapping_q.cc
mapping_q_cache.cc
mapping_manifold.cc
+ mapping_related_data.cc
)
set(_separate_src
mapping_q_eulerian.inst.in
mapping_q.inst.in
mapping_manifold.inst.in
+ mapping_related_data.inst.in
)
file(GLOB _header
{
namespace FEValuesImplementation
{
- template <int dim, int spacedim>
- void
- MappingRelatedData<dim, spacedim>::initialize(
- const unsigned int n_quadrature_points,
- const UpdateFlags flags)
- {
- if (flags & update_quadrature_points)
- this->quadrature_points.resize(
- n_quadrature_points,
- Point<spacedim>(numbers::signaling_nan<Tensor<1, spacedim>>()));
-
- if (flags & update_JxW_values)
- this->JxW_values.resize(n_quadrature_points,
- numbers::signaling_nan<double>());
-
- if (flags & update_jacobians)
- this->jacobians.resize(
- n_quadrature_points,
- numbers::signaling_nan<DerivativeForm<1, dim, spacedim>>());
-
- if (flags & update_jacobian_grads)
- this->jacobian_grads.resize(
- n_quadrature_points,
- numbers::signaling_nan<DerivativeForm<2, dim, spacedim>>());
-
- if (flags & update_jacobian_pushed_forward_grads)
- this->jacobian_pushed_forward_grads.resize(
- n_quadrature_points, numbers::signaling_nan<Tensor<3, spacedim>>());
-
- if (flags & update_jacobian_2nd_derivatives)
- this->jacobian_2nd_derivatives.resize(
- n_quadrature_points,
- numbers::signaling_nan<DerivativeForm<3, dim, spacedim>>());
-
- if (flags & update_jacobian_pushed_forward_2nd_derivatives)
- this->jacobian_pushed_forward_2nd_derivatives.resize(
- n_quadrature_points, numbers::signaling_nan<Tensor<4, spacedim>>());
-
- if (flags & update_jacobian_3rd_derivatives)
- this->jacobian_3rd_derivatives.resize(n_quadrature_points);
-
- if (flags & update_jacobian_pushed_forward_3rd_derivatives)
- this->jacobian_pushed_forward_3rd_derivatives.resize(
- n_quadrature_points, numbers::signaling_nan<Tensor<5, spacedim>>());
-
- if (flags & update_inverse_jacobians)
- this->inverse_jacobians.resize(
- n_quadrature_points,
- numbers::signaling_nan<DerivativeForm<1, spacedim, dim>>());
-
- if (flags & update_boundary_forms)
- this->boundary_forms.resize(
- n_quadrature_points, numbers::signaling_nan<Tensor<1, spacedim>>());
-
- if (flags & update_normal_vectors)
- this->normal_vectors.resize(
- n_quadrature_points, numbers::signaling_nan<Tensor<1, spacedim>>());
- }
-
-
-
- template <int dim, int spacedim>
- std::size_t
- MappingRelatedData<dim, spacedim>::memory_consumption() const
- {
- return (
- MemoryConsumption::memory_consumption(JxW_values) +
- MemoryConsumption::memory_consumption(jacobians) +
- MemoryConsumption::memory_consumption(jacobian_grads) +
- MemoryConsumption::memory_consumption(jacobian_pushed_forward_grads) +
- MemoryConsumption::memory_consumption(jacobian_2nd_derivatives) +
- MemoryConsumption::memory_consumption(
- jacobian_pushed_forward_2nd_derivatives) +
- MemoryConsumption::memory_consumption(jacobian_3rd_derivatives) +
- MemoryConsumption::memory_consumption(
- jacobian_pushed_forward_3rd_derivatives) +
- MemoryConsumption::memory_consumption(inverse_jacobians) +
- MemoryConsumption::memory_consumption(quadrature_points) +
- MemoryConsumption::memory_consumption(normal_vectors) +
- MemoryConsumption::memory_consumption(boundary_forms));
- }
-
-
-
template <int dim, int spacedim>
void
FiniteElementRelatedData<dim, spacedim>::initialize(
\{
namespace FEValuesImplementation
\{
- template class MappingRelatedData<deal_II_dimension,
- deal_II_space_dimension>;
template class FiniteElementRelatedData<deal_II_dimension,
deal_II_space_dimension>;
\}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/base/memory_consumption.h>
+#include <deal.II/base/numbers.h>
+#include <deal.II/base/signaling_nan.h>
+
+#include <deal.II/fe/mapping_related_data.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace internal
+{
+ namespace FEValuesImplementation
+ {
+ template <int dim, int spacedim>
+ void
+ MappingRelatedData<dim, spacedim>::initialize(
+ const unsigned int n_quadrature_points,
+ const UpdateFlags flags)
+ {
+ if (flags & update_quadrature_points)
+ this->quadrature_points.resize(
+ n_quadrature_points,
+ Point<spacedim>(numbers::signaling_nan<Tensor<1, spacedim>>()));
+
+ if (flags & update_JxW_values)
+ this->JxW_values.resize(n_quadrature_points,
+ numbers::signaling_nan<double>());
+
+ if (flags & update_jacobians)
+ this->jacobians.resize(
+ n_quadrature_points,
+ numbers::signaling_nan<DerivativeForm<1, dim, spacedim>>());
+
+ if (flags & update_jacobian_grads)
+ this->jacobian_grads.resize(
+ n_quadrature_points,
+ numbers::signaling_nan<DerivativeForm<2, dim, spacedim>>());
+
+ if (flags & update_jacobian_pushed_forward_grads)
+ this->jacobian_pushed_forward_grads.resize(
+ n_quadrature_points, numbers::signaling_nan<Tensor<3, spacedim>>());
+
+ if (flags & update_jacobian_2nd_derivatives)
+ this->jacobian_2nd_derivatives.resize(
+ n_quadrature_points,
+ numbers::signaling_nan<DerivativeForm<3, dim, spacedim>>());
+
+ if (flags & update_jacobian_pushed_forward_2nd_derivatives)
+ this->jacobian_pushed_forward_2nd_derivatives.resize(
+ n_quadrature_points, numbers::signaling_nan<Tensor<4, spacedim>>());
+
+ if (flags & update_jacobian_3rd_derivatives)
+ this->jacobian_3rd_derivatives.resize(n_quadrature_points);
+
+ if (flags & update_jacobian_pushed_forward_3rd_derivatives)
+ this->jacobian_pushed_forward_3rd_derivatives.resize(
+ n_quadrature_points, numbers::signaling_nan<Tensor<5, spacedim>>());
+
+ if (flags & update_inverse_jacobians)
+ this->inverse_jacobians.resize(
+ n_quadrature_points,
+ numbers::signaling_nan<DerivativeForm<1, spacedim, dim>>());
+
+ if (flags & update_boundary_forms)
+ this->boundary_forms.resize(
+ n_quadrature_points, numbers::signaling_nan<Tensor<1, spacedim>>());
+
+ if (flags & update_normal_vectors)
+ this->normal_vectors.resize(
+ n_quadrature_points, numbers::signaling_nan<Tensor<1, spacedim>>());
+ }
+
+
+
+ template <int dim, int spacedim>
+ std::size_t
+ MappingRelatedData<dim, spacedim>::memory_consumption() const
+ {
+ return (
+ MemoryConsumption::memory_consumption(JxW_values) +
+ MemoryConsumption::memory_consumption(jacobians) +
+ MemoryConsumption::memory_consumption(jacobian_grads) +
+ MemoryConsumption::memory_consumption(jacobian_pushed_forward_grads) +
+ MemoryConsumption::memory_consumption(jacobian_2nd_derivatives) +
+ MemoryConsumption::memory_consumption(
+ jacobian_pushed_forward_2nd_derivatives) +
+ MemoryConsumption::memory_consumption(jacobian_3rd_derivatives) +
+ MemoryConsumption::memory_consumption(
+ jacobian_pushed_forward_3rd_derivatives) +
+ MemoryConsumption::memory_consumption(inverse_jacobians) +
+ MemoryConsumption::memory_consumption(quadrature_points) +
+ MemoryConsumption::memory_consumption(normal_vectors) +
+ MemoryConsumption::memory_consumption(boundary_forms));
+ }
+ } // namespace FEValuesImplementation
+} // namespace internal
+
+#include "mapping_related_data.inst"
+
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 - 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef DOXYGEN
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
+ {
+# if deal_II_dimension <= deal_II_space_dimension
+ namespace internal
+ \{
+ namespace FEValuesImplementation
+ \{
+ template class MappingRelatedData<deal_II_dimension,
+ deal_II_space_dimension>;
+ \}
+ \}
+# endif
+ }
+#endif