36#include <boost/container/small_vector.hpp>
49template <
int dim,
int spacedim>
61template <
int dim,
int spacedim>
79template <
int dim,
int spacedim>
88 const unsigned int n_q_points = quadrature.
size();
106 const std::array<Quadrature<1>, dim> &quad_array =
110 if (quad_array[i - 1].size() != quad_array[i].size())
117 const std::vector<Point<1>> &points_1 =
118 quad_array[i - 1].get_points();
119 const std::vector<Point<1>> &points_2 =
120 quad_array[i].get_points();
121 const std::vector<double> &weights_1 =
122 quad_array[i - 1].get_weights();
123 const std::vector<double> &weights_2 =
124 quad_array[i].get_weights();
125 for (
unsigned int j = 0; j < quad_array[i].size(); ++j)
127 if (
std::abs(points_1[j][0] - points_2[j][0]) > 1.e-10 ||
128 std::abs(weights_1[j] - weights_2[j]) > 1.e-10)
130 tensor_product_quadrature =
false;
137 if (tensor_product_quadrature)
144 shape_info.lexicographic_numbering =
147 shape_info.n_q_points = n_q_points;
148 shape_info.dofs_per_component_on_cell =
157template <
int dim,
int spacedim>
160 const UpdateFlags update_flags,
162 const unsigned int n_original_q_points)
164 reinit(update_flags, quadrature);
168 if (dim > 1 && tensor_product_quadrature)
170 constexpr unsigned int facedim = dim - 1;
173 shape_info.lexicographic_numbering =
176 shape_info.n_q_points = n_original_q_points;
177 shape_info.dofs_per_component_on_cell =
183 if (this->update_each &
187 aux[0].resize(n_original_q_points);
189 aux[1].resize(n_original_q_points);
194 unit_tangentials[i].resize(n_original_q_points);
195 std::fill(unit_tangentials[i].begin(),
196 unit_tangentials[i].end(),
201 .resize(n_original_q_points);
216template <
int dim,
int spacedim>
224 FETools::lexicographic_to_hierarchic_numbering<dim>(p))
226 internal::MappingQImplementation::unit_support_points<dim>(
231 compute_support_point_weights_perimeter_to_interior(
235 internal::MappingQImplementation::compute_support_point_weights_cell<dim>(
239 ExcMessage(
"It only makes sense to create polynomial mappings "
240 "with a polynomial degree greater or equal to one."));
245template <
int dim,
int spacedim>
260template <
int dim,
int spacedim>
261std::unique_ptr<Mapping<dim, spacedim>>
264 return std::make_unique<MappingQ<dim, spacedim>>(*this);
269template <
int dim,
int spacedim>
278template <
int dim,
int spacedim>
319template <
int dim,
int spacedim>
338 const Point<1> &initial_p_unit)
const
370 const Point<2> &initial_p_unit)
const
400 const Point<3> &initial_p_unit)
const
430 const Point<1> &initial_p_unit)
const
433 const int spacedim = 2;
441 get_data(update_flags, point_quadrature));
463 const Point<2> &initial_p_unit)
const
466 const int spacedim = 3;
474 get_data(update_flags, point_quadrature));
504template <
int dim,
int spacedim>
513 ((dim == 1) || ((dim == 2) && (dim == spacedim))))
540 for (
unsigned int i = 0; i < vertices.size(); ++i)
541 vertices[i] = vertices_[i];
566 const double eps = 1e-15;
567 if (-eps <= point[1] && point[1] <= 1 + eps &&
568 -eps <= point[0] && point[0] <= 1 + eps)
605 return initial_p_unit;
610 for (
unsigned int d = 0; d < dim; ++d)
611 initial_p_unit[d] = 0.5;
618 AssertThrow(p_unit[0] != std::numeric_limits<double>::lowest(),
625template <
int dim,
int spacedim>
643 std::vector<Point<spacedim>> support_points_higher_order;
644 boost::container::small_vector<Point<spacedim>,
658 support_points_higher_order.data(),
666 const unsigned int n_points = real_points.size();
671 for (
unsigned int i = 0; i < n_points; i += n_lanes)
672 if (n_points - i > 1)
675 for (
unsigned int j = 0; j < n_lanes; ++j)
676 if (i + j < n_points)
677 for (
unsigned int d = 0; d < spacedim; ++d)
678 p_vec[d][j] = real_points[i + j][d];
680 for (
unsigned int d = 0; d < spacedim; ++d)
681 p_vec[d][j] = real_points[i][d];
687 inverse_approximation.compute(p_vec),
697 for (
unsigned int j = 0; j < n_lanes && i + j < n_points; ++j)
698 if (unit_point[0][j] != std::numeric_limits<double>::lowest())
699 for (
unsigned int d = 0; d < dim; ++d)
700 unit_points[i + j][d] = unit_point[d][j];
705 inverse_approximation.compute(real_points[i + j]),
714 inverse_approximation.compute(real_points[i]),
717 renumber_lexicographic_to_hierarchic);
722template <
int dim,
int spacedim>
732 UpdateFlags out = in;
733 for (
unsigned int i = 0; i < 5; ++i)
778template <
int dim,
int spacedim>
779std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
783 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
785 data_ptr->reinit(update_flags, q);
791template <
int dim,
int spacedim>
792std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
794 const UpdateFlags update_flags,
799 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
805 quadrature[0].size());
812template <
int dim,
int spacedim>
813std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
815 const UpdateFlags update_flags,
818 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
831template <
int dim,
int spacedim>
847 const unsigned int n_q_points = quadrature.
size();
863 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
884 computed_cell_similarity,
894 computed_cell_similarity,
905 computed_cell_similarity,
915 spacedim>(computed_cell_similarity,
924 spacedim>(computed_cell_similarity,
933 computed_cell_similarity,
942 spacedim>(computed_cell_similarity,
951 computed_cell_similarity,
959 const std::vector<double> &weights = quadrature.
get_weights();
967 (output_data.
JxW_values.size() == n_q_points),
977 for (
unsigned int point = 0; point < n_q_points; ++point)
992 cell->
center(), det, point)));
994 output_data.
JxW_values[point] = weights[point] * det;
1002 for (
unsigned int i = 0; i < spacedim; ++i)
1003 for (
unsigned int j = 0; j < dim; ++j)
1004 DX_t[j][i] = output_data.
jacobians[point][i][j];
1007 for (
unsigned int i = 0; i < dim; ++i)
1008 for (
unsigned int j = 0; j < dim; ++j)
1009 G[i][j] = DX_t[i] * DX_t[j];
1015 if (computed_cell_similarity ==
1026 Assert(spacedim == dim + 1,
1028 "There is no (unique) cell normal for " +
1030 "-dimensional cells in " +
1032 "-dimensional space. This only works if the "
1033 "space dimension is one greater than the "
1034 "dimensionality of the mesh cells."));
1054 return computed_cell_similarity;
1059template <
int dim,
int spacedim>
1063 const unsigned int face_no,
1091 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell;
1110 quadrature[0].
size()),
1120template <
int dim,
int spacedim>
1124 const unsigned int face_no,
1125 const unsigned int subface_no,
1151 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell;
1182template <
int dim,
int spacedim>
1199 const unsigned int n_q_points = quadrature.size();
1205 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
1275 const UpdateFlags update_flags = data.
update_each;
1276 const std::vector<double> &weights = quadrature.get_weights();
1288 for (
unsigned int point = 0; point < n_q_points; ++point)
1300 cell->
center(), det, point)));
1304 for (
unsigned int d = 0; d < spacedim; d++)
1308 output_data.
JxW_values[point] = weights[point] * det * normal.
norm();
1312 normal /= normal.
norm();
1321template <
int dim,
int spacedim>
1326 const UpdateFlags update_flags,
1338 output_data.
initialize(unit_points.size(), update_flags);
1340 auto internal_data =
1343 unit_points.end())));
1350 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
1369template <
int dim,
int spacedim>
1373 const unsigned int face_no,
1391 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
1413template <
int dim,
int spacedim>
1417 const MappingKind mapping_kind,
1429template <
int dim,
int spacedim>
1433 const MappingKind mapping_kind,
1445template <
int dim,
int spacedim>
1449 const MappingKind mapping_kind,
1453 switch (mapping_kind)
1477template <
int dim,
int spacedim>
1481 const MappingKind mapping_kind,
1489 &data = *
static_cast<const InternalData &
>(mapping_data).output_data;
1491 switch (mapping_kind)
1497 "update_covariant_transformation"));
1499 for (
unsigned int q = 0; q < output.size(); ++q)
1516template <
int dim,
int spacedim>
1520 const MappingKind mapping_kind,
1524 switch (mapping_kind)
1541template <
int dim,
int spacedim>
1550 for (
unsigned int line_no = 0;
1551 line_no < GeometryInfo<dim>::lines_per_cell;
1558 cell->
line(line_no));
1572 std::vector<Point<spacedim>> tmp_points;
1573 for (
unsigned int line_no = 0;
1574 line_no < GeometryInfo<dim>::lines_per_cell;
1581 cell->
line(line_no));
1590 const std::array<Point<spacedim>, 2> vertices{
1591 {cell->
vertex(reference_cell.line_to_cell_vertices(line_no, 0)),
1592 cell->
vertex(reference_cell.line_to_cell_vertices(line_no, 1))}};
1594 const std::size_t n_rows =
1596 a.resize(a.size() + n_rows);
1617 std::vector<Point<3>> tmp_points;
1620 for (
unsigned int face_no = 0; face_no <
faces_per_cell; ++face_no)
1635 Assert(face->vertex_index(i) ==
1637 face_no, i, face_orientation, face_flip, face_rotation)),
1645 face_no, i, face_orientation, face_flip, face_rotation)),
1651 boost::container::small_vector<Point<3>, 200> tmp_points(
1657 for (
unsigned int line = 0; line < GeometryInfo<2>::lines_per_cell;
1666 const std::size_t n_rows =
1668 a.resize(a.size() + n_rows);
1670 face->get_manifold().get_new_points(
1687 vertices[i] = cell->
vertex(i);
1700 const std::size_t n_rows = weights.size(0);
1701 a.resize(a.size() + n_rows);
1709template <
int dim,
int spacedim>
1720template <
int dim,
int spacedim>
1721std::vector<Point<spacedim>>
1726 std::vector<Point<spacedim>> a;
1729 a.push_back(cell->
vertex(i));
1738 bool all_manifold_ids_are_equal = (dim == spacedim);
1739 if (all_manifold_ids_are_equal &&
1745 all_manifold_ids_are_equal =
false;
1748 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
1750 all_manifold_ids_are_equal =
false;
1753 if (all_manifold_ids_are_equal)
1756 a.resize(a.size() + n_rows);
1775 if (dim != spacedim)
1779 const std::size_t n_rows =
1781 a.resize(a.size() + n_rows);
1797 const std::size_t n_rows =
1799 a.resize(a.size() + n_rows);
1819template <
int dim,
int spacedim>
1829template <
int dim,
int spacedim>
1834 Assert(dim == reference_cell.get_dimension(),
1835 ExcMessage(
"The dimension of your mapping (" +
1837 ") and the reference cell cell_type (" +
1839 " ) do not agree."));
1841 return reference_cell.is_hyper_cube();
1847#include "fe/mapping_q.inst"
auto make_const_array_view(const Container &container) -> decltype(make_array_view(container))
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
double diameter(const Mapping< dim, spacedim > &mapping) const
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
bool direction_flag() const
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
virtual void get_new_points(const ArrayView< const Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, ArrayView< Point< spacedim > > new_points) const
QGaussLobatto< 1 > line_support_points
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
virtual std::size_t memory_consumption() const override
std::vector< Point< spacedim > > mapping_support_points
const unsigned int polynomial_degree
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
bool tensor_product_quadrature
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > * output_data
std::vector< AlignedVector< Tensor< 1, spacedim > > > aux
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
InternalData(const unsigned int polynomial_degree)
AlignedVector< double > volume_elements
std::vector< Point< dim > > quadrature_points
const unsigned int n_shape_functions
const std::vector< unsigned int > renumber_lexicographic_to_hierarchic
const Table< 2, double > support_point_weights_cell
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void fill_mapping_data_for_face_quadrature(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_number, const Quadrature< dim - 1 > &face_quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const override
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual void fill_fe_immersed_surface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
const unsigned int polynomial_degree
virtual bool preserves_vertex_locations() const override
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
const std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
void fill_mapping_data_for_generic_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points, const UpdateFlags update_flags, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
const std::vector< Point< 1 > > line_support_points
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
const std::vector< Polynomials::Polynomial< double > > polynomials_1d
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
const std::vector< Point< dim > > unit_cell_support_points
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
MappingQ(const unsigned int polynomial_degree)
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
unsigned int get_degree() const
virtual std::size_t memory_consumption() const
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const
virtual boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
static DataSetDescriptor subface(const ReferenceCell &reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static DataSetDescriptor cell()
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
static Quadrature< dim > project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
bool is_tensor_product() const
const std::vector< double > & get_weights() const
const std::array< Quadrature< 1 >, dim > & get_tensor_basis() const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
numbers::NumberTraits< Number >::real_type norm() const
const Triangulation< dim, spacedim > & get_triangulation() const
Point< structdim > real_to_unit_cell_affine_approximation(const Point< spacedim > &point) const
const Manifold< dim, spacedim > & get_manifold() const
unsigned int vertex_index(const unsigned int i) const
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
Point< spacedim > & vertex(const unsigned int i) const
typename::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line(const unsigned int i) const
Triangulation< dim, spacedim > & get_triangulation()
static constexpr std::size_t size()
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
constexpr bool running_in_debug_mode()
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
constexpr unsigned int GeometryInfo< dim >::faces_per_cell
constexpr unsigned int GeometryInfo< dim >::vertices_per_face
constexpr unsigned int GeometryInfo< dim >::lines_per_face
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
static ::ExceptionBase & ExcTransformationFailed()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename IteratorSelector::line_iterator line_iterator
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face_iterator
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.
@ update_default
No update.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
@ mapping_covariant_gradient
@ mapping_contravariant_hessian
@ mapping_covariant_hessian
@ mapping_contravariant_gradient
bool face_rotation(const unsigned int face) const
bool face_orientation(const unsigned int face) const
types::geometric_orientation combined_face_orientation(const unsigned int face) const
bool face_flip(const unsigned int face) const
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
constexpr const ReferenceCell & get_hypercube()
constexpr unsigned int max_n_vertices()
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
constexpr T fixed_power(const T t)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
constexpr T pow(const T base, const int iexp)
Point< 1 > transform_real_to_unit_cell(const std::array< Point< spacedim >, GeometryInfo< 1 >::vertices_per_cell > &vertices, const Point< spacedim > &p)
void maybe_update_jacobian_pushed_forward_2nd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 4, spacedim > > &jacobian_pushed_forward_2nd_derivatives)
Point< dim, Number > do_transform_real_to_unit_cell_internal(const Point< spacedim, Number > &p, const Point< dim, Number > &initial_p_unit, const ArrayView< const Point< spacedim > > &points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber, const bool print_iterations_to_deallog=false)
void transform_differential_forms(const ArrayView< const DerivativeForm< rank, dim, spacedim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank+1, spacedim > > &output)
void maybe_update_jacobian_grads(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads)
void do_fill_fe_face_values(const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const typename QProjector< dim >::DataSetDescriptor data_set, const Quadrature< dim - 1 > &quadrature, const typename ::MappingQ< dim, spacedim >::InternalData &data, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
void transform_fields(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void maybe_update_jacobian_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 4, dim, spacedim > > &jacobian_3rd_derivatives)
void maybe_update_jacobian_pushed_forward_grads(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 3, spacedim > > &jacobian_pushed_forward_grads)
void maybe_update_q_points_Jacobians_generic(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Point< spacedim > > &quadrature_points, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians)
void transform_gradients(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void maybe_update_q_points_Jacobians_and_grads_tensor(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Point< spacedim > > &quadrature_points, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads)
void transform_hessians(const ArrayView< const Tensor< 3, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< 3, spacedim > > &output)
void maybe_update_jacobian_pushed_forward_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< Tensor< 5, spacedim > > &jacobian_pushed_forward_3rd_derivatives)
Point< dim > do_transform_real_to_unit_cell_internal_codim1(const Point< dim+1 > &p, const Point< dim > &initial_p_unit, const ArrayView< const Point< dim+1 > > &points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber)
void maybe_update_jacobian_2nd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::MappingQ< dim, spacedim >::InternalData &data, const ArrayView< const Point< dim > > &unit_points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, std::vector< DerivativeForm< 3, dim, spacedim > > &jacobian_2nd_derivatives)
Tensor< 3, spacedim, Number > apply_covariant_gradient(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 2, dim, spacedim, Number > &input)
ProductTypeNoPoint< Number, Number2 >::type evaluate_tensor_product_value_linear(const Number *values, const Point< dim, Number2 > &p)
ProductTypeNoPoint< Number, Number2 >::type evaluate_tensor_product_value(const std::vector< Polynomials::Polynomial< double > > &poly, const ArrayView< const Number > &values, const Point< dim, Number2 > &p, const bool d_linear=false, const std::vector< unsigned int > &renumber={})
constexpr unsigned int invalid_unsigned_int
constexpr types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static constexpr ndarray< Tensor< 1, dim >, faces_per_cell, dim - 1 > unit_tangential_vectors
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static constexpr unsigned int vertices_per_cell
static constexpr unsigned int lines_per_cell
static constexpr unsigned int faces_per_cell
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static constexpr unsigned int vertices_per_face
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
static unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static constexpr unsigned int lines_per_face
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)