|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
42 template <
int dim,
int spacedim>
52 template <
int dim,
int spacedim>
64 template <
int dim,
int spacedim>
73 template <
int dim,
int spacedim>
92 template <
int dim,
int spacedim>
93 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
97 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
98 std_cxx14::make_unique<InternalData>(q);
111 template <
int dim,
int spacedim>
112 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
117 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
118 std_cxx14::make_unique<InternalData>(
129 data.update_each = update_flags;
136 template <
int dim,
int spacedim>
137 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
142 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
143 std_cxx14::make_unique<InternalData>(
154 data.update_each = update_flags;
161 template <
int dim,
int spacedim>
197 template <
int dim,
int spacedim>
215 template <
int dim,
int spacedim>
219 const unsigned int face_no,
229 cell->face_orientation(
231 cell->face_flip(face_no),
232 cell->face_rotation(face_no),
242 template <
int dim,
int spacedim>
246 const unsigned int face_no,
247 const unsigned int sub_no,
253 if (cell->face(face_no)->has_children())
264 cell->face_orientation(face_no),
265 cell->face_flip(face_no),
266 cell->face_rotation(face_no),
268 cell->subface_case(face_no));
276 template <
int dim,
int spacedim>
291 for (
unsigned int d = 0;
d < dim; ++
d)
299 template <
int dim,
int spacedim>
302 const unsigned int face_no,
310 std::fill(normal_vectors.begin(),
311 normal_vectors.end(),
318 template <
int dim,
int spacedim>
329 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
333 for (
unsigned int i = 0;
339 for (
unsigned int i = 0;
346 for (
unsigned int i = 0;
353 for (
unsigned int i = 0;
360 for (
unsigned int i = 0;
370 template <
int dim,
int spacedim>
399 for (
unsigned int d = 1;
d < dim; ++
d)
403 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
410 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
413 for (
unsigned int j = 0; j < dim; ++j)
426 for (
unsigned int j = 0; j < dim; ++j)
430 return cell_similarity;
435 template <
int dim,
int spacedim>
439 const unsigned int face_no,
465 for (
unsigned int d = 0;
d < dim; ++
d)
470 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
474 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
480 for (
unsigned int d = 1;
d < dim; ++
d)
486 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
489 for (
unsigned int d = 0;
d < dim; ++
d)
499 for (
unsigned int d = 0;
d < dim; ++
d)
506 template <
int dim,
int spacedim>
510 const unsigned int face_no,
511 const unsigned int subface_no,
533 for (
unsigned int d = 0;
d < dim; ++
d)
543 const unsigned int n_subfaces =
544 cell->face(face_no)->has_children() ?
545 cell->face(face_no)->n_children() :
547 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
552 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
558 for (
unsigned int d = 1;
d < dim; ++
d)
564 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
567 for (
unsigned int d = 0;
d < dim; ++
d)
577 for (
unsigned int d = 0;
d < dim; ++
d)
584 template <
int dim,
int spacedim>
597 switch (mapping_kind)
603 "update_covariant_transformation"));
605 for (
unsigned int i = 0; i < output.size(); ++i)
606 for (
unsigned int d = 0;
d < dim; ++
d)
615 "update_contravariant_transformation"));
617 for (
unsigned int i = 0; i < output.size(); ++i)
618 for (
unsigned int d = 0;
d < dim; ++
d)
626 "update_contravariant_transformation"));
629 "update_volume_elements"));
631 for (
unsigned int i = 0; i < output.size(); ++i)
632 for (
unsigned int d = 0;
d < dim; ++
d)
644 template <
int dim,
int spacedim>
657 switch (mapping_kind)
663 "update_covariant_transformation"));
665 for (
unsigned int i = 0; i < output.size(); ++i)
666 for (
unsigned int d1 = 0; d1 < dim; ++d1)
667 for (
unsigned int d2 = 0; d2 < dim; ++d2)
668 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
676 "update_contravariant_transformation"));
678 for (
unsigned int i = 0; i < output.size(); ++i)
679 for (
unsigned int d1 = 0; d1 < dim; ++d1)
680 for (
unsigned int d2 = 0; d2 < dim; ++d2)
681 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
689 "update_covariant_transformation"));
691 for (
unsigned int i = 0; i < output.size(); ++i)
692 for (
unsigned int d1 = 0; d1 < dim; ++d1)
693 for (
unsigned int d2 = 0; d2 < dim; ++d2)
694 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
703 "update_contravariant_transformation"));
705 for (
unsigned int i = 0; i < output.size(); ++i)
706 for (
unsigned int d1 = 0; d1 < dim; ++d1)
707 for (
unsigned int d2 = 0; d2 < dim; ++d2)
708 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
717 "update_contravariant_transformation"));
720 "update_volume_elements"));
722 for (
unsigned int i = 0; i < output.size(); ++i)
723 for (
unsigned int d1 = 0; d1 < dim; ++d1)
724 for (
unsigned int d2 = 0; d2 < dim; ++d2)
725 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
734 "update_contravariant_transformation"));
737 "update_volume_elements"));
739 for (
unsigned int i = 0; i < output.size(); ++i)
740 for (
unsigned int d1 = 0; d1 < dim; ++d1)
741 for (
unsigned int d2 = 0; d2 < dim; ++d2)
742 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
754 template <
int dim,
int spacedim>
767 switch (mapping_kind)
773 "update_covariant_transformation"));
775 for (
unsigned int i = 0; i < output.size(); ++i)
776 for (
unsigned int d1 = 0; d1 < dim; ++d1)
777 for (
unsigned int d2 = 0; d2 < dim; ++d2)
778 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
786 "update_contravariant_transformation"));
788 for (
unsigned int i = 0; i < output.size(); ++i)
789 for (
unsigned int d1 = 0; d1 < dim; ++d1)
790 for (
unsigned int d2 = 0; d2 < dim; ++d2)
791 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
799 "update_covariant_transformation"));
801 for (
unsigned int i = 0; i < output.size(); ++i)
802 for (
unsigned int d1 = 0; d1 < dim; ++d1)
803 for (
unsigned int d2 = 0; d2 < dim; ++d2)
804 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
813 "update_contravariant_transformation"));
815 for (
unsigned int i = 0; i < output.size(); ++i)
816 for (
unsigned int d1 = 0; d1 < dim; ++d1)
817 for (
unsigned int d2 = 0; d2 < dim; ++d2)
818 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
827 "update_contravariant_transformation"));
830 "update_volume_elements"));
832 for (
unsigned int i = 0; i < output.size(); ++i)
833 for (
unsigned int d1 = 0; d1 < dim; ++d1)
834 for (
unsigned int d2 = 0; d2 < dim; ++d2)
835 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
844 "update_contravariant_transformation"));
847 "update_volume_elements"));
849 for (
unsigned int i = 0; i < output.size(); ++i)
850 for (
unsigned int d1 = 0; d1 < dim; ++d1)
851 for (
unsigned int d2 = 0; d2 < dim; ++d2)
852 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
864 template <
int dim,
int spacedim>
877 switch (mapping_kind)
883 "update_covariant_transformation"));
885 for (
unsigned int q = 0; q < output.size(); ++q)
886 for (
unsigned int i = 0; i < spacedim; ++i)
887 for (
unsigned int j = 0; j < spacedim; ++j)
888 for (
unsigned int k = 0; k < spacedim; ++k)
890 output[q][i][j][k] = input[q][i][j][k] /
903 template <
int dim,
int spacedim>
916 switch (mapping_kind)
922 "update_covariant_transformation"));
925 "update_contravariant_transformation"));
927 for (
unsigned int q = 0; q < output.size(); ++q)
928 for (
unsigned int i = 0; i < spacedim; ++i)
929 for (
unsigned int j = 0; j < spacedim; ++j)
930 for (
unsigned int k = 0; k < spacedim; ++k)
943 "update_covariant_transformation"));
945 for (
unsigned int q = 0; q < output.size(); ++q)
946 for (
unsigned int i = 0; i < spacedim; ++i)
947 for (
unsigned int j = 0; j < spacedim; ++j)
948 for (
unsigned int k = 0; k < spacedim; ++k)
962 "update_covariant_transformation"));
965 "update_contravariant_transformation"));
968 "update_volume_elements"));
970 for (
unsigned int q = 0; q < output.size(); ++q)
971 for (
unsigned int i = 0; i < spacedim; ++i)
972 for (
unsigned int j = 0; j < spacedim; ++j)
973 for (
unsigned int k = 0; k < spacedim; ++k)
991 template <
int dim,
int spacedim>
1002 length[0] = cell->vertex(1)(0) - start(0);
1005 length[0] = cell->vertex(1)(0) - start(0);
1006 length[1] = cell->vertex(2)(1) - start(1);
1009 length[0] = cell->vertex(1)(0) - start(0);
1010 length[1] = cell->vertex(2)(1) - start(1);
1011 length[2] = cell->vertex(4)(2) - start(2);
1018 for (
unsigned int d = 0;
d < dim; ++
d)
1019 p_real(
d) += length[
d] * p(
d);
1026 template <
int dim,
int spacedim>
1032 if (dim != spacedim)
1041 real(0) /= cell->vertex(1)(0) - start(0);
1044 real(0) /= cell->vertex(1)(0) - start(0);
1045 real(1) /= cell->vertex(2)(1) - start(1);
1048 real(0) /= cell->vertex(1)(0) - start(0);
1049 real(1) /= cell->vertex(2)(1) - start(1);
1050 real(2) /= cell->vertex(4)(2) - start(2);
1060 template <
int dim,
int spacedim>
1061 std::unique_ptr<Mapping<dim, spacedim>>
1064 return std_cxx14::make_unique<MappingCartesian<dim, spacedim>>(*this);
1070 #include "mapping_cartesian.inst"
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
@ update_quadrature_points
Transformed quadrature points.
@ update_jacobian_pushed_forward_3rd_derivatives
std::vector< Point< dim > > quadrature_points
void maybe_update_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, std::vector< Point< dim >> &quadrature_points) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
@ mapping_covariant_gradient
InternalData(const Quadrature< dim > &quadrature)
static ::ExceptionBase & ExcNotImplemented()
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
@ mapping_contravariant_gradient
#define AssertIndexRange(index, range)
double weight(const unsigned int i) const
void maybe_update_subface_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const InternalData &data, std::vector< Point< dim >> &quadrature_points) const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
virtual bool preserves_vertex_locations() const override
@ update_jacobian_pushed_forward_grads
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim >> &normal_vectors) const
Abstract base class for mapping classes.
@ mapping_covariant_hessian
Tensor< 1, dim > cell_extents
@ update_covariant_transformation
Covariant transformation.
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) 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
#define DEAL_II_NAMESPACE_OPEN
static DataSetDescriptor subface(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)
@ update_jacobians
Volume element.
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
@ update_jacobian_pushed_forward_2nd_derivatives
#define AssertDimension(dim1, dim2)
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
@ update_jacobian_3rd_derivatives
@ update_jacobian_grads
Gradient of volume element.
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
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
static ::ExceptionBase & ExcInternalError()
@ update_JxW_values
Transformed quadrature weights.
@ update_normal_vectors
Normal vectors.
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim >>> &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim, spacedim >::mapping)
#define Assert(cond, exc)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void transform_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const typename QProjector< dim >::DataSetDescriptor &offset, std::vector< Point< dim >> &quadrature_points) const
static DataSetDescriptor cell()
@ update_inverse_jacobians
Volume element.
void maybe_update_face_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const InternalData &data, std::vector< Point< dim >> &quadrature_points) const
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
#define DEAL_II_NAMESPACE_CLOSE
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
virtual std::size_t memory_consumption() const override
@ update_boundary_forms
Outer normal vector, not normalized.
@ mapping_contravariant_hessian
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
@ update_jacobian_2nd_derivatives