PolyDEAL
 
Loading...
Searching...
No Matches
mapping_box.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2024 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
20#include <deal.II/base/tensor.h>
21
23
25
26#include <deal.II/grid/tria.h>
28
30
31#include <mapping_box.h>
32
33#include <algorithm>
34
35
36
38
40 ExcCellNotAssociatedWithBox,
41 "You are using MappingBox, but the incoming element is not associated with a"
42 "Bounding Box Cartesian.");
43
44
45
50template <typename CellType>
51bool
52has_box(const CellType &cell,
53 const std::map<types::global_cell_index, types::global_cell_index>
54 &translator)
55{
56 Assert((cell->reference_cell().is_hyper_cube() ||
57 cell->reference_cell().is_simplex()),
58 ExcNotImplemented());
59 Assert((translator.find(cell->active_cell_index()) != translator.cend()),
60 ExcCellNotAssociatedWithBox());
61
62 return true;
63}
64
65
66
67template <int dim, int spacedim>
69 const std::vector<BoundingBox<dim>> &input_boxes,
70 const std::map<types::global_cell_index, types::global_cell_index>
71 &global_to_polytope)
72{
73 Assert(input_boxes.size() > 0,
74 ExcMessage("Invalid number of bounding boxes."));
75
76 // copy boxes and map
77 boxes.resize(input_boxes.size());
78 for (unsigned int i = 0; i < input_boxes.size(); ++i)
79 boxes[i] = input_boxes[i];
80 polytope_translator = global_to_polytope;
81}
82
83
84
85template <int dim, int spacedim>
87 : cell_extents(numbers::signaling_nan<Tensor<1, dim>>())
88 , traslation(numbers::signaling_nan<Tensor<1, dim>>())
89 , inverse_cell_extents(numbers::signaling_nan<Tensor<1, dim>>())
90 , volume_element(numbers::signaling_nan<double>())
91 , quadrature_points(q.get_points())
92{}
93
94
95
96template <int dim, int spacedim>
97void
99 const Quadrature<dim> &)
100{
101 // store the flags in the internal data object so we can access them
102 // in fill_fe_*_values(). use the transitive hull of the required
103 // flags
104 this->update_each = update_flags;
105}
106
107
108
109template <int dim, int spacedim>
110std::size_t
119
120
121
122template <int dim, int spacedim>
123bool
128
129
130
131template <int dim, int spacedim>
132bool
134 const ReferenceCell &reference_cell) const
135{
136 Assert(dim == reference_cell.get_dimension(),
137 ExcMessage("The dimension of your mapping (" +
139 ") and the reference cell cell_type (" +
140 Utilities::to_string(reference_cell.get_dimension()) +
141 " ) do not agree."));
142
143 return reference_cell.is_hyper_cube() || reference_cell.is_simplex();
144}
145
146
147
148template <int dim, int spacedim>
151{
152 // this mapping is pretty simple in that it can basically compute
153 // every piece of information wanted by FEValues without requiring
154 // computing any other quantities. boundary forms are one exception
155 // since they can be computed from the normal vectors without much
156 // further ado
157 UpdateFlags out = in;
158 if (out & update_boundary_forms)
160
161 return out;
162}
163
164
165
166template <int dim, int spacedim>
167std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
169 const Quadrature<dim> &q) const
170{
171 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
172 std::make_unique<InternalData>();
173 data_ptr->reinit(requires_update_flags(update_flags), q);
174
175 return data_ptr;
176}
177
178
179
180template <int dim, int spacedim>
181std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
183 const UpdateFlags update_flags,
184 const Quadrature<dim - 1> &quadrature) const
185{
186 (void)update_flags;
187 (void)quadrature;
189 return {};
190}
191
192
193
194template <int dim, int spacedim>
195void
198 const CellSimilarity::Similarity cell_similarity,
199 const InternalData &data) const
200{
201 // Compute start point and sizes along axes. The vertices to be looked at
202 // are 1, 2, 4 compared to the base vertex 0.
203 if (cell_similarity != CellSimilarity::translation)
204 {
205 const BoundingBox<dim> &current_box =
206 boxes[polytope_translator.at(cell->active_cell_index())];
207 const std::pair<Point<dim>, Point<dim>> &bdary_points =
208 current_box.get_boundary_points();
209
210 for (unsigned int d = 0; d < dim; ++d)
211 {
212 const double cell_extent_d = current_box.side_length(d);
213 data.cell_extents[d] = cell_extent_d;
214
215 data.traslation[d] =
216 .5 * (bdary_points.first[d] +
217 bdary_points.second[d]); // midpoint of each interval
218
219 Assert(cell_extent_d != 0.,
220 ExcMessage("Cell does not appear to be Cartesian!"));
221 data.inverse_cell_extents[d] = 1. / cell_extent_d;
222 }
223 }
224}
225
226
227
228namespace
229{
230 template <int dim>
231 void
232 transform_quadrature_points(
233 const BoundingBox<dim> &box,
234 const ArrayView<const Point<dim>> &unit_quadrature_points,
235 std::vector<Point<dim>> &quadrature_points)
236 {
237 for (unsigned int i = 0; i < quadrature_points.size(); ++i)
238 quadrature_points[i] = box.unit_to_real(unit_quadrature_points[i]);
239 }
240} // namespace
241
242
243
244template <int dim, int spacedim>
245void
248 const InternalData &data,
249 const ArrayView<const Point<dim>> &unit_quadrature_points,
250 std::vector<Point<dim>> &quadrature_points) const
251{
253 transform_quadrature_points(
254 boxes[polytope_translator.at(cell->active_cell_index())],
255 unit_quadrature_points,
256 quadrature_points);
257}
258
259
260
261template <int dim, int spacedim>
262void
264 const unsigned int face_no,
265 const InternalData &data,
266 std::vector<Tensor<1, dim>> &normal_vectors) const
267{
268 // compute normal vectors. All normals on a face have the same value.
270 {
271 Assert(face_no < GeometryInfo<dim>::faces_per_cell, ExcInternalError());
272 std::fill(normal_vectors.begin(),
273 normal_vectors.end(),
275 }
276}
277
278
279
280template <int dim, int spacedim>
281void
283 const InternalData &data,
284 const CellSimilarity::Similarity cell_similarity,
286 &output_data) const
287{
288 if (cell_similarity != CellSimilarity::translation)
289 {
291 for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
293
295 for (unsigned int i = 0;
296 i < output_data.jacobian_pushed_forward_grads.size();
297 ++i)
299
301 for (unsigned int i = 0;
302 i < output_data.jacobian_2nd_derivatives.size();
303 ++i)
304 output_data.jacobian_2nd_derivatives[i] =
306
308 for (unsigned int i = 0;
309 i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
310 ++i)
313
315 for (unsigned int i = 0;
316 i < output_data.jacobian_3rd_derivatives.size();
317 ++i)
318 output_data.jacobian_3rd_derivatives[i] =
320
322 for (unsigned int i = 0;
323 i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
324 ++i)
327 }
328}
329
330
331
332template <int dim, int spacedim>
333void
335 const InternalData &data) const
336{
338 {
339 double volume = data.cell_extents[0];
340 for (unsigned int d = 1; d < dim; ++d)
341 volume *= data.cell_extents[d];
342 data.volume_element = volume;
343 }
344}
345
346
347
348template <int dim, int spacedim>
349void
351 const InternalData &data,
352 const CellSimilarity::Similarity cell_similarity,
354 &output_data) const
355{
356 // "compute" Jacobian at the quadrature points, which are all the
357 // same
358 if (data.update_each & update_jacobians)
359 if (cell_similarity != CellSimilarity::translation)
360 for (unsigned int i = 0; i < output_data.jacobians.size(); ++i)
361 {
363 for (unsigned int j = 0; j < dim; ++j)
364 output_data.jacobians[i][j][j] = data.cell_extents[j];
365 }
366}
367
368
369
370template <int dim, int spacedim>
371void
373 const InternalData &data,
374 const CellSimilarity::Similarity cell_similarity,
376 &output_data) const
377{
378 // "compute" inverse Jacobian at the quadrature points, which are
379 // all the same
381 if (cell_similarity != CellSimilarity::translation)
382 for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
383 {
384 output_data.inverse_jacobians[i] = Tensor<2, dim>();
385 for (unsigned int j = 0; j < dim; ++j)
386 output_data.inverse_jacobians[i][j][j] =
387 data.inverse_cell_extents[j];
388 }
389}
390
391
392
393template <int dim, int spacedim>
397 const CellSimilarity::Similarity cell_similarity,
398 const Quadrature<dim> &quadrature,
399 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
401 &output_data) const
402{
403 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
404
405 // convert data object to internal data for this class. fails with
406 // an exception if that is not possible
407 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
408 ExcInternalError());
409 const InternalData &data = static_cast<const InternalData &>(internal_data);
410
411
412 update_cell_extents(cell, cell_similarity, data);
413
415 data,
416 quadrature.get_points(),
417 output_data.quadrature_points);
418
419 // compute Jacobian determinant. all values are equal and are the
420 // product of the local lengths in each coordinate direction
422 if (cell_similarity != CellSimilarity::translation)
423 {
424 double J = data.cell_extents[0];
425 for (unsigned int d = 1; d < dim; ++d)
426 J *= data.cell_extents[d];
427 data.volume_element = J;
429 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
430 output_data.JxW_values[i] = quadrature.weight(i);
431 }
432
433
434 maybe_update_jacobians(data, cell_similarity, output_data);
435 maybe_update_jacobian_derivatives(data, cell_similarity, output_data);
436 maybe_update_inverse_jacobians(data, cell_similarity, output_data);
437
438 return cell_similarity;
439}
440
441
442
443template <int dim, int spacedim>
444void
447 const unsigned int face_no,
448 const unsigned int subface_no,
449 const Quadrature<dim - 1> &quadrature,
450 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
452 &output_data) const
453{
454 (void)cell;
455 (void)face_no;
456 (void)subface_no;
457 (void)quadrature;
458 (void)internal_data;
459 (void)output_data;
461}
462
463
464
465template <int dim, int spacedim>
466void
470 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
472 &output_data) const
473{
474 AssertDimension(dim, spacedim);
475 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
476
477 // Convert data object to internal data for this class. Fails with an
478 // exception if that is not possible.
479 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
480 ExcInternalError());
481 const InternalData &data = static_cast<const InternalData &>(internal_data);
482
483
485
487 data,
488 quadrature.get_points(),
489 output_data.quadrature_points);
490
492 for (unsigned int i = 0; i < output_data.normal_vectors.size(); ++i)
493 output_data.normal_vectors[i] = quadrature.normal_vector(i);
494
496 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
497 output_data.JxW_values[i] = quadrature.weight(i);
498
503}
504
505
506
507template <int dim, int spacedim>
508void
510 const ArrayView<const Tensor<1, dim>> &input,
511 const MappingKind mapping_kind,
512 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
513 const ArrayView<Tensor<1, spacedim>> &output) const
514{
515 AssertDimension(input.size(), output.size());
516 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
517 ExcInternalError());
518 const InternalData &data = static_cast<const InternalData &>(mapping_data);
519
520 switch (mapping_kind)
521 {
523 {
526 "update_covariant_transformation"));
527
528 for (unsigned int i = 0; i < output.size(); ++i)
529 for (unsigned int d = 0; d < dim; ++d)
530 output[i][d] = input[i][d] * data.inverse_cell_extents[d];
531 return;
532 }
533
535 {
538 "update_contravariant_transformation"));
539
540 for (unsigned int i = 0; i < output.size(); ++i)
541 for (unsigned int d = 0; d < dim; ++d)
542 output[i][d] = input[i][d] * data.cell_extents[d];
543 return;
544 }
545 case mapping_piola:
546 {
549 "update_contravariant_transformation"));
552 "update_volume_elements"));
553
554 for (unsigned int i = 0; i < output.size(); ++i)
555 for (unsigned int d = 0; d < dim; ++d)
556 output[i][d] =
557 input[i][d] * data.cell_extents[d] / data.volume_element;
558 return;
559 }
560 default:
562 }
563}
564
565
566
567template <int dim, int spacedim>
568void
571 const MappingKind mapping_kind,
572 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
573 const ArrayView<Tensor<2, spacedim>> &output) const
574{
575 AssertDimension(input.size(), output.size());
576 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
577 ExcInternalError());
578 const InternalData &data = static_cast<const InternalData &>(mapping_data);
579
580 switch (mapping_kind)
581 {
583 {
586 "update_covariant_transformation"));
587
588 for (unsigned int i = 0; i < output.size(); ++i)
589 for (unsigned int d1 = 0; d1 < dim; ++d1)
590 for (unsigned int d2 = 0; d2 < dim; ++d2)
591 output[i][d1][d2] =
592 input[i][d1][d2] * data.inverse_cell_extents[d2];
593 return;
594 }
595
597 {
600 "update_contravariant_transformation"));
601
602 for (unsigned int i = 0; i < output.size(); ++i)
603 for (unsigned int d1 = 0; d1 < dim; ++d1)
604 for (unsigned int d2 = 0; d2 < dim; ++d2)
605 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
606 return;
607 }
608
610 {
613 "update_covariant_transformation"));
614
615 for (unsigned int i = 0; i < output.size(); ++i)
616 for (unsigned int d1 = 0; d1 < dim; ++d1)
617 for (unsigned int d2 = 0; d2 < dim; ++d2)
618 output[i][d1][d2] = input[i][d1][d2] *
619 data.inverse_cell_extents[d2] *
620 data.inverse_cell_extents[d1];
621 return;
622 }
623
625 {
628 "update_contravariant_transformation"));
629
630 for (unsigned int i = 0; i < output.size(); ++i)
631 for (unsigned int d1 = 0; d1 < dim; ++d1)
632 for (unsigned int d2 = 0; d2 < dim; ++d2)
633 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
634 data.inverse_cell_extents[d1];
635 return;
636 }
637
638 case mapping_piola:
639 {
642 "update_contravariant_transformation"));
645 "update_volume_elements"));
646
647 for (unsigned int i = 0; i < output.size(); ++i)
648 for (unsigned int d1 = 0; d1 < dim; ++d1)
649 for (unsigned int d2 = 0; d2 < dim; ++d2)
650 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
651 data.volume_element;
652 return;
653 }
654
656 {
659 "update_contravariant_transformation"));
662 "update_volume_elements"));
663
664 for (unsigned int i = 0; i < output.size(); ++i)
665 for (unsigned int d1 = 0; d1 < dim; ++d1)
666 for (unsigned int d2 = 0; d2 < dim; ++d2)
667 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
668 data.inverse_cell_extents[d1] /
669 data.volume_element;
670 return;
671 }
672
673 default:
675 }
676}
677
678
679
680template <int dim, int spacedim>
681void
683 const ArrayView<const Tensor<2, dim>> &input,
684 const MappingKind mapping_kind,
685 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
686 const ArrayView<Tensor<2, spacedim>> &output) const
687{
688 AssertDimension(input.size(), output.size());
689 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
690 ExcInternalError());
691 const InternalData &data = static_cast<const InternalData &>(mapping_data);
692
693 switch (mapping_kind)
694 {
696 {
699 "update_covariant_transformation"));
700
701 for (unsigned int i = 0; i < output.size(); ++i)
702 for (unsigned int d1 = 0; d1 < dim; ++d1)
703 for (unsigned int d2 = 0; d2 < dim; ++d2)
704 output[i][d1][d2] =
705 input[i][d1][d2] * data.inverse_cell_extents[d2];
706 return;
707 }
708
710 {
713 "update_contravariant_transformation"));
714
715 for (unsigned int i = 0; i < output.size(); ++i)
716 for (unsigned int d1 = 0; d1 < dim; ++d1)
717 for (unsigned int d2 = 0; d2 < dim; ++d2)
718 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
719 return;
720 }
721
723 {
726 "update_covariant_transformation"));
727
728 for (unsigned int i = 0; i < output.size(); ++i)
729 for (unsigned int d1 = 0; d1 < dim; ++d1)
730 for (unsigned int d2 = 0; d2 < dim; ++d2)
731 output[i][d1][d2] = input[i][d1][d2] *
732 data.inverse_cell_extents[d2] *
733 data.inverse_cell_extents[d1];
734 return;
735 }
736
738 {
741 "update_contravariant_transformation"));
742
743 for (unsigned int i = 0; i < output.size(); ++i)
744 for (unsigned int d1 = 0; d1 < dim; ++d1)
745 for (unsigned int d2 = 0; d2 < dim; ++d2)
746 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
747 data.inverse_cell_extents[d1];
748 return;
749 }
750
751 case mapping_piola:
752 {
755 "update_contravariant_transformation"));
758 "update_volume_elements"));
759
760 for (unsigned int i = 0; i < output.size(); ++i)
761 for (unsigned int d1 = 0; d1 < dim; ++d1)
762 for (unsigned int d2 = 0; d2 < dim; ++d2)
763 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
764 data.volume_element;
765 return;
766 }
767
769 {
772 "update_contravariant_transformation"));
775 "update_volume_elements"));
776
777 for (unsigned int i = 0; i < output.size(); ++i)
778 for (unsigned int d1 = 0; d1 < dim; ++d1)
779 for (unsigned int d2 = 0; d2 < dim; ++d2)
780 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
781 data.inverse_cell_extents[d1] /
782 data.volume_element;
783 return;
784 }
785
786 default:
788 }
789}
790
791
792
793template <int dim, int spacedim>
794void
797 const MappingKind mapping_kind,
798 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
799 const ArrayView<Tensor<3, spacedim>> &output) const
800{
801 AssertDimension(input.size(), output.size());
802 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
803 ExcInternalError());
804 const InternalData &data = static_cast<const InternalData &>(mapping_data);
805
806 switch (mapping_kind)
807 {
809 {
812 "update_covariant_transformation"));
813
814 for (unsigned int q = 0; q < output.size(); ++q)
815 for (unsigned int i = 0; i < spacedim; ++i)
816 for (unsigned int j = 0; j < spacedim; ++j)
817 for (unsigned int k = 0; k < spacedim; ++k)
818 {
819 output[q][i][j][k] = input[q][i][j][k] *
820 data.inverse_cell_extents[j] *
821 data.inverse_cell_extents[k];
822 }
823 return;
824 }
825 default:
827 }
828}
829
830
831
832template <int dim, int spacedim>
833void
835 const ArrayView<const Tensor<3, dim>> &input,
836 const MappingKind mapping_kind,
837 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
838 const ArrayView<Tensor<3, spacedim>> &output) const
839{
840 AssertDimension(input.size(), output.size());
841 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
842 ExcInternalError());
843 const InternalData &data = static_cast<const InternalData &>(mapping_data);
844
845 switch (mapping_kind)
846 {
848 {
851 "update_covariant_transformation"));
854 "update_contravariant_transformation"));
855
856 for (unsigned int q = 0; q < output.size(); ++q)
857 for (unsigned int i = 0; i < spacedim; ++i)
858 for (unsigned int j = 0; j < spacedim; ++j)
859 for (unsigned int k = 0; k < spacedim; ++k)
860 {
861 output[q][i][j][k] = input[q][i][j][k] *
862 data.cell_extents[i] *
863 data.inverse_cell_extents[j] *
864 data.inverse_cell_extents[k];
865 }
866 return;
867 }
868
870 {
873 "update_covariant_transformation"));
874
875 for (unsigned int q = 0; q < output.size(); ++q)
876 for (unsigned int i = 0; i < spacedim; ++i)
877 for (unsigned int j = 0; j < spacedim; ++j)
878 for (unsigned int k = 0; k < spacedim; ++k)
879 {
880 output[q][i][j][k] = input[q][i][j][k] *
881 (data.inverse_cell_extents[i] *
882 data.inverse_cell_extents[j]) *
883 data.inverse_cell_extents[k];
884 }
885
886 return;
887 }
888
890 {
893 "update_covariant_transformation"));
896 "update_contravariant_transformation"));
899 "update_volume_elements"));
900
901 for (unsigned int q = 0; q < output.size(); ++q)
902 for (unsigned int i = 0; i < spacedim; ++i)
903 for (unsigned int j = 0; j < spacedim; ++j)
904 for (unsigned int k = 0; k < spacedim; ++k)
905 {
906 output[q][i][j][k] =
907 input[q][i][j][k] *
908 (data.cell_extents[i] / data.volume_element *
909 data.inverse_cell_extents[j]) *
910 data.inverse_cell_extents[k];
911 }
912
913 return;
914 }
915
916 default:
918 }
919}
920
921
922
923template <int dim, int spacedim>
927 const Point<dim> &p) const
928{
929 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
930 Assert(dim == spacedim, ExcNotImplemented());
931
932 return boxes[polytope_translator.at(cell->active_cell_index())].unit_to_real(
933 p);
934}
935
936
937
938template <int dim, int spacedim>
942 const Point<spacedim> &p) const
943{
944 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
945 Assert(dim == spacedim, ExcNotImplemented());
946
947 return boxes[polytope_translator.at(cell->active_cell_index())].real_to_unit(
948 p);
949}
950
951
952
953template <int dim, int spacedim>
954void
957 const ArrayView<const Point<spacedim>> &real_points,
958 const ArrayView<Point<dim>> &unit_points) const
959{
960 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
961 AssertDimension(real_points.size(), unit_points.size());
962
963 if (dim != spacedim)
965 for (unsigned int i = 0; i < real_points.size(); ++i)
966 unit_points[i] =
967 boxes[polytope_translator.at(cell->active_cell_index())].real_to_unit(
968 real_points[i]);
969}
970
971
972
973template <int dim, int spacedim>
974std::unique_ptr<Mapping<dim, spacedim>>
976{
977 return std::make_unique<MappingBox<dim, spacedim>>(*this);
978}
979
980
981//---------------------------------------------------------------------------
982// explicit instantiations
983template class MappingBox<1>;
984template class MappingBox<2>;
985template class MappingBox<3>;
986
987
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
Number side_length(const unsigned int direction) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
double weight(const unsigned int i) const
const std::vector< Point< dim > > & get_points() const
virtual std::size_t memory_consumption() const override
Tensor< 1, dim > traslation
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Tensor< 1, dim > cell_extents
Tensor< 1, dim > inverse_cell_extents
void maybe_update_volume_elements(const InternalData &data) const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
MappingBox(const std::vector< BoundingBox< dim > > &local_boxes, const std::map< types::global_cell_index, types::global_cell_index > &polytope_translator)
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual bool preserves_vertex_locations() 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 Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) 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
void update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) 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
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 std::unique_ptr< Mapping< dim, spacedim > > clone() const override
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
std::vector< BoundingBox< dim > > boxes
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void maybe_update_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const ArrayView< const Point< dim > > &unit_quadrature_points, std::vector< Point< dim > > &quadrature_points) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const
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
std::map< types::global_cell_index, types::global_cell_index > polytope_translator
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
std::vector< Tensor< 1, spacedim > > normal_vectors
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< Point< spacedim > > quadrature_points
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
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
UpdateFlags
update_jacobian_pushed_forward_2nd_derivatives
update_volume_elements
update_contravariant_transformation
update_jacobian_pushed_forward_grads
update_jacobian_3rd_derivatives
update_jacobian_grads
update_normal_vectors
update_JxW_values
update_covariant_transformation
update_jacobians
update_inverse_jacobians
update_quadrature_points
update_jacobian_pushed_forward_3rd_derivatives
update_boundary_forms
update_jacobian_2nd_derivatives
MappingKind
mapping_piola
mapping_covariant_gradient
mapping_covariant
mapping_contravariant
mapping_contravariant_hessian
mapping_covariant_hessian
mapping_contravariant_gradient
mapping_piola_gradient
mapping_piola_hessian
#define DEAL_II_NOT_IMPLEMENTED()
DEAL_II_NAMESPACE_OPEN DeclExceptionMsg(ExcCellNotAssociatedWithBox, "You are using MappingBox, but the incoming element is not associated with a" "Bounding Box Cartesian.")
bool has_box(const CellType &cell, const std::map< types::global_cell_index, types::global_cell_index > &translator)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
double volume(const Triangulation< dim, spacedim > &tria)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)