Reference documentation for deal.II version 9.4.1
|
Classes | |
class | InverseQuadraticApproximation |
Functions | |
template<int dim> | |
std::vector< Point< dim > > | unit_support_points (const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering) |
inline ::Table< 2, double > | compute_support_point_weights_on_quad (const unsigned int polynomial_degree) |
inline ::Table< 2, double > | compute_support_point_weights_on_hex (const unsigned int polynomial_degree) |
std::vector<::Table< 2, double > > | compute_support_point_weights_perimeter_to_interior (const unsigned int polynomial_degree, const unsigned int dim) |
template<int dim> | |
inline ::Table< 2, double > | compute_support_point_weights_cell (const unsigned int polynomial_degree) |
template<int dim, int spacedim> | |
Point< spacedim > | compute_mapped_location_of_point (const typename ::MappingQ< dim, spacedim >::InternalData &data) |
template<int dim, int spacedim, typename Number > | |
Point< dim, Number > | do_transform_real_to_unit_cell_internal (const Point< spacedim, Number > &p, const Point< dim, Number > &initial_p_unit, const std::vector< Point< spacedim > > &points, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber, const bool print_iterations_to_deallog=false) |
template<int dim> | |
Point< dim > | do_transform_real_to_unit_cell_internal_codim1 (const typename ::Triangulation< dim, dim+1 >::cell_iterator &cell, const Point< dim+1 > &p, const Point< dim > &initial_p_unit, typename ::MappingQ< dim, dim+1 >::InternalData &mdata) |
template<int dim, int spacedim> | |
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< 2, dim, spacedim > > &jacobian_grads) |
template<int dim, int spacedim> | |
void | maybe_compute_q_points (const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Point< spacedim > > &quadrature_points) |
template<int dim, int spacedim> | |
void | maybe_update_Jacobians (const CellSimilarity::Similarity cell_similarity, const typename ::QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_grads (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_pushed_forward_grads (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Tensor< 3, spacedim > > &jacobian_pushed_forward_grads) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_2nd_derivatives (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< DerivativeForm< 3, dim, spacedim > > &jacobian_2nd_derivatives) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_pushed_forward_2nd_derivatives (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Tensor< 4, spacedim > > &jacobian_pushed_forward_2nd_derivatives) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_3rd_derivatives (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< DerivativeForm< 4, dim, spacedim > > &jacobian_3rd_derivatives) |
template<int dim, int spacedim> | |
void | maybe_update_jacobian_pushed_forward_3rd_derivatives (const CellSimilarity::Similarity cell_similarity, const typename QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingQ< dim, spacedim >::InternalData &data, std::vector< Tensor< 5, spacedim > > &jacobian_pushed_forward_3rd_derivatives) |
template<int dim, int spacedim> | |
void | maybe_compute_face_data (const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int n_q_points, const std::vector< double > &weights, const typename ::MappingQ< dim, spacedim >::InternalData &data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) |
template<int dim, int spacedim> | |
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, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) |
template<int dim, int spacedim, int rank> | |
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) |
template<int dim, int spacedim, int rank> | |
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) |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim, int rank> | |
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) |
Internal namespace to implement methods of MappingQ, such as the evaluation of the mapping and the transformation between real and unit cell.
std::vector< Point< dim > > internal::MappingQImplementation::unit_support_points | ( | const std::vector< Point< 1 > > & | line_support_points, |
const std::vector< unsigned int > & | renumbering | ||
) |
This function generates the reference cell support points from the 1d support points by expanding the tensor product.
Definition at line 217 of file mapping_q_internal.h.
inline ::Table< 2, double > internal::MappingQImplementation::compute_support_point_weights_on_quad | ( | const unsigned int | polynomial_degree | ) |
This function is needed by the constructor of MappingQ<dim,spacedim>
for dim=
2 and 3.
For the definition of the support_point_weights_on_quad
please refer to the description of TransfiniteInterpolationManifold.
Definition at line 247 of file mapping_q_internal.h.
inline ::Table< 2, double > internal::MappingQImplementation::compute_support_point_weights_on_hex | ( | const unsigned int | polynomial_degree | ) |
This function is needed by the constructor of MappingQ<3>
.
For the definition of the support_point_weights_on_quad
please refer to the description of TransfiniteInterpolationManifold.
Definition at line 299 of file mapping_q_internal.h.
|
inline |
This function collects the output of compute_support_point_weights_on_{quad,hex} in a single data structure.
Definition at line 388 of file mapping_q_internal.h.
inline ::Table< 2, double > internal::MappingQImplementation::compute_support_point_weights_cell | ( | const unsigned int | polynomial_degree | ) |
Collects all interior points for the various dimensions.
Definition at line 423 of file mapping_q_internal.h.
|
inline |
Using the relative weights of the shape functions evaluated at one point on the reference cell (and stored in data.shape_values and accessed via data.shape(0,i)) and the locations of mapping support points (stored in data.mapping_support_points), compute the mapped location of that point in real space.
Definition at line 455 of file mapping_q_internal.h.
|
inline |
Implementation of transform_real_to_unit_cell for either type double or VectorizedArray<double>
Definition at line 477 of file mapping_q_internal.h.
|
inline |
Implementation of transform_real_to_unit_cell for dim==spacedim-1
Definition at line 734 of file mapping_q_internal.h.
|
inline |
In case the quadrature formula is a tensor product, this is a replacement for maybe_compute_q_points(), maybe_update_Jacobians() and maybe_update_jacobian_grads()
Definition at line 1114 of file mapping_q_internal.h.
|
inline |
Compute the locations of quadrature points on the object described by the first argument (and the cell for which the mapping support points have already been set), but only if the update_flags of the data
argument indicate so.
Definition at line 1296 of file mapping_q_internal.h.
|
inline |
Update the co- and contravariant matrices as well as their determinant, for the cell described stored in the data object, but only if the update_flags of the data
argument indicate so.
Skip the computation if possible as indicated by the first argument.
Definition at line 1328 of file mapping_q_internal.h.
|
inline |
Update the Hessian of the transformation from unit to real cell, the Jacobian gradients.
Skip the computation if possible as indicated by the first argument.
Definition at line 1407 of file mapping_q_internal.h.
|
inline |
Update the Hessian of the transformation from unit to real cell, the Jacobian gradients, pushed forward to the real cell coordinates.
Skip the computation if possible as indicated by the first argument.
Definition at line 1454 of file mapping_q_internal.h.
|
inline |
Update the third derivatives of the transformation from unit to real cell, the Jacobian hessians.
Skip the computation if possible as indicated by the first argument.
Definition at line 1528 of file mapping_q_internal.h.
|
inline |
Update the Hessian of the Hessian of the transformation from unit to real cell, the Jacobian Hessian gradients, pushed forward to the real cell coordinates.
Skip the computation if possible as indicated by the first argument.
Definition at line 1584 of file mapping_q_internal.h.
|
inline |
Update the fourth derivatives of the transformation from unit to real cell, the Jacobian hessian gradients.
Skip the computation if possible as indicated by the first argument.
Definition at line 1685 of file mapping_q_internal.h.
|
inline |
Update the Hessian gradient of the transformation from unit to real cell, the Jacobian Hessians, pushed forward to the real cell coordinates.
Skip the computation if possible as indicated by the first argument.
Definition at line 1744 of file mapping_q_internal.h.
|
inline |
Depending on what information is called for in the update flags of the data
object, compute the various pieces of information that is required by the fill_fe_face_values() and fill_fe_subface_values() functions. This function simply unifies the work that would be done by those two functions.
The resulting data is put into the output_data
argument.
Definition at line 1871 of file mapping_q_internal.h.
|
inline |
Do the work of MappingQ::fill_fe_face_values() and MappingQ::fill_fe_subface_values() in a generic way, using the 'data_set' to differentiate whether we will work on a face (and if so, which one) or subface.
Definition at line 2035 of file mapping_q_internal.h.
|
inline |
Implementation of MappingQ::transform() for generic tensors.
Definition at line 2108 of file mapping_q_internal.h.
|
inline |
Implementation of MappingQ::transform() for gradients.
Definition at line 2186 of file mapping_q_internal.h.
|
inline |
Implementation of MappingQ::transform() for hessians.
Definition at line 2284 of file mapping_q_internal.h.
|
inline |
Implementation of MappingQ::transform() for DerivativeForm arguments.
Definition at line 2451 of file mapping_q_internal.h.