deal.II version 9.7.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\}}\)
Loading...
Searching...
No Matches
mapping_p1.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
19#include <deal.II/base/tensor.h>
20
24
26
27#include <algorithm>
28#include <cmath>
29#include <memory>
30
31
33
34namespace
35{
36 template <int dim, int spacedim>
38 compute_linear_transformation(
40 {
42 for (unsigned int j = 0; j < spacedim; ++j)
43 for (unsigned int i = 1; i < dim + 1; ++i)
44 result[j][i - 1] = cell->vertex(i)[j] - cell->vertex(0)[j];
45 return result;
46 }
47} // namespace
48
49
50
51template <int dim, int spacedim>
53 const ArrayView<const Point<dim>> &quadrature_points)
54{
55 quadrature.initialize(quadrature_points);
56}
57
58
59
60template <int dim, int spacedim>
65
66
67
68template <int dim, int spacedim>
69void
71 const UpdateFlags update_flags,
73{
74 this->quadrature = quadrature;
75 this->update_each = update_flags;
76}
77
78
79
80template <int dim, int spacedim>
81std::size_t
91
92
93
94template <int dim, int spacedim>
95bool
100
101
102
103template <int dim, int spacedim>
104bool
106 const ReferenceCell &reference_cell) const
107{
108 return reference_cell.is_simplex();
109}
110
111
112
113template <int dim, int spacedim>
116{
117 // Like MappingCartesian, this mapping is simple and has minimal update
118 // interdependencies
119 UpdateFlags out = in;
120 if (out & update_boundary_forms)
122
123 return out;
124}
125
126
127
128template <int dim, int spacedim>
129std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
130MappingP1<dim, spacedim>::get_data(const UpdateFlags update_flags,
131 const Quadrature<dim> &quadrature) const
132{
133 auto data_ptr = std::make_unique<InternalData>(quadrature);
134
135 // verify that we have computed the transitive hull of the required
136 // flags and that FEValues has faithfully passed them on to us
137 Assert(update_flags == requires_update_flags(update_flags),
139
140 // store the flags in the internal data object so we can access them
141 // in fill_fe_*_values(). use the transitive hull of the required
142 // flags
143 data_ptr->update_each = requires_update_flags(update_flags);
144
145 return data_ptr;
146}
147
148
149
150template <int dim, int spacedim>
151std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
153 const UpdateFlags update_flags,
154 const hp::QCollection<dim - 1> &quadrature) const
155{
156 auto data_ptr = std::make_unique<InternalData>(
158 quadrature));
159
160 // verify that we have computed the transitive hull of the required
161 // flags and that FEValues has faithfully passed them on to us
162 Assert(update_flags == requires_update_flags(update_flags),
164
165 // store the flags in the internal data object so we can access them
166 // in fill_fe_*_values()
167 data_ptr->update_each = update_flags;
168
169 return data_ptr;
170}
171
172
173
174template <int dim, int spacedim>
175std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
177 const UpdateFlags update_flags,
178 const Quadrature<dim - 1> &quadrature) const
179{
180 auto data_ptr = std::make_unique<InternalData>(
182 quadrature));
183
184 // verify that we have computed the transitive hull of the required
185 // flags and that FEValues has faithfully passed them on to us
186 Assert(update_flags == requires_update_flags(update_flags),
188
189 // store the flags in the internal data object so we can access them
190 // in fill_fe_*_values()
191 data_ptr->update_each = update_flags;
192
193 return data_ptr;
194}
195
196
197
198template <int dim, int spacedim>
199void
202 const InternalData &data) const
203{
204 data.affine_component = cell->vertex(0);
205 data.linear_component = compute_linear_transformation<dim, spacedim>(cell);
208}
209
210
211
212template <int dim, int spacedim>
213void
216 const InternalData &data,
217 const typename QProjector<dim>::DataSetDescriptor &offset,
218 std::vector<Point<spacedim>> &quadrature_points) const
219{
220 Assert(cell->vertex(0) == data.affine_component, ExcInternalError());
221 for (unsigned int i = 0; i < quadrature_points.size(); ++i)
222 quadrature_points[i] =
223 data.affine_component +
225 data.quadrature.point(offset + i));
226}
227
228
229
230template <int dim, int spacedim>
231void
233 const unsigned int face_no,
234 const InternalData &data,
235 std::vector<Tensor<1, spacedim>> &normal_vectors) const
236{
237 const Tensor<1, dim> ref_normal_vector =
238 ReferenceCells::get_simplex<dim>().template unit_normal_vectors<dim>(
239 face_no);
240 Tensor<1, spacedim> normal_vector =
241 apply_transformation(data.covariant, ref_normal_vector);
242 normal_vector /= normal_vector.norm();
243
244 std::fill(normal_vectors.begin(), normal_vectors.end(), normal_vector);
245}
246
247
248
249template <int dim, int spacedim>
250void
252 const InternalData &data,
253 const CellSimilarity::Similarity cell_similarity,
255 &output_data) const
256{
257 // The Jacobian is constant so its derivatives are zero
258 if (cell_similarity != CellSimilarity::translation)
259 {
261 std::fill(output_data.jacobian_grads.begin(),
262 output_data.jacobian_grads.end(),
264
266 std::fill(output_data.jacobian_pushed_forward_grads.begin(),
267 output_data.jacobian_pushed_forward_grads.end(),
269
271 std::fill(output_data.jacobian_2nd_derivatives.begin(),
272 output_data.jacobian_2nd_derivatives.end(),
274
276 std::fill(output_data.jacobian_pushed_forward_2nd_derivatives.begin(),
279
281 std::fill(output_data.jacobian_3rd_derivatives.begin(),
282 output_data.jacobian_3rd_derivatives.end(),
284
286 std::fill(output_data.jacobian_pushed_forward_3rd_derivatives.begin(),
289 }
290}
291
292
293
294template <int dim, int spacedim>
295void
297 const InternalData &data,
298 const CellSimilarity::Similarity cell_similarity,
300 &output_data) const
301{
302 // "compute" Jacobian at the quadrature points, which are all the
303 // same
304 if (data.update_each & update_jacobians)
305 if (cell_similarity != CellSimilarity::translation)
306 std::fill(output_data.jacobians.begin(),
307 output_data.jacobians.end(),
308 data.linear_component);
309}
310
311
312
313template <int dim, int spacedim>
314void
316 const InternalData &data,
317 const CellSimilarity::Similarity cell_similarity,
319 &output_data) const
320{
322 if (cell_similarity != CellSimilarity::translation)
323 {
324 const auto inverse = data.covariant.transpose();
325 std::fill(output_data.inverse_jacobians.begin(),
326 output_data.inverse_jacobians.end(),
327 inverse);
328 }
329}
330
331
332
333template <int dim, int spacedim>
337 const CellSimilarity::Similarity cell_similarity,
338 const Quadrature<dim> &quadrature,
339 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
341 &output_data) const
342{
343 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
345 const InternalData &data = static_cast<const InternalData &>(internal_data);
346
347 update_transformation(cell, data);
350 data,
352 output_data.quadrature_points);
353
354 // A special property of this mapping is that the Jacobian is constant over
355 // the cell
357 if (cell_similarity != CellSimilarity::translation)
358 {
359 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
360 output_data.JxW_values[i] = data.determinant * quadrature.weight(i);
361 }
362
363 maybe_update_jacobians(data, cell_similarity, output_data);
364 maybe_update_jacobian_derivatives(data, cell_similarity, output_data);
365 maybe_update_inverse_jacobians(data, cell_similarity, output_data);
366
368 {
369 Assert(spacedim == dim + 1,
370 ExcMessage("There is no (unique) cell normal for " +
372 "-dimensional cells in " +
373 Utilities::int_to_string(spacedim) +
374 "-dimensional space. This only works if the "
375 "space dimension is one greater than the "
376 "dimensionality of the mesh cells."));
377
378 Tensor<1, spacedim> normal;
379 // avoid warnings by only computing cross products in supported dimensions
380 if constexpr (dim == 1 && spacedim == 2)
381 normal = cross_product_2d(-data.linear_component.transpose()[0]);
382 else if constexpr (dim == 2 && spacedim == 3)
383 {
384 const auto transpose = data.linear_component.transpose();
385 normal = cross_product_3d(transpose[0], transpose[1]);
386 }
387 else
389 normal /= normal.norm();
390
391 if (cell->direction_flag() == false)
392 normal *= -1.0;
393
394 std::fill(output_data.normal_vectors.begin(),
395 output_data.normal_vectors.end(),
396 normal);
397 }
398
399 return cell_similarity;
400}
401
402
403
404template <int dim, int spacedim>
405void
408 const ArrayView<const Point<dim>> &unit_points,
409 const UpdateFlags update_flags,
411 &output_data) const
412{
413 if (update_flags == update_default)
414 return;
415
416 Assert(update_flags & update_inverse_jacobians ||
417 update_flags & update_jacobians ||
418 update_flags & update_quadrature_points,
420
421 output_data.initialize(unit_points.size(), update_flags);
422
423 InternalData data(unit_points);
424 data.update_each = update_flags;
425 update_transformation(cell, data);
426
429 data,
431 output_data.quadrature_points);
432
435}
436
437
438
439template <int dim, int spacedim>
440void
443 const unsigned int face_no,
444 const hp::QCollection<dim - 1> &quadrature_collection,
445 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
447 &output_data) const
448{
449 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
451 const InternalData &data = static_cast<const InternalData &>(internal_data);
452
453 update_transformation(cell, data);
454 const auto offset =
456 face_no,
458 face_no),
459 quadrature_collection);
462 data,
463 offset,
464 output_data.quadrature_points);
465
466 maybe_update_normal_vectors(face_no, data, output_data.normal_vectors);
467
469 {
470 // Since the quadrature weights presently sum to
471 // cell->reference_cell().face_measure(face_no), we have to rescale so
472 // they sum to the area of the face
473 const double J = cell->face(face_no)->measure() /
474 cell->reference_cell().face_measure(face_no);
476 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
477 output_data.JxW_values[i] = J * data.quadrature.weight(i + offset);
478
480 for (unsigned int i = 0; i < output_data.boundary_forms.size(); ++i)
481 output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
482 }
483
487}
488
489
490
491template <int dim, int spacedim>
492void
495 const unsigned int face_no,
496 const unsigned int subface_no,
497 const Quadrature<dim - 1> &quadrature,
498 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
500 &output_data) const
501{
502 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
504 const InternalData &data = static_cast<const InternalData &>(internal_data);
505
506 update_transformation(cell, data);
507 const auto offset =
509 face_no,
510 subface_no,
512 face_no),
513 quadrature.size(),
514 cell->subface_case(face_no));
517 data,
518 offset,
519 output_data.quadrature_points);
520
521 maybe_update_normal_vectors(face_no, data, output_data.normal_vectors);
522
524 {
525 // Same as fill_fe_face_values()
526 const double J =
527 cell->face(face_no)->measure() /
528 cell->face(face_no)->reference_cell().volume() /
529 // TODO: once we support 3d refinement this should be updated to the
530 // simplex version of GeometryInfo::subface_ratio()
531 (dim == 2 ? 2 : 4);
533 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
534 output_data.JxW_values[i] = J * quadrature.weight(i);
535
537 for (unsigned int i = 0; i < output_data.boundary_forms.size(); ++i)
538 output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
539 }
540
544}
545
546
547
548template <int dim, int spacedim>
549void
551 const ArrayView<const Tensor<1, dim>> &input,
552 const MappingKind mapping_kind,
553 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
554 const ArrayView<Tensor<1, spacedim>> &output) const
555{
556 AssertDimension(input.size(), output.size());
557 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
559 const InternalData &data = static_cast<const InternalData &>(mapping_data);
560
561 switch (mapping_kind)
562 {
564 {
565 for (unsigned int i = 0; i < output.size(); ++i)
566 output[i] = apply_transformation(data.covariant, input[i]);
567 return;
568 }
569
571 {
572 for (unsigned int i = 0; i < output.size(); ++i)
573 output[i] = apply_transformation(data.linear_component, input[i]);
574 return;
575 }
576 case mapping_piola:
577 {
578 auto transformation = data.linear_component;
579 Assert(data.determinant > 0.0, ExcDivideByZero());
580 for (unsigned int d = 0; d < spacedim; ++d)
581 transformation[d] *= 1.0 / data.determinant;
582 for (unsigned int i = 0; i < output.size(); ++i)
583 output[i] = apply_transformation(transformation, input[i]);
584 return;
585 }
586 default:
588 }
589}
590
591
592
593template <int dim, int spacedim>
594void
597 const MappingKind mapping_kind,
598 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
599 const ArrayView<Tensor<2, spacedim>> &output) const
600{
601 AssertDimension(input.size(), output.size());
602 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
604 const InternalData &data = static_cast<const InternalData &>(mapping_data);
605
606 switch (mapping_kind)
607 {
609 {
610 for (unsigned int i = 0; i < output.size(); ++i)
611 output[i] = apply_transformation(data.covariant, input[i]);
612
613 return;
614 }
615 default:
617 }
618}
619
620
621
622template <int dim, int spacedim>
623void
625 const ArrayView<const Tensor<2, dim>> &input,
626 const MappingKind mapping_kind,
627 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
628 const ArrayView<Tensor<2, spacedim>> &output) const
629{
630 AssertDimension(input.size(), output.size());
631 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
633 const InternalData &data = static_cast<const InternalData &>(mapping_data);
634
635 switch (mapping_kind)
636 {
638 {
641 "update_covariant_transformation"));
642
643 for (unsigned int i = 0; i < output.size(); ++i)
644 output[i] = apply_transformation(data.covariant, input[i]);
645 return;
646 }
647
649 {
652 "update_contravariant_transformation"));
653
654 for (unsigned int i = 0; i < output.size(); ++i)
655 output[i] = apply_transformation(data.linear_component, input[i]);
656 return;
657 }
658
660 {
663 "update_covariant_transformation"));
664
665 for (unsigned int i = 0; i < output.size(); ++i)
666 {
669 output[i] = apply_transformation(data.covariant, A.transpose());
670 }
671 return;
672 }
673
675 {
678 "update_covariant_transformation"));
679
680 for (unsigned int i = 0; i < output.size(); ++i)
681 {
684 transpose(input[i]));
685 output[i] = apply_transformation(data.covariant, A.transpose());
686 }
687
688 return;
689 }
690
692 {
695 "update_contravariant_transformation"));
696
697 auto scaled_contravariant = data.linear_component;
698 for (unsigned int d = 0; d < spacedim; ++d)
699 scaled_contravariant[d] *= 1.0 / data.determinant;
700 for (unsigned int i = 0; i < output.size(); ++i)
701 {
703 apply_transformation(data.covariant, input[i]);
704 const Tensor<2, spacedim> T =
705 apply_transformation(scaled_contravariant, A.transpose());
706
707 output[i] = transpose(T);
708 }
709
710 return;
711 }
712
713 default:
714 Assert(false, ExcNotImplemented());
715 }
716}
717
718
719
720template <int dim, int spacedim>
721void
724 const MappingKind mapping_kind,
725 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
726 const ArrayView<Tensor<3, spacedim>> &output) const
727{
728 AssertDimension(input.size(), output.size());
729 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
731 const InternalData &data = static_cast<const InternalData &>(mapping_data);
732
733 switch (mapping_kind)
734 {
736 {
739 "update_covariant_transformation"));
740
741 for (unsigned int i = 0; i < output.size(); ++i)
742 output[i] =
744
745 return;
746 }
747 default:
749 }
750}
751
752
753
754template <int dim, int spacedim>
755void
757 const ArrayView<const Tensor<3, dim>> &input,
758 const MappingKind mapping_kind,
759 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
760 const ArrayView<Tensor<3, spacedim>> &output) const
761{
762 AssertDimension(input.size(), output.size());
763 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
765 const InternalData &data = static_cast<const InternalData &>(mapping_data);
766
767 switch (mapping_kind)
768 {
770 {
773 "update_covariant_transformation"));
776 "update_contravariant_transformation"));
777
778 for (unsigned int i = 0; i < output.size(); ++i)
779 output[i] =
781 data.linear_component,
782 input[i]);
783
784 return;
785 }
786
788 {
791 "update_covariant_transformation"));
792
793 for (unsigned int i = 0; i < output.size(); ++i)
794 output[i] =
796
797 return;
798 }
799
801 {
804 "update_covariant_transformation"));
807 "update_contravariant_transformation"));
808
809 for (unsigned int i = 0; i < output.size(); ++i)
811 data.linear_component,
812 data.determinant,
813 input[i]);
814
815 return;
816 }
817
818 default:
820 }
821}
822
823
824
825template <int dim, int spacedim>
829 const Point<dim> &p) const
830{
831 const DerivativeForm<1, dim, spacedim> linear_component =
832 compute_linear_transformation<dim, spacedim>(cell);
833 const Tensor<1, spacedim> sheared = apply_transformation(linear_component, p);
834 return cell->vertex(0) + sheared;
835}
836
837
838
839template <int dim, int spacedim>
843 const Point<spacedim> &p) const
844{
845 const DerivativeForm<1, spacedim, dim> linear_component =
846 compute_linear_transformation<dim, spacedim>(cell)
848 .transpose();
849 const Tensor<1, spacedim> offset = cell->vertex(0);
850 return Point<dim>(apply_transformation(linear_component, p - offset));
851}
852
853
854
855template <int dim, int spacedim>
856void
859 const ArrayView<const Point<spacedim>> &real_points,
860 const ArrayView<Point<dim>> &unit_points) const
861{
862 const DerivativeForm<1, spacedim, dim> linear_component =
863 compute_linear_transformation<dim, spacedim>(cell)
865 .transpose();
866 const Tensor<1, spacedim> offset = cell->vertex(0);
867 for (unsigned int i = 0; i < real_points.size(); ++i)
868 unit_points[i] = Point<dim>(
869 apply_transformation(linear_component, real_points[i] - offset));
870}
871
872
873
874template <int dim, int spacedim>
875std::unique_ptr<Mapping<dim, spacedim>>
877{
878 return std::make_unique<MappingP1<dim, spacedim>>(*this);
879}
880
881
882//---------------------------------------------------------------------------
883// explicit instantiations
884#include "fe/mapping_p1.inst"
885
886
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
bool direction_flag() const
Number determinant() const
DerivativeForm< 1, dim, spacedim, Number > covariant_form() const
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
DerivativeForm< 1, dim, spacedim > covariant
Definition mapping_p1.h:231
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Definition mapping_p1.cc:70
DerivativeForm< 1, dim, spacedim > linear_component
Definition mapping_p1.h:226
Tensor< 1, spacedim > affine_component
Definition mapping_p1.h:221
Quadrature< dim > quadrature
Definition mapping_p1.h:242
InternalData(const ArrayView< const Point< dim > > &quadrature_points)
Definition mapping_p1.cc:52
virtual std::size_t memory_consumption() const override
Definition mapping_p1.cc:82
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) 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 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
void transform_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const typename QProjector< dim >::DataSetDescriptor &offset, std::vector< Point< spacedim > > &quadrature_points) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual bool preserves_vertex_locations() const override
Definition mapping_p1.cc:96
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 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 Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
void update_transformation(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, spacedim > > &normal_vectors) const
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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
UpdateFlags update_each
Definition mapping.h:704
virtual std::size_t memory_consumption() const
Definition point.h:113
Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to...
Definition qprojector.h:334
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()
Definition qprojector.h:516
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)
double weight(const unsigned int i) const
unsigned int size() const
double face_measure(const unsigned int face_no) const
numbers::NumberTraits< Number >::real_type norm() const
Point< spacedim > & vertex(const unsigned int i) const
ReferenceCell reference_cell() const
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
void initialize(const unsigned int n_quadrature_points, const UpdateFlags flags)
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
Tensor< 1, spacedim, typename ProductType< Number1, Number2 >::type > apply_transformation(const DerivativeForm< 1, dim, spacedim, Number1 > &grad_F, const Tensor< 1, dim, Number2 > &d_x)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
UpdateFlags
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_jacobian_3rd_derivatives
@ 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.
@ update_jacobian_2nd_derivatives
@ mapping_piola
Definition mapping.h:116
@ mapping_covariant_gradient
Definition mapping.h:102
@ mapping_covariant
Definition mapping.h:91
@ mapping_contravariant
Definition mapping.h:96
@ mapping_contravariant_hessian
Definition mapping.h:158
@ mapping_covariant_hessian
Definition mapping.h:152
@ mapping_contravariant_gradient
Definition mapping.h:108
@ mapping_piola_gradient
Definition mapping.h:122
@ mapping_piola_hessian
Definition mapping.h:164
types::geometric_orientation combined_face_orientation(const unsigned int face) const
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
constexpr const ReferenceCell & get_simplex()
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:466
Tensor< 3, spacedim, Number > apply_contravariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_piola_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Number &volume_element, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_gradient(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 2, dim, spacedim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const Tensor< 3, dim, Number > &input)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition tensor.h:2647
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)
Definition tensor.h:2672