Reference documentation for deal.II version 9.6.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
scratch_data.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2023 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
16
17#include <memory>
18
20
21namespace MeshWorker
22{
23 template <int dim, int spacedim>
25 const Mapping<dim, spacedim> &mapping,
27 const Quadrature<dim> &quadrature,
28 const UpdateFlags &update_flags,
29 const Quadrature<dim - 1> &face_quadrature,
30 const UpdateFlags &face_update_flags)
31 : mapping(&mapping)
32 , fe(&fe)
33 , cell_quadrature(quadrature)
34 , face_quadrature(face_quadrature)
35 , hp_capability_enabled(false)
36 , cell_update_flags(update_flags)
37 , neighbor_cell_update_flags(update_flags)
38 , face_update_flags(face_update_flags)
39 , neighbor_face_update_flags(face_update_flags)
40 , local_dof_indices(fe.n_dofs_per_cell())
41 , neighbor_dof_indices(fe.n_dofs_per_cell())
42 {}
43
44
45
46 template <int dim, int spacedim>
48 const Mapping<dim, spacedim> &mapping,
50 const Quadrature<dim> &quadrature,
51 const UpdateFlags &update_flags,
52 const UpdateFlags &neighbor_update_flags,
53 const Quadrature<dim - 1> &face_quadrature,
54 const UpdateFlags &face_update_flags,
55 const UpdateFlags &neighbor_face_update_flags)
56 : mapping(&mapping)
57 , fe(&fe)
58 , cell_quadrature(quadrature)
59 , face_quadrature(face_quadrature)
60 , hp_capability_enabled(false)
61 , cell_update_flags(update_flags)
62 , neighbor_cell_update_flags(neighbor_update_flags)
63 , face_update_flags(face_update_flags)
64 , neighbor_face_update_flags(neighbor_face_update_flags)
65 , local_dof_indices(fe.n_dofs_per_cell())
66 , neighbor_dof_indices(fe.n_dofs_per_cell())
67 {}
68
69
70
71 template <int dim, int spacedim>
74 const Quadrature<dim> &quadrature,
75 const UpdateFlags &update_flags,
76 const Quadrature<dim - 1> &face_quadrature,
77 const UpdateFlags &face_update_flags)
78 : ScratchData(fe.reference_cell()
79 .template get_default_linear_mapping<dim, spacedim>(),
80 fe,
81 quadrature,
82 update_flags,
83 face_quadrature,
84 face_update_flags)
85 {}
86
87
88
89 template <int dim, int spacedim>
92 const Quadrature<dim> &quadrature,
93 const UpdateFlags &update_flags,
94 const UpdateFlags &neighbor_update_flags,
95 const Quadrature<dim - 1> &face_quadrature,
96 const UpdateFlags &face_update_flags,
97 const UpdateFlags &neighbor_face_update_flags)
98 : ScratchData(fe.reference_cell()
99 .template get_default_linear_mapping<dim, spacedim>(),
100 fe,
101 quadrature,
102 update_flags,
103 neighbor_update_flags,
104 face_quadrature,
105 face_update_flags,
106 neighbor_face_update_flags)
107 {}
108
109
110
111 template <int dim, int spacedim>
113 const hp::MappingCollection<dim, spacedim> &mapping_collection,
114 const hp::FECollection<dim, spacedim> &fe_collection,
115 const hp::QCollection<dim> &cell_quadrature_collection,
116 const UpdateFlags &cell_update_flags,
117 const hp::QCollection<dim - 1> &face_quadrature_collection,
118 const UpdateFlags &face_update_flags)
119 : mapping_collection(&mapping_collection)
120 , fe_collection(&fe_collection)
121 , cell_quadrature_collection(cell_quadrature_collection)
122 , face_quadrature_collection(face_quadrature_collection)
123 , hp_capability_enabled(true)
124 , cell_update_flags(cell_update_flags)
125 , neighbor_cell_update_flags(cell_update_flags)
126 , face_update_flags(face_update_flags)
127 , neighbor_face_update_flags(face_update_flags)
128 {
129 local_dof_indices.reserve(fe_collection.max_dofs_per_cell());
130 neighbor_dof_indices.reserve(fe_collection.max_dofs_per_cell());
131 }
132
133
134
135 template <int dim, int spacedim>
137 const hp::MappingCollection<dim, spacedim> &mapping_collection,
138 const hp::FECollection<dim, spacedim> &fe_collection,
139 const hp::QCollection<dim> &cell_quadrature_collection,
140 const UpdateFlags &cell_update_flags,
141 const UpdateFlags &neighbor_cell_update_flags,
142 const hp::QCollection<dim - 1> &face_quadrature_collection,
143 const UpdateFlags &face_update_flags,
144 const UpdateFlags &neighbor_face_update_flags)
145 : mapping_collection(&mapping_collection)
146 , fe_collection(&fe_collection)
147 , cell_quadrature_collection(cell_quadrature_collection)
148 , face_quadrature_collection(face_quadrature_collection)
149 , hp_capability_enabled(true)
150 , cell_update_flags(cell_update_flags)
151 , neighbor_cell_update_flags(neighbor_cell_update_flags)
152 , face_update_flags(face_update_flags)
153 , neighbor_face_update_flags(neighbor_face_update_flags)
154 {
155 local_dof_indices.reserve(fe_collection.max_dofs_per_cell());
156 neighbor_dof_indices.reserve(fe_collection.max_dofs_per_cell());
157 }
158
159
160
161 template <int dim, int spacedim>
163 const hp::FECollection<dim, spacedim> &fe_collection,
164 const hp::QCollection<dim> &cell_quadrature_collection,
165 const UpdateFlags &cell_update_flags,
166 const hp::QCollection<dim - 1> &face_quadrature_collection,
167 const UpdateFlags &face_update_flags)
168 : ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
169 fe_collection,
170 cell_quadrature_collection,
171 cell_update_flags,
172 face_quadrature_collection,
173 face_update_flags)
174 {}
175
176
177
178 template <int dim, int spacedim>
180 const hp::FECollection<dim, spacedim> &fe_collection,
181 const hp::QCollection<dim> &cell_quadrature_collection,
182 const UpdateFlags &cell_update_flags,
183 const UpdateFlags &neighbor_cell_update_flags,
184 const hp::QCollection<dim - 1> &face_quadrature_collection,
185 const UpdateFlags &face_update_flags,
186 const UpdateFlags &neighbor_face_update_flags)
187 : ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
188 fe_collection,
189 cell_quadrature_collection,
190 cell_update_flags,
191 neighbor_cell_update_flags,
192 face_quadrature_collection,
193 face_update_flags,
194 neighbor_face_update_flags)
195 {}
196
197
198
199 template <int dim, int spacedim>
201 const ScratchData<dim, spacedim> &scratch)
202 : mapping(scratch.mapping)
203 , fe(scratch.fe)
204 , cell_quadrature(scratch.cell_quadrature)
205 , face_quadrature(scratch.face_quadrature)
206 , mapping_collection(scratch.mapping_collection)
207 , fe_collection(scratch.fe_collection)
208 , cell_quadrature_collection(scratch.cell_quadrature_collection)
209 , face_quadrature_collection(scratch.face_quadrature_collection)
210 , hp_capability_enabled(scratch.hp_capability_enabled)
211 , cell_update_flags(scratch.cell_update_flags)
212 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
213 , face_update_flags(scratch.face_update_flags)
214 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
215 , local_dof_indices(scratch.local_dof_indices)
216 , neighbor_dof_indices(scratch.neighbor_dof_indices)
217 , user_data_storage(scratch.user_data_storage)
218 , internal_data_storage(scratch.internal_data_storage)
219 {}
220
221
222
223 template <int dim, int spacedim>
227 {
228 if (hp_capability_enabled == false)
229 {
230 if (!fe_values)
231 fe_values = std::make_unique<FEValues<dim, spacedim>>(
232 *mapping, *fe, cell_quadrature, cell_update_flags);
233
234 fe_values->reinit(cell);
235 local_dof_indices.resize(fe_values->dofs_per_cell);
236 cell->get_dof_indices(local_dof_indices);
237 current_fe_values = fe_values.get();
238 return *fe_values;
239 }
240 else
241 {
242 if (!hp_fe_values)
243 hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
244 *mapping_collection,
245 *fe_collection,
246 cell_quadrature_collection,
247 cell_update_flags);
248
249 hp_fe_values->reinit(cell);
250 const auto &fe_values = hp_fe_values->get_present_fe_values();
251
253 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
254 fe_values.dofs_per_cell);
255 local_dof_indices.resize(fe_values.dofs_per_cell);
256 cell->get_dof_indices(local_dof_indices);
257
258 current_fe_values = &fe_values;
259 return fe_values;
260 }
261 }
262
263
264
265 template <int dim, int spacedim>
269 const unsigned int face_no)
270 {
271 if (hp_capability_enabled == false)
272 {
273 if (!fe_face_values)
274 fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
275 *mapping, *fe, face_quadrature, face_update_flags);
276
277 fe_face_values->reinit(cell, face_no);
278 local_dof_indices.resize(fe->n_dofs_per_cell());
279 cell->get_dof_indices(local_dof_indices);
280 current_fe_values = fe_face_values.get();
281 return *fe_face_values;
282 }
283 else
284 {
285 return reinit(cell, cell, face_no);
286 }
287 }
288
289
290
291 template <int dim, int spacedim>
296 &neighbor_cell,
297 const unsigned int face_no)
298 {
299 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
300
301 if (!hp_fe_face_values)
302 hp_fe_face_values = std::make_unique<hp::FEFaceValues<dim, spacedim>>(
303 *mapping_collection,
304 *fe_collection,
305 face_quadrature_collection,
306 face_update_flags);
307
308 if (neighbor_cell == cell)
309 {
310 hp_fe_face_values->reinit(cell, face_no);
311 }
312 else
313 {
314 // When we want to ensure some agreement between the cell face and its
315 // neighbor on the quadrature order and mapping to use on this face,
316 // then we defer to the dominance of one FE over another. This should
317 // ensure that the optimal integration order and mapping order are
318 // selected for this situation.
319 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
320 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
321
322 // TODO: find_dominated_fe returns invalid_fe_index when no dominated FE
323 // has been found. We want to pass this value to FEFaceValues, but it
324 // expects an invalid_unsigned_int in this case. We need to match the
325 // interfaces in the future.
326 if (dominated_fe_index == numbers::invalid_fe_index)
327 dominated_fe_index = numbers::invalid_unsigned_int;
328
329 hp_fe_face_values->reinit(cell,
330 face_no,
331 dominated_fe_index,
332 dominated_fe_index);
333 }
334
335 const auto &fe_face_values = hp_fe_face_values->get_present_fe_values();
336 const auto &fe = (*fe_collection)[cell->active_fe_index()];
337
338 local_dof_indices.resize(fe.n_dofs_per_cell());
339 cell->get_dof_indices(local_dof_indices);
340
341 current_fe_values = &fe_face_values;
342 return fe_face_values;
343 }
344
345
346
347 template <int dim, int spacedim>
351 const unsigned int face_no,
352 const unsigned int subface_no)
353 {
354 if (subface_no != numbers::invalid_unsigned_int)
355 {
356 if (hp_capability_enabled == false)
357 {
358 if (!fe_subface_values)
359 fe_subface_values =
360 std::make_unique<FESubfaceValues<dim, spacedim>>(
361 *mapping, *fe, face_quadrature, face_update_flags);
362 fe_subface_values->reinit(cell, face_no, subface_no);
363 local_dof_indices.resize(fe->n_dofs_per_cell());
364 cell->get_dof_indices(local_dof_indices);
365
366 current_fe_values = fe_subface_values.get();
367 return *fe_subface_values;
368 }
369 else
370 {
371 return reinit(cell, cell, face_no, subface_no);
372 }
373 }
374 else
375 return reinit(cell, face_no);
376 }
377
378
379
380 template <int dim, int spacedim>
385 &neighbor_cell,
386 const unsigned int face_no,
387 const unsigned int subface_no)
388 {
389 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
390
391 if (subface_no != numbers::invalid_unsigned_int)
392 {
393 if (!hp_fe_subface_values)
394 hp_fe_subface_values =
395 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
396 *mapping_collection,
397 *fe_collection,
398 face_quadrature_collection,
399 face_update_flags);
400
401 if (neighbor_cell == cell)
402 {
403 hp_fe_subface_values->reinit(cell, face_no, subface_no);
404 }
405 else
406 {
407 // When we want to ensure some agreement between the cell face and
408 // its neighbor on the quadrature order and mapping to use on this
409 // face, then we defer to the dominance of one FE over another. This
410 // should ensure that the optimal integration order and mapping
411 // order are selected for this situation.
412 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
413 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
414
415 // TODO: find_dominated_fe returns invalid_fe_index when no
416 // dominated FE has been found. We want to pass this value to
417 // FEFaceValues, but it expects an invalid_unsigned_int in this
418 // case. We need to match the interfaces in the future.
419 if (dominated_fe_index == numbers::invalid_fe_index)
420 dominated_fe_index = numbers::invalid_unsigned_int;
421
422 hp_fe_subface_values->reinit(cell,
423 face_no,
424 subface_no,
425 dominated_fe_index,
426 dominated_fe_index);
427 }
428
429 const auto &fe_subface_values =
430 hp_fe_subface_values->get_present_fe_values();
431 const auto &fe = (*fe_collection)[cell->active_fe_index()];
432
433 local_dof_indices.resize(fe.n_dofs_per_cell());
434 cell->get_dof_indices(local_dof_indices);
435
436 current_fe_values = &fe_subface_values;
437 return fe_subface_values;
438 }
439 else
440 {
441 return reinit(cell, neighbor_cell, face_no);
442 }
443 }
444
445
446
447 template <int dim, int spacedim>
451 const unsigned int face_no,
453 &cell_neighbor,
454 const unsigned int face_no_neighbor)
455 {
456 return reinit(cell,
457 face_no,
459 cell_neighbor,
460 face_no_neighbor,
462 }
463
464
465
466 template <int dim, int spacedim>
470 const unsigned int face_no,
471 const unsigned int sub_face_no,
473 &cell_neighbor,
474 const unsigned int face_no_neighbor,
475 const unsigned int sub_face_no_neighbor)
476 {
477 if (hp_capability_enabled == false)
478 {
479 if (!interface_fe_values)
480 interface_fe_values =
481 std::make_unique<FEInterfaceValues<dim, spacedim>>(
482 *mapping, *fe, face_quadrature, face_update_flags);
483
484 interface_fe_values->reinit(cell,
485 face_no,
486 sub_face_no,
487 cell_neighbor,
488 face_no_neighbor,
489 sub_face_no_neighbor);
490 }
491 else
492 {
493 if (!interface_fe_values)
494 interface_fe_values =
495 std::make_unique<FEInterfaceValues<dim, spacedim>>(
496 *mapping_collection,
497 *fe_collection,
498 face_quadrature_collection,
499 face_update_flags);
500
501 // When we want to ensure some agreement between the cell face and
502 // its neighbor on the quadrature order and mapping to use on this
503 // face, then we defer to the dominance of one FE over another. This
504 // should ensure that the optimal integration order and mapping
505 // order are selected for this situation.
506 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
507 {cell->active_fe_index(), cell_neighbor->active_fe_index()});
508
509 // TODO: find_dominated_fe returns invalid_fe_index when no dominated FE
510 // has been found. We want to pass this value to FEFaceValues, but it
511 // expects an invalid_unsigned_int in this case. We need to match the
512 // interfaces in the future.
513 if (dominated_fe_index == numbers::invalid_fe_index)
514 dominated_fe_index = numbers::invalid_unsigned_int;
515
516 interface_fe_values->reinit(cell,
517 face_no,
518 sub_face_no,
519 cell_neighbor,
520 face_no_neighbor,
521 sub_face_no_neighbor,
522 dominated_fe_index,
523 dominated_fe_index);
524 }
525
526 current_fe_values = &interface_fe_values->get_fe_face_values(0);
527 current_neighbor_fe_values = &interface_fe_values->get_fe_face_values(1);
528
529 local_dof_indices = interface_fe_values->get_interface_dof_indices();
530 return *interface_fe_values;
531 }
532
533
534
535 template <int dim, int spacedim>
539 {
540 if (hp_capability_enabled == false)
541 {
542 if (!neighbor_fe_values)
543 neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
544 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
545
546 neighbor_fe_values->reinit(cell);
547 cell->get_dof_indices(neighbor_dof_indices);
548 current_neighbor_fe_values = neighbor_fe_values.get();
549 return *neighbor_fe_values;
550 }
551 else
552 {
553 if (!neighbor_hp_fe_values)
554 neighbor_hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
555 *mapping_collection,
556 *fe_collection,
557 cell_quadrature_collection,
558 neighbor_cell_update_flags);
559
560 neighbor_hp_fe_values->reinit(cell);
561 const auto &neighbor_fe_values =
562 neighbor_hp_fe_values->get_present_fe_values();
563
565 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
566 neighbor_fe_values.dofs_per_cell);
567 neighbor_dof_indices.resize(neighbor_fe_values.dofs_per_cell);
568 cell->get_dof_indices(neighbor_dof_indices);
569
570 current_neighbor_fe_values = &neighbor_fe_values;
571 return neighbor_fe_values;
572 }
573 }
574
575
576
577 template <int dim, int spacedim>
581 const unsigned int face_no)
582 {
583 if (hp_capability_enabled == false)
584 {
585 if (!neighbor_fe_face_values)
586 neighbor_fe_face_values =
587 std::make_unique<FEFaceValues<dim, spacedim>>(
588 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
589 neighbor_fe_face_values->reinit(cell, face_no);
590 cell->get_dof_indices(neighbor_dof_indices);
591 current_neighbor_fe_values = neighbor_fe_face_values.get();
592 return *neighbor_fe_face_values;
593 }
594 else
595 {
596 return reinit_neighbor(cell, cell, face_no);
597 }
598 }
599
600
601
602 template <int dim, int spacedim>
607 &neighbor_cell,
608 const unsigned int face_no)
609 {
610 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
611
612 if (!neighbor_hp_fe_face_values)
613 neighbor_hp_fe_face_values =
614 std::make_unique<hp::FEFaceValues<dim, spacedim>>(
615 *mapping_collection,
616 *fe_collection,
617 face_quadrature_collection,
618 neighbor_face_update_flags);
619
620 if (neighbor_cell == cell)
621 {
622 neighbor_hp_fe_face_values->reinit(neighbor_cell, face_no);
623 }
624 else
625 {
626 // When we want to ensure some agreement between the cell face and its
627 // neighbor on the quadrature order and mapping to use on this face,
628 // then we defer to the dominance of one FE over another. This should
629 // ensure that the optimal integration order and mapping order are
630 // selected for this situation.
631 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
632 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
633
634 // TODO: find_dominated_fe returns invalid_fe_index when no dominated FE
635 // has been found. We want to pass this value to FEFaceValues, but it
636 // expects an invalid_unsigned_int in this case. We need to match the
637 // interfaces in the future.
638 if (dominated_fe_index == numbers::invalid_fe_index)
639 dominated_fe_index = numbers::invalid_unsigned_int;
640
641 neighbor_hp_fe_face_values->reinit(neighbor_cell,
642 face_no,
643 dominated_fe_index,
644 dominated_fe_index);
645 }
646
647 const auto &neighbor_fe_face_values =
648 neighbor_hp_fe_face_values->get_present_fe_values();
649 const auto &neighbor_fe =
650 (*fe_collection)[neighbor_cell->active_fe_index()];
651
652 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
653 neighbor_cell->get_dof_indices(neighbor_dof_indices);
654
655 current_neighbor_fe_values = &neighbor_fe_face_values;
656 return neighbor_fe_face_values;
657 }
658
659
660
661 template <int dim, int spacedim>
665 const unsigned int face_no,
666 const unsigned int subface_no)
667 {
668 if (subface_no != numbers::invalid_unsigned_int)
669 {
670 if (hp_capability_enabled == false)
671 {
672 if (!neighbor_fe_subface_values)
673 neighbor_fe_subface_values =
674 std::make_unique<FESubfaceValues<dim, spacedim>>(
675 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
676 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
677 cell->get_dof_indices(neighbor_dof_indices);
678 current_neighbor_fe_values = neighbor_fe_subface_values.get();
679 return *neighbor_fe_subface_values;
680 }
681 else
682 {
683 return reinit_neighbor(cell, cell, face_no, subface_no);
684 }
685 }
686 else
687 return reinit_neighbor(cell, face_no);
688 }
689
690
691
692 template <int dim, int spacedim>
697 &neighbor_cell,
698 const unsigned int face_no,
699 const unsigned int subface_no)
700 {
701 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
702
703 if (subface_no != numbers::invalid_unsigned_int)
704 {
705 if (!neighbor_hp_fe_subface_values)
706 neighbor_hp_fe_subface_values =
707 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
708 *mapping_collection,
709 *fe_collection,
710 face_quadrature_collection,
711 neighbor_face_update_flags);
712
713 if (neighbor_cell == cell)
714 {
715 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
716 face_no,
717 subface_no);
718 }
719 else
720 {
721 // When we want to ensure some agreement between the cell face and
722 // its neighbor on the quadrature order and mapping to use on this
723 // face, then we defer to the dominance of one FE over another. This
724 // should ensure that the optimal integration order and mapping
725 // order are selected for this situation.
726 unsigned int dominated_fe_index = fe_collection->find_dominated_fe(
727 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
728
729 // TODO: find_dominated_fe returns invalid_fe_index when no
730 // dominated FE has been found. We want to pass this value to
731 // FEFaceValues, but it expects an invalid_unsigned_int in this
732 // case. We need to match the interfaces in the future.
733 if (dominated_fe_index == numbers::invalid_fe_index)
734 dominated_fe_index = numbers::invalid_unsigned_int;
735
736 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
737 face_no,
738 subface_no,
739 dominated_fe_index,
740 dominated_fe_index);
741 }
742
743 const auto &neighbor_fe_subface_values =
744 neighbor_hp_fe_subface_values->get_present_fe_values();
745 const auto &neighbor_fe =
746 (*fe_collection)[neighbor_cell->active_fe_index()];
747
748 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
749 neighbor_cell->get_dof_indices(neighbor_dof_indices);
750
751 current_neighbor_fe_values = &neighbor_fe_subface_values;
752 return neighbor_fe_subface_values;
753 }
754 else
755 {
756 return reinit_neighbor(cell, neighbor_cell, face_no);
757 }
758 }
759
760
761
762 template <int dim, int spacedim>
765 {
766 Assert(current_fe_values != nullptr,
767 ExcMessage("You have to initialize the cache using one of the "
768 "reinit functions first!"));
769 return *current_fe_values;
770 }
771
772
773
774 template <int dim, int spacedim>
777 {
778 Assert(interface_fe_values != nullptr,
779 ExcMessage("You have to initialize the cache using one of the "
780 "reinit functions first!"));
781 return *interface_fe_values;
782 }
783
784
785
786 template <int dim, int spacedim>
789 {
790 Assert(current_neighbor_fe_values != nullptr,
791 ExcMessage("You have to initialize the cache using one of the "
792 "reinit functions first!"));
793 return *current_neighbor_fe_values;
794 }
795
796
797
798 template <int dim, int spacedim>
799 const std::vector<Point<spacedim>> &
801 {
802 return get_current_fe_values().get_quadrature_points();
803 }
804
805
806
807 template <int dim, int spacedim>
808 const std::vector<double> &
810 {
811 return get_current_fe_values().get_JxW_values();
812 }
813
814
815
816 template <int dim, int spacedim>
817 const std::vector<double> &
819 {
820 return get_current_neighbor_fe_values().get_JxW_values();
821 }
822
823
824
825 template <int dim, int spacedim>
826 const std::vector<Tensor<1, spacedim>> &
828 {
829 return get_current_fe_values().get_normal_vectors();
830 }
831
832
833
834 template <int dim, int spacedim>
835 const std::vector<Tensor<1, spacedim>> &
837 {
838 return get_current_neighbor_fe_values().get_normal_vectors();
839 }
840
841
842
843 template <int dim, int spacedim>
844 const std::vector<types::global_dof_index> &
846 {
847 return local_dof_indices;
848 }
849
850
851
852 template <int dim, int spacedim>
853 unsigned int
855 {
856 return local_dof_indices.size();
857 }
858
859
860
861 template <int dim, int spacedim>
862 const std::vector<types::global_dof_index> &
864 {
865 return neighbor_dof_indices;
866 }
867
868
869
870 template <int dim, int spacedim>
871 unsigned int
873 {
874 return neighbor_dof_indices.size();
875 }
876
877
878
879 template <int dim, int spacedim>
882 {
883 return user_data_storage;
884 }
885
886
887
888 template <int dim, int spacedim>
889 const GeneralDataStorage &
891 {
892 return user_data_storage;
893 }
894
895
896
897 template <int dim, int spacedim>
900 {
901 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
902 Assert(mapping, ExcNotInitialized());
903 return *mapping;
904 }
905
906
907
908 template <int dim, int spacedim>
911 {
912 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
914 return *fe;
915 }
916
917
918
919 template <int dim, int spacedim>
920 const Quadrature<dim> &
922 {
923 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
924 return cell_quadrature;
925 }
926
927
928
929 template <int dim, int spacedim>
930 const Quadrature<dim - 1> &
932 {
933 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
934 return face_quadrature;
935 }
936
937
938
939 template <int dim, int spacedim>
942 {
943 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
944 Assert(mapping_collection, ExcNotInitialized());
945 return *mapping_collection;
946 }
947
948
949
950 template <int dim, int spacedim>
953 {
954 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
955 Assert(fe_collection, ExcNotInitialized());
956 return *fe_collection;
957 }
958
959
960
961 template <int dim, int spacedim>
964 {
965 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
966 return cell_quadrature_collection;
967 }
968
969
970
971 template <int dim, int spacedim>
972 const hp::QCollection<dim - 1> &
974 {
975 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
976 return face_quadrature_collection;
977 }
978
979
980
981 template <int dim, int spacedim>
982 bool
984 {
985 return hp_capability_enabled;
986 }
987
988
989
990 template <int dim, int spacedim>
993 {
994 return cell_update_flags;
995 }
996
997
998
999 template <int dim, int spacedim>
1002 {
1003 return neighbor_cell_update_flags;
1004 }
1005
1006
1007
1008 template <int dim, int spacedim>
1011 {
1012 return face_update_flags;
1013 }
1014
1015
1016
1017 template <int dim, int spacedim>
1020 {
1021 return neighbor_face_update_flags;
1022 }
1023
1024} // namespace MeshWorker
1026
1027// Explicit instantiations
1029namespace MeshWorker
1030{
1031#include "scratch_data.inst"
1032}
void reinit(const Triangulation< dim, spacedim > &tria)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
Abstract base class for mapping classes.
Definition mapping.h:318
bool has_hp_capabilities() const
const hp::MappingCollection< dim, spacedim > & get_mapping_collection() const
const FEValuesBase< dim, spacedim > & get_current_fe_values() const
GeneralDataStorage & get_general_data_storage()
const Quadrature< dim - 1 > & get_face_quadrature() const
const std::vector< Tensor< 1, spacedim > > & get_neighbor_normal_vectors()
unsigned int n_dofs_per_cell() const
unsigned int n_neighbor_dofs_per_cell() const
const std::vector< double > & get_JxW_values() const
const Mapping< dim, spacedim > & get_mapping() const
const std::vector< double > & get_neighbor_JxW_values() const
const FiniteElement< dim, spacedim > & get_fe() const
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
const hp::QCollection< dim > & get_cell_quadrature_collection() const
std::vector< types::global_dof_index > neighbor_dof_indices
const std::vector< types::global_dof_index > & get_neighbor_dof_indices() const
UpdateFlags get_face_update_flags() const
UpdateFlags get_neighbor_cell_update_flags() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors() const
const FEValues< dim, spacedim > & reinit_neighbor(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
std::vector< types::global_dof_index > local_dof_indices
const std::vector< Point< spacedim > > & get_quadrature_points() const
const FEValuesBase< dim, spacedim > & get_current_neighbor_fe_values() const
const std::vector< types::global_dof_index > & get_local_dof_indices() const
const FEInterfaceValues< dim, spacedim > & get_current_interface_fe_values() const
UpdateFlags get_cell_update_flags() const
const Quadrature< dim > & get_cell_quadrature() const
const hp::QCollection< dim - 1 > & get_face_quadrature_collection() const
ScratchData(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags &update_flags, const Quadrature< dim - 1 > &face_quadrature=Quadrature< dim - 1 >(), const UpdateFlags &face_update_flags=update_default)
const FEValues< dim, spacedim > & reinit(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell)
UpdateFlags get_neighbor_face_update_flags() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::active_cell_iterator active_cell_iterator
UpdateFlags
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
const types::fe_index invalid_fe_index
Definition types.h:243
static const unsigned int invalid_unsigned_int
Definition types.h:220