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
133#if DEAL_II_VERSION_GTE(9, 8, 0)
136#else
138 const ReferenceCell &reference_cell) const
139#endif
140{
141 Assert(dim == reference_cell.get_dimension(),
142 ExcMessage("The dimension of your mapping (" +
144 ") and the reference cell cell_type (" +
145 Utilities::to_string(reference_cell.get_dimension()) +
146 " ) do not agree."));
147
148 return reference_cell.is_hyper_cube() || reference_cell.is_simplex();
149}
150
151
152
153template <int dim, int spacedim>
156{
157 // this mapping is pretty simple in that it can basically compute
158 // every piece of information wanted by FEValues without requiring
159 // computing any other quantities. boundary forms are one exception
160 // since they can be computed from the normal vectors without much
161 // further ado
162 UpdateFlags out = in;
163 if (out & update_boundary_forms)
165
166 return out;
167}
168
169
170
171template <int dim, int spacedim>
172std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
174 const Quadrature<dim> &q) const
175{
176 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
177 std::make_unique<InternalData>();
178 data_ptr->reinit(requires_update_flags(update_flags), q);
179
180 return data_ptr;
181}
182
183
184
185template <int dim, int spacedim>
186std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
188 const UpdateFlags update_flags,
189 const Quadrature<dim - 1> &quadrature) const
190{
191 (void)update_flags;
192 (void)quadrature;
194 return {};
195}
196
197
198
199template <int dim, int spacedim>
200void
203 const CellSimilarity::Similarity cell_similarity,
204 const InternalData &data) const
205{
206 // Compute start point and sizes along axes. The vertices to be looked at
207 // are 1, 2, 4 compared to the base vertex 0.
208 if (cell_similarity != CellSimilarity::translation)
209 {
210 const BoundingBox<dim> &current_box =
211 boxes[polytope_translator.at(cell->active_cell_index())];
212 const std::pair<Point<dim>, Point<dim>> &bdary_points =
213 current_box.get_boundary_points();
214
215 for (unsigned int d = 0; d < dim; ++d)
216 {
217 const double cell_extent_d = current_box.side_length(d);
218 data.cell_extents[d] = cell_extent_d;
219
220 data.traslation[d] =
221 .5 * (bdary_points.first[d] +
222 bdary_points.second[d]); // midpoint of each interval
223
224 Assert(cell_extent_d != 0.,
225 ExcMessage("Cell does not appear to be Cartesian!"));
226 data.inverse_cell_extents[d] = 1. / cell_extent_d;
227 }
228 }
229}
230
231
232
233namespace
234{
235 template <int dim>
236 void
237 transform_quadrature_points(
238 const BoundingBox<dim> &box,
239 const ArrayView<const Point<dim>> &unit_quadrature_points,
240 std::vector<Point<dim>> &quadrature_points)
241 {
242 for (unsigned int i = 0; i < quadrature_points.size(); ++i)
243 quadrature_points[i] = box.unit_to_real(unit_quadrature_points[i]);
244 }
245} // namespace
246
247
248
249template <int dim, int spacedim>
250void
253 const InternalData &data,
254 const ArrayView<const Point<dim>> &unit_quadrature_points,
255 std::vector<Point<dim>> &quadrature_points) const
256{
257 if (data.update_each & update_quadrature_points)
258 transform_quadrature_points(
259 boxes[polytope_translator.at(cell->active_cell_index())],
260 unit_quadrature_points,
261 quadrature_points);
262}
263
264
265
266template <int dim, int spacedim>
267void
269 const unsigned int face_no,
270 const InternalData &data,
271 std::vector<Tensor<1, dim>> &normal_vectors) const
272{
273 // compute normal vectors. All normals on a face have the same value.
274 if (data.update_each & update_normal_vectors)
275 {
276 Assert(face_no < GeometryInfo<dim>::faces_per_cell, ExcInternalError());
277 std::fill(normal_vectors.begin(),
278 normal_vectors.end(),
280 }
281}
282
283
284
285template <int dim, int spacedim>
286void
288 const InternalData &data,
289 const CellSimilarity::Similarity cell_similarity,
291 &output_data) const
292{
293 if (cell_similarity != CellSimilarity::translation)
294 {
295 if (data.update_each & update_jacobian_grads)
296 for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
298
300 for (unsigned int i = 0;
301 i < output_data.jacobian_pushed_forward_grads.size();
302 ++i)
304
305 if (data.update_each & update_jacobian_2nd_derivatives)
306 for (unsigned int i = 0;
307 i < output_data.jacobian_2nd_derivatives.size();
308 ++i)
309 output_data.jacobian_2nd_derivatives[i] =
311
313 for (unsigned int i = 0;
314 i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
315 ++i)
318
319 if (data.update_each & update_jacobian_3rd_derivatives)
320 for (unsigned int i = 0;
321 i < output_data.jacobian_3rd_derivatives.size();
322 ++i)
323 output_data.jacobian_3rd_derivatives[i] =
325
327 for (unsigned int i = 0;
328 i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
329 ++i)
332 }
333}
334
335
336
337template <int dim, int spacedim>
338void
340 const InternalData &data) const
341{
342 if (data.update_each & update_volume_elements)
343 {
344 double volume = data.cell_extents[0];
345 for (unsigned int d = 1; d < dim; ++d)
346 volume *= data.cell_extents[d];
347 data.volume_element = volume;
348 }
349}
350
351
352
353template <int dim, int spacedim>
354void
356 const InternalData &data,
357 const CellSimilarity::Similarity cell_similarity,
359 &output_data) const
360{
361 // "compute" Jacobian at the quadrature points, which are all the
362 // same
363 if (data.update_each & update_jacobians)
364 if (cell_similarity != CellSimilarity::translation)
365 for (unsigned int i = 0; i < output_data.jacobians.size(); ++i)
366 {
368 for (unsigned int j = 0; j < dim; ++j)
369 output_data.jacobians[i][j][j] = data.cell_extents[j];
370 }
371}
372
373
374
375template <int dim, int spacedim>
376void
378 const InternalData &data,
379 const CellSimilarity::Similarity cell_similarity,
381 &output_data) const
382{
383 // "compute" inverse Jacobian at the quadrature points, which are
384 // all the same
385 if (data.update_each & update_inverse_jacobians)
386 if (cell_similarity != CellSimilarity::translation)
387 for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
388 {
389 output_data.inverse_jacobians[i] = Tensor<2, dim>();
390 for (unsigned int j = 0; j < dim; ++j)
391 output_data.inverse_jacobians[i][j][j] =
392 data.inverse_cell_extents[j];
393 }
394}
395
396
397
398template <int dim, int spacedim>
402 const CellSimilarity::Similarity cell_similarity,
403 const Quadrature<dim> &quadrature,
404 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
406 &output_data) const
407{
408 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
409
410 // convert data object to internal data for this class. fails with
411 // an exception if that is not possible
412 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
413 ExcInternalError());
414 const InternalData &data = static_cast<const InternalData &>(internal_data);
415
416
417 update_cell_extents(cell, cell_similarity, data);
418
420 data,
421 quadrature.get_points(),
422 output_data.quadrature_points);
423
424 // compute Jacobian determinant. all values are equal and are the
425 // product of the local lengths in each coordinate direction
426 if (data.update_each & (update_JxW_values | update_volume_elements))
427 if (cell_similarity != CellSimilarity::translation)
428 {
429 double J = data.cell_extents[0];
430 for (unsigned int d = 1; d < dim; ++d)
431 J *= data.cell_extents[d];
432 data.volume_element = J;
433 if (data.update_each & update_JxW_values)
434 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
435 output_data.JxW_values[i] = quadrature.weight(i);
436 }
437
438
439 maybe_update_jacobians(data, cell_similarity, output_data);
440 maybe_update_jacobian_derivatives(data, cell_similarity, output_data);
441 maybe_update_inverse_jacobians(data, cell_similarity, output_data);
442
443 return cell_similarity;
444}
445
446
447
448template <int dim, int spacedim>
449void
452 const unsigned int face_no,
453 const unsigned int subface_no,
454 const Quadrature<dim - 1> &quadrature,
455 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
457 &output_data) const
458{
459 (void)cell;
460 (void)face_no;
461 (void)subface_no;
462 (void)quadrature;
463 (void)internal_data;
464 (void)output_data;
466}
467
468
469
470template <int dim, int spacedim>
471void
475 const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
477 &output_data) const
478{
479 AssertDimension(dim, spacedim);
480 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
481
482 // Convert data object to internal data for this class. Fails with an
483 // exception if that is not possible.
484 Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
485 ExcInternalError());
486 const InternalData &data = static_cast<const InternalData &>(internal_data);
487
488
490
492 data,
493 quadrature.get_points(),
494 output_data.quadrature_points);
495
496 if (data.update_each & update_normal_vectors)
497 for (unsigned int i = 0; i < output_data.normal_vectors.size(); ++i)
498 output_data.normal_vectors[i] = quadrature.normal_vector(i);
499
500 if (data.update_each & update_JxW_values)
501 for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
502 output_data.JxW_values[i] = quadrature.weight(i);
503
508}
509
510
511
512template <int dim, int spacedim>
513void
515 const ArrayView<const Tensor<1, dim>> &input,
516 const MappingKind mapping_kind,
517 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
518 const ArrayView<Tensor<1, spacedim>> &output) const
519{
520 AssertDimension(input.size(), output.size());
521 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
522 ExcInternalError());
523 const InternalData &data = static_cast<const InternalData &>(mapping_data);
524
525 switch (mapping_kind)
526 {
528 {
531 "update_covariant_transformation"));
532
533 for (unsigned int i = 0; i < output.size(); ++i)
534 for (unsigned int d = 0; d < dim; ++d)
535 output[i][d] = input[i][d] * data.inverse_cell_extents[d];
536 return;
537 }
538
540 {
543 "update_contravariant_transformation"));
544
545 for (unsigned int i = 0; i < output.size(); ++i)
546 for (unsigned int d = 0; d < dim; ++d)
547 output[i][d] = input[i][d] * data.cell_extents[d];
548 return;
549 }
550 case mapping_piola:
551 {
554 "update_contravariant_transformation"));
555 Assert(data.update_each & update_volume_elements,
557 "update_volume_elements"));
558
559 for (unsigned int i = 0; i < output.size(); ++i)
560 for (unsigned int d = 0; d < dim; ++d)
561 output[i][d] =
562 input[i][d] * data.cell_extents[d] / data.volume_element;
563 return;
564 }
565 default:
567 }
568}
569
570
571
572template <int dim, int spacedim>
573void
576 const MappingKind mapping_kind,
577 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
578 const ArrayView<Tensor<2, spacedim>> &output) const
579{
580 AssertDimension(input.size(), output.size());
581 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
582 ExcInternalError());
583 const InternalData &data = static_cast<const InternalData &>(mapping_data);
584
585 switch (mapping_kind)
586 {
588 {
591 "update_covariant_transformation"));
592
593 for (unsigned int i = 0; i < output.size(); ++i)
594 for (unsigned int d1 = 0; d1 < dim; ++d1)
595 for (unsigned int d2 = 0; d2 < dim; ++d2)
596 output[i][d1][d2] =
597 input[i][d1][d2] * data.inverse_cell_extents[d2];
598 return;
599 }
600
602 {
605 "update_contravariant_transformation"));
606
607 for (unsigned int i = 0; i < output.size(); ++i)
608 for (unsigned int d1 = 0; d1 < dim; ++d1)
609 for (unsigned int d2 = 0; d2 < dim; ++d2)
610 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
611 return;
612 }
613
615 {
618 "update_covariant_transformation"));
619
620 for (unsigned int i = 0; i < output.size(); ++i)
621 for (unsigned int d1 = 0; d1 < dim; ++d1)
622 for (unsigned int d2 = 0; d2 < dim; ++d2)
623 output[i][d1][d2] = input[i][d1][d2] *
624 data.inverse_cell_extents[d2] *
625 data.inverse_cell_extents[d1];
626 return;
627 }
628
630 {
633 "update_contravariant_transformation"));
634
635 for (unsigned int i = 0; i < output.size(); ++i)
636 for (unsigned int d1 = 0; d1 < dim; ++d1)
637 for (unsigned int d2 = 0; d2 < dim; ++d2)
638 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
639 data.inverse_cell_extents[d1];
640 return;
641 }
642
643 case mapping_piola:
644 {
647 "update_contravariant_transformation"));
648 Assert(data.update_each & update_volume_elements,
650 "update_volume_elements"));
651
652 for (unsigned int i = 0; i < output.size(); ++i)
653 for (unsigned int d1 = 0; d1 < dim; ++d1)
654 for (unsigned int d2 = 0; d2 < dim; ++d2)
655 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
656 data.volume_element;
657 return;
658 }
659
661 {
664 "update_contravariant_transformation"));
665 Assert(data.update_each & update_volume_elements,
667 "update_volume_elements"));
668
669 for (unsigned int i = 0; i < output.size(); ++i)
670 for (unsigned int d1 = 0; d1 < dim; ++d1)
671 for (unsigned int d2 = 0; d2 < dim; ++d2)
672 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
673 data.inverse_cell_extents[d1] /
674 data.volume_element;
675 return;
676 }
677
678 default:
680 }
681}
682
683
684
685template <int dim, int spacedim>
686void
688 const ArrayView<const Tensor<2, dim>> &input,
689 const MappingKind mapping_kind,
690 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
691 const ArrayView<Tensor<2, spacedim>> &output) const
692{
693 AssertDimension(input.size(), output.size());
694 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
695 ExcInternalError());
696 const InternalData &data = static_cast<const InternalData &>(mapping_data);
697
698 switch (mapping_kind)
699 {
701 {
704 "update_covariant_transformation"));
705
706 for (unsigned int i = 0; i < output.size(); ++i)
707 for (unsigned int d1 = 0; d1 < dim; ++d1)
708 for (unsigned int d2 = 0; d2 < dim; ++d2)
709 output[i][d1][d2] =
710 input[i][d1][d2] * data.inverse_cell_extents[d2];
711 return;
712 }
713
715 {
718 "update_contravariant_transformation"));
719
720 for (unsigned int i = 0; i < output.size(); ++i)
721 for (unsigned int d1 = 0; d1 < dim; ++d1)
722 for (unsigned int d2 = 0; d2 < dim; ++d2)
723 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
724 return;
725 }
726
728 {
731 "update_covariant_transformation"));
732
733 for (unsigned int i = 0; i < output.size(); ++i)
734 for (unsigned int d1 = 0; d1 < dim; ++d1)
735 for (unsigned int d2 = 0; d2 < dim; ++d2)
736 output[i][d1][d2] = input[i][d1][d2] *
737 data.inverse_cell_extents[d2] *
738 data.inverse_cell_extents[d1];
739 return;
740 }
741
743 {
746 "update_contravariant_transformation"));
747
748 for (unsigned int i = 0; i < output.size(); ++i)
749 for (unsigned int d1 = 0; d1 < dim; ++d1)
750 for (unsigned int d2 = 0; d2 < dim; ++d2)
751 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
752 data.inverse_cell_extents[d1];
753 return;
754 }
755
756 case mapping_piola:
757 {
760 "update_contravariant_transformation"));
761 Assert(data.update_each & update_volume_elements,
763 "update_volume_elements"));
764
765 for (unsigned int i = 0; i < output.size(); ++i)
766 for (unsigned int d1 = 0; d1 < dim; ++d1)
767 for (unsigned int d2 = 0; d2 < dim; ++d2)
768 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
769 data.volume_element;
770 return;
771 }
772
774 {
777 "update_contravariant_transformation"));
778 Assert(data.update_each & update_volume_elements,
780 "update_volume_elements"));
781
782 for (unsigned int i = 0; i < output.size(); ++i)
783 for (unsigned int d1 = 0; d1 < dim; ++d1)
784 for (unsigned int d2 = 0; d2 < dim; ++d2)
785 output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] *
786 data.inverse_cell_extents[d1] /
787 data.volume_element;
788 return;
789 }
790
791 default:
793 }
794}
795
796
797
798template <int dim, int spacedim>
799void
802 const MappingKind mapping_kind,
803 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
804 const ArrayView<Tensor<3, spacedim>> &output) const
805{
806 AssertDimension(input.size(), output.size());
807 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
808 ExcInternalError());
809 const InternalData &data = static_cast<const InternalData &>(mapping_data);
810
811 switch (mapping_kind)
812 {
814 {
817 "update_covariant_transformation"));
818
819 for (unsigned int q = 0; q < output.size(); ++q)
820 for (unsigned int i = 0; i < spacedim; ++i)
821 for (unsigned int j = 0; j < spacedim; ++j)
822 for (unsigned int k = 0; k < spacedim; ++k)
823 {
824 output[q][i][j][k] = input[q][i][j][k] *
825 data.inverse_cell_extents[j] *
826 data.inverse_cell_extents[k];
827 }
828 return;
829 }
830 default:
832 }
833}
834
835
836
837template <int dim, int spacedim>
838void
840 const ArrayView<const Tensor<3, dim>> &input,
841 const MappingKind mapping_kind,
842 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
843 const ArrayView<Tensor<3, spacedim>> &output) const
844{
845 AssertDimension(input.size(), output.size());
846 Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
847 ExcInternalError());
848 const InternalData &data = static_cast<const InternalData &>(mapping_data);
849
850 switch (mapping_kind)
851 {
853 {
856 "update_covariant_transformation"));
859 "update_contravariant_transformation"));
860
861 for (unsigned int q = 0; q < output.size(); ++q)
862 for (unsigned int i = 0; i < spacedim; ++i)
863 for (unsigned int j = 0; j < spacedim; ++j)
864 for (unsigned int k = 0; k < spacedim; ++k)
865 {
866 output[q][i][j][k] = input[q][i][j][k] *
867 data.cell_extents[i] *
868 data.inverse_cell_extents[j] *
869 data.inverse_cell_extents[k];
870 }
871 return;
872 }
873
875 {
878 "update_covariant_transformation"));
879
880 for (unsigned int q = 0; q < output.size(); ++q)
881 for (unsigned int i = 0; i < spacedim; ++i)
882 for (unsigned int j = 0; j < spacedim; ++j)
883 for (unsigned int k = 0; k < spacedim; ++k)
884 {
885 output[q][i][j][k] = input[q][i][j][k] *
886 (data.inverse_cell_extents[i] *
887 data.inverse_cell_extents[j]) *
888 data.inverse_cell_extents[k];
889 }
890
891 return;
892 }
893
895 {
898 "update_covariant_transformation"));
901 "update_contravariant_transformation"));
902 Assert(data.update_each & update_volume_elements,
904 "update_volume_elements"));
905
906 for (unsigned int q = 0; q < output.size(); ++q)
907 for (unsigned int i = 0; i < spacedim; ++i)
908 for (unsigned int j = 0; j < spacedim; ++j)
909 for (unsigned int k = 0; k < spacedim; ++k)
910 {
911 output[q][i][j][k] =
912 input[q][i][j][k] *
913 (data.cell_extents[i] / data.volume_element *
914 data.inverse_cell_extents[j]) *
915 data.inverse_cell_extents[k];
916 }
917
918 return;
919 }
920
921 default:
923 }
924}
925
926
927
928template <int dim, int spacedim>
932 const Point<dim> &p) const
933{
934 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
935 Assert(dim == spacedim, ExcNotImplemented());
936
937 return boxes[polytope_translator.at(cell->active_cell_index())].unit_to_real(
938 p);
939}
940
941
942
943template <int dim, int spacedim>
947 const Point<spacedim> &p) const
948{
949 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
950 Assert(dim == spacedim, ExcNotImplemented());
951
952 return boxes[polytope_translator.at(cell->active_cell_index())].real_to_unit(
953 p);
954}
955
956
957
958template <int dim, int spacedim>
959void
962 const ArrayView<const Point<spacedim>> &real_points,
963 const ArrayView<Point<dim>> &unit_points) const
964{
965 Assert(has_box(cell, polytope_translator), ExcCellNotAssociatedWithBox());
966 AssertDimension(real_points.size(), unit_points.size());
967
968 if (dim != spacedim)
970 for (unsigned int i = 0; i < real_points.size(); ++i)
971 unit_points[i] =
972 boxes[polytope_translator.at(cell->active_cell_index())].real_to_unit(
973 real_points[i]);
974}
975
976
977
978template <int dim, int spacedim>
979std::unique_ptr<Mapping<dim, spacedim>>
981{
982 return std::make_unique<MappingBox<dim, spacedim>>(*this);
983}
984
985
986//---------------------------------------------------------------------------
987// explicit instantiations
988template class MappingBox<1>;
989template class MappingBox<2>;
990template class MappingBox<3>;
991
992
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
virtual std::size_t memory_consumption() 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
std::vector< Point< dim > > quadrature_points
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 DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
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
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)
std::vector< index_type > data
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell< dim > &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)
static constexpr unsigned int faces_per_cell
static constexpr std::array< Tensor< 1, dim >, faces_per_cell > unit_normal_vector