Reference documentation for deal.II version 9.4.1
\(\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// Copyright (C) 2019 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
17
18#include <memory>
19
21
22namespace MeshWorker
23{
24 template <int dim, int spacedim>
26 const Mapping<dim, spacedim> & mapping,
28 const Quadrature<dim> & quadrature,
29 const UpdateFlags & update_flags,
30 const Quadrature<dim - 1> & face_quadrature,
31 const UpdateFlags & face_update_flags)
32 : mapping(&mapping)
33 , fe(&fe)
34 , cell_quadrature(quadrature)
35 , face_quadrature(face_quadrature)
36 , hp_capability_enabled(false)
37 , cell_update_flags(update_flags)
38 , neighbor_cell_update_flags(update_flags)
39 , face_update_flags(face_update_flags)
40 , neighbor_face_update_flags(face_update_flags)
41 , local_dof_indices(fe.n_dofs_per_cell())
42 , neighbor_dof_indices(fe.n_dofs_per_cell())
43 {}
44
45
46
47 template <int dim, int spacedim>
49 const Mapping<dim, spacedim> & mapping,
51 const Quadrature<dim> & quadrature,
52 const UpdateFlags & update_flags,
53 const UpdateFlags & neighbor_update_flags,
54 const Quadrature<dim - 1> & face_quadrature,
55 const UpdateFlags & face_update_flags,
56 const UpdateFlags & neighbor_face_update_flags)
57 : mapping(&mapping)
58 , fe(&fe)
59 , cell_quadrature(quadrature)
60 , face_quadrature(face_quadrature)
61 , hp_capability_enabled(false)
62 , cell_update_flags(update_flags)
63 , neighbor_cell_update_flags(neighbor_update_flags)
64 , face_update_flags(face_update_flags)
65 , neighbor_face_update_flags(neighbor_face_update_flags)
66 , local_dof_indices(fe.n_dofs_per_cell())
67 , neighbor_dof_indices(fe.n_dofs_per_cell())
68 {}
69
70
71
72 template <int dim, int spacedim>
75 const Quadrature<dim> & quadrature,
76 const UpdateFlags & update_flags,
77 const Quadrature<dim - 1> & face_quadrature,
78 const UpdateFlags & face_update_flags)
79 : ScratchData(fe.reference_cell()
80 .template get_default_linear_mapping<dim, spacedim>(),
81 fe,
82 quadrature,
83 update_flags,
84 face_quadrature,
85 face_update_flags)
86 {}
87
88
89
90 template <int dim, int spacedim>
93 const Quadrature<dim> & quadrature,
94 const UpdateFlags & update_flags,
95 const UpdateFlags & neighbor_update_flags,
96 const Quadrature<dim - 1> & face_quadrature,
97 const UpdateFlags & face_update_flags,
98 const UpdateFlags & neighbor_face_update_flags)
99 : ScratchData(fe.reference_cell()
100 .template get_default_linear_mapping<dim, spacedim>(),
101 fe,
102 quadrature,
103 update_flags,
104 neighbor_update_flags,
105 face_quadrature,
106 face_update_flags,
107 neighbor_face_update_flags)
108 {}
109
110
111
112 template <int dim, int spacedim>
114 const hp::MappingCollection<dim, spacedim> &mapping_collection,
115 const hp::FECollection<dim, spacedim> & fe_collection,
116 const hp::QCollection<dim> & cell_quadrature_collection,
117 const UpdateFlags & cell_update_flags,
118 const hp::QCollection<dim - 1> & face_quadrature_collection,
119 const UpdateFlags & face_update_flags)
120 : mapping_collection(&mapping_collection)
121 , fe_collection(&fe_collection)
122 , cell_quadrature_collection(cell_quadrature_collection)
123 , face_quadrature_collection(face_quadrature_collection)
124 , hp_capability_enabled(true)
125 , cell_update_flags(cell_update_flags)
126 , neighbor_cell_update_flags(cell_update_flags)
127 , face_update_flags(face_update_flags)
128 , neighbor_face_update_flags(face_update_flags)
129 {
130 local_dof_indices.reserve(fe_collection.max_dofs_per_cell());
131 neighbor_dof_indices.reserve(fe_collection.max_dofs_per_cell());
132 }
133
134
135
136 template <int dim, int spacedim>
138 const hp::MappingCollection<dim, spacedim> &mapping_collection,
139 const hp::FECollection<dim, spacedim> & fe_collection,
140 const hp::QCollection<dim> & cell_quadrature_collection,
141 const UpdateFlags & cell_update_flags,
142 const UpdateFlags & neighbor_cell_update_flags,
143 const hp::QCollection<dim - 1> & face_quadrature_collection,
144 const UpdateFlags & face_update_flags,
145 const UpdateFlags & neighbor_face_update_flags)
146 : mapping_collection(&mapping_collection)
147 , fe_collection(&fe_collection)
148 , cell_quadrature_collection(cell_quadrature_collection)
149 , face_quadrature_collection(face_quadrature_collection)
150 , hp_capability_enabled(true)
151 , cell_update_flags(cell_update_flags)
152 , neighbor_cell_update_flags(neighbor_cell_update_flags)
153 , face_update_flags(face_update_flags)
154 , neighbor_face_update_flags(neighbor_face_update_flags)
155 {
156 local_dof_indices.reserve(fe_collection.max_dofs_per_cell());
157 neighbor_dof_indices.reserve(fe_collection.max_dofs_per_cell());
158 }
159
160
161
162 template <int dim, int spacedim>
164 const hp::FECollection<dim, spacedim> &fe_collection,
165 const hp::QCollection<dim> & cell_quadrature_collection,
166 const UpdateFlags & cell_update_flags,
167 const hp::QCollection<dim - 1> & face_quadrature_collection,
168 const UpdateFlags & face_update_flags)
169 : ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
170 fe_collection,
171 cell_quadrature_collection,
172 cell_update_flags,
173 face_quadrature_collection,
174 face_update_flags)
175 {}
176
177
178
179 template <int dim, int spacedim>
181 const hp::FECollection<dim, spacedim> &fe_collection,
182 const hp::QCollection<dim> & cell_quadrature_collection,
183 const UpdateFlags & cell_update_flags,
184 const UpdateFlags & neighbor_cell_update_flags,
185 const hp::QCollection<dim - 1> & face_quadrature_collection,
186 const UpdateFlags & face_update_flags,
187 const UpdateFlags & neighbor_face_update_flags)
188 : ScratchData(fe_collection.get_reference_cell_default_linear_mapping(),
189 fe_collection,
190 cell_quadrature_collection,
191 cell_update_flags,
192 neighbor_cell_update_flags,
193 face_quadrature_collection,
194 face_update_flags,
195 neighbor_face_update_flags)
196 {}
197
198
199
200 template <int dim, int spacedim>
202 const ScratchData<dim, spacedim> &scratch)
203 : mapping(scratch.mapping)
204 , fe(scratch.fe)
205 , cell_quadrature(scratch.cell_quadrature)
206 , face_quadrature(scratch.face_quadrature)
207 , mapping_collection(scratch.mapping_collection)
208 , fe_collection(scratch.fe_collection)
209 , cell_quadrature_collection(scratch.cell_quadrature_collection)
210 , face_quadrature_collection(scratch.face_quadrature_collection)
211 , hp_capability_enabled(scratch.hp_capability_enabled)
212 , cell_update_flags(scratch.cell_update_flags)
213 , neighbor_cell_update_flags(scratch.neighbor_cell_update_flags)
214 , face_update_flags(scratch.face_update_flags)
215 , neighbor_face_update_flags(scratch.neighbor_face_update_flags)
216 , local_dof_indices(scratch.local_dof_indices)
217 , neighbor_dof_indices(scratch.neighbor_dof_indices)
218 , user_data_storage(scratch.user_data_storage)
219 , internal_data_storage(scratch.internal_data_storage)
220 {}
221
222
223
224 template <int dim, int spacedim>
228 {
229 if (hp_capability_enabled == false)
230 {
231 if (!fe_values)
232 fe_values = std::make_unique<FEValues<dim, spacedim>>(
233 *mapping, *fe, cell_quadrature, cell_update_flags);
234
235 fe_values->reinit(cell);
236 local_dof_indices.resize(fe_values->dofs_per_cell);
237 cell->get_dof_indices(local_dof_indices);
238 current_fe_values = fe_values.get();
239 return *fe_values;
240 }
241 else
242 {
243 if (!hp_fe_values)
244 hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
245 *mapping_collection,
246 *fe_collection,
247 cell_quadrature_collection,
248 cell_update_flags);
249
250 hp_fe_values->reinit(cell);
251 const auto &fe_values = hp_fe_values->get_present_fe_values();
252
254 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
255 fe_values.dofs_per_cell);
256 local_dof_indices.resize(fe_values.dofs_per_cell);
257 cell->get_dof_indices(local_dof_indices);
258
259 current_fe_values = &fe_values;
260 return fe_values;
261 }
262 }
263
264
265
266 template <int dim, int spacedim>
270 const unsigned int face_no)
271 {
272 if (hp_capability_enabled == false)
273 {
274 if (!fe_face_values)
275 fe_face_values = std::make_unique<FEFaceValues<dim, spacedim>>(
276 *mapping, *fe, face_quadrature, face_update_flags);
277
278 fe_face_values->reinit(cell, face_no);
279 local_dof_indices.resize(fe->n_dofs_per_cell());
280 cell->get_dof_indices(local_dof_indices);
281 current_fe_values = fe_face_values.get();
282 return *fe_face_values;
283 }
284 else
285 {
286 return reinit(cell, cell, face_no);
287 }
288 }
289
290
291
292 template <int dim, int spacedim>
297 & neighbor_cell,
298 const unsigned int face_no)
299 {
300 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
301
302 if (!hp_fe_face_values)
303 hp_fe_face_values = std::make_unique<hp::FEFaceValues<dim, spacedim>>(
304 *mapping_collection,
305 *fe_collection,
306 face_quadrature_collection,
307 face_update_flags);
308
309 if (neighbor_cell == cell)
310 {
311 hp_fe_face_values->reinit(cell, face_no);
312 }
313 else
314 {
315 // When we want to ensure some agreement between the cell face and its
316 // neighbor on the quadrature order and mapping to use on this face,
317 // then we defer to the dominance of one FE over another. This should
318 // ensure that the optimal integration order and mapping order are
319 // selected for this situation.
320 const unsigned int dominated_fe_index =
321 fe_collection->find_dominated_fe(
322 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
323
324 hp_fe_face_values->reinit(cell,
325 face_no,
326 dominated_fe_index,
327 dominated_fe_index);
328 }
329
330 const auto &fe_face_values = hp_fe_face_values->get_present_fe_values();
331 const auto &fe = (*fe_collection)[cell->active_fe_index()];
332
333 local_dof_indices.resize(fe.n_dofs_per_cell());
334 cell->get_dof_indices(local_dof_indices);
335
336 current_fe_values = &fe_face_values;
337 return fe_face_values;
338 }
339
340
341
342 template <int dim, int spacedim>
346 const unsigned int face_no,
347 const unsigned int subface_no)
348 {
349 if (subface_no != numbers::invalid_unsigned_int)
350 {
351 if (hp_capability_enabled == false)
352 {
353 if (!fe_subface_values)
354 fe_subface_values =
355 std::make_unique<FESubfaceValues<dim, spacedim>>(
356 *mapping, *fe, face_quadrature, face_update_flags);
357 fe_subface_values->reinit(cell, face_no, subface_no);
358 local_dof_indices.resize(fe->n_dofs_per_cell());
359 cell->get_dof_indices(local_dof_indices);
360
361 current_fe_values = fe_subface_values.get();
362 return *fe_subface_values;
363 }
364 else
365 {
366 return reinit(cell, cell, face_no, subface_no);
367 }
368 }
369 else
370 return reinit(cell, face_no);
371 }
372
373
374
375 template <int dim, int spacedim>
380 & neighbor_cell,
381 const unsigned int face_no,
382 const unsigned int subface_no)
383 {
384 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
385
386 if (subface_no != numbers::invalid_unsigned_int)
387 {
388 if (!hp_fe_subface_values)
389 hp_fe_subface_values =
390 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
391 *mapping_collection,
392 *fe_collection,
393 face_quadrature_collection,
394 face_update_flags);
395
396 if (neighbor_cell == cell)
397 {
398 hp_fe_subface_values->reinit(cell, face_no, subface_no);
399 }
400 else
401 {
402 // When we want to ensure some agreement between the cell face and
403 // its neighbor on the quadrature order and mapping to use on this
404 // face, then we defer to the dominance of one FE over another. This
405 // should ensure that the optimal integration order and mapping
406 // order are selected for this situation.
407 const unsigned int dominated_fe_index =
408 fe_collection->find_dominated_fe(
409 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
410
411 hp_fe_subface_values->reinit(cell,
412 face_no,
413 subface_no,
414 dominated_fe_index,
415 dominated_fe_index);
416 }
417
418 const auto &fe_subface_values =
419 hp_fe_subface_values->get_present_fe_values();
420 const auto &fe = (*fe_collection)[cell->active_fe_index()];
421
422 local_dof_indices.resize(fe.n_dofs_per_cell());
423 cell->get_dof_indices(local_dof_indices);
424
425 current_fe_values = &fe_subface_values;
426 return fe_subface_values;
427 }
428 else
429 {
430 return reinit(cell, neighbor_cell, face_no);
431 }
432 }
433
434
435
436 template <int dim, int spacedim>
440 const unsigned int face_no,
441 const unsigned int sub_face_no,
443 & cell_neighbor,
444 const unsigned int face_no_neighbor,
445 const unsigned int sub_face_no_neighbor)
446 {
447 if (hp_capability_enabled == false)
448 {
449 if (!interface_fe_values)
450 interface_fe_values =
451 std::make_unique<FEInterfaceValues<dim, spacedim>>(
452 *mapping, *fe, face_quadrature, face_update_flags);
453 interface_fe_values->reinit(cell,
454 face_no,
455 sub_face_no,
456 cell_neighbor,
457 face_no_neighbor,
458 sub_face_no_neighbor);
459
460 current_fe_values = &interface_fe_values->get_fe_face_values(0);
461 current_neighbor_fe_values =
462 &interface_fe_values->get_fe_face_values(1);
463
464 cell_neighbor->get_dof_indices(neighbor_dof_indices);
465 local_dof_indices = interface_fe_values->get_interface_dof_indices();
466 return *interface_fe_values;
467 }
468 else
469 {
470 AssertThrow(false, ExcOnlyAvailableWithoutHP());
471 return *interface_fe_values;
472 }
473 }
474
475
476
477 template <int dim, int spacedim>
481 {
482 if (hp_capability_enabled == false)
483 {
484 if (!neighbor_fe_values)
485 neighbor_fe_values = std::make_unique<FEValues<dim, spacedim>>(
486 *mapping, *fe, cell_quadrature, neighbor_cell_update_flags);
487
488 neighbor_fe_values->reinit(cell);
489 cell->get_dof_indices(neighbor_dof_indices);
490 current_neighbor_fe_values = neighbor_fe_values.get();
491 return *neighbor_fe_values;
492 }
493 else
494 {
495 if (!neighbor_hp_fe_values)
496 neighbor_hp_fe_values = std::make_unique<hp::FEValues<dim, spacedim>>(
497 *mapping_collection,
498 *fe_collection,
499 cell_quadrature_collection,
500 neighbor_cell_update_flags);
501
502 neighbor_hp_fe_values->reinit(cell);
503 const auto &neighbor_fe_values =
504 neighbor_hp_fe_values->get_present_fe_values();
505
507 (*fe_collection)[cell->active_fe_index()].n_dofs_per_cell(),
508 neighbor_fe_values.dofs_per_cell);
509 neighbor_dof_indices.resize(neighbor_fe_values.dofs_per_cell);
510 cell->get_dof_indices(neighbor_dof_indices);
511
512 current_neighbor_fe_values = &neighbor_fe_values;
513 return neighbor_fe_values;
514 }
515 }
516
517
518
519 template <int dim, int spacedim>
523 const unsigned int face_no)
524 {
525 if (hp_capability_enabled == false)
526 {
527 if (!neighbor_fe_face_values)
528 neighbor_fe_face_values =
529 std::make_unique<FEFaceValues<dim, spacedim>>(
530 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
531 neighbor_fe_face_values->reinit(cell, face_no);
532 cell->get_dof_indices(neighbor_dof_indices);
533 current_neighbor_fe_values = neighbor_fe_face_values.get();
534 return *neighbor_fe_face_values;
535 }
536 else
537 {
538 return reinit_neighbor(cell, cell, face_no);
539 }
540 }
541
542
543
544 template <int dim, int spacedim>
549 & neighbor_cell,
550 const unsigned int face_no)
551 {
552 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
553
554 if (!neighbor_hp_fe_face_values)
555 neighbor_hp_fe_face_values =
556 std::make_unique<hp::FEFaceValues<dim, spacedim>>(
557 *mapping_collection,
558 *fe_collection,
559 face_quadrature_collection,
560 neighbor_face_update_flags);
561
562 if (neighbor_cell == cell)
563 {
564 neighbor_hp_fe_face_values->reinit(neighbor_cell, face_no);
565 }
566 else
567 {
568 // When we want to ensure some agreement between the cell face and its
569 // neighbor on the quadrature order and mapping to use on this face,
570 // then we defer to the dominance of one FE over another. This should
571 // ensure that the optimal integration order and mapping order are
572 // selected for this situation.
573 const unsigned int dominated_fe_index =
574 fe_collection->find_dominated_fe(
575 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
576
577 neighbor_hp_fe_face_values->reinit(neighbor_cell,
578 face_no,
579 dominated_fe_index,
580 dominated_fe_index);
581 }
582
583 const auto &neighbor_fe_face_values =
584 neighbor_hp_fe_face_values->get_present_fe_values();
585 const auto &neighbor_fe =
586 (*fe_collection)[neighbor_cell->active_fe_index()];
587
588 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
589 neighbor_cell->get_dof_indices(neighbor_dof_indices);
590
591 current_neighbor_fe_values = &neighbor_fe_face_values;
592 return neighbor_fe_face_values;
593 }
594
595
596
597 template <int dim, int spacedim>
601 const unsigned int face_no,
602 const unsigned int subface_no)
603 {
604 if (subface_no != numbers::invalid_unsigned_int)
605 {
606 if (hp_capability_enabled == false)
607 {
608 if (!neighbor_fe_subface_values)
609 neighbor_fe_subface_values =
610 std::make_unique<FESubfaceValues<dim, spacedim>>(
611 *mapping, *fe, face_quadrature, neighbor_face_update_flags);
612 neighbor_fe_subface_values->reinit(cell, face_no, subface_no);
613 cell->get_dof_indices(neighbor_dof_indices);
614 current_neighbor_fe_values = neighbor_fe_subface_values.get();
615 return *neighbor_fe_subface_values;
616 }
617 else
618 {
619 return reinit_neighbor(cell, cell, face_no, subface_no);
620 }
621 }
622 else
623 return reinit_neighbor(cell, face_no);
624 }
625
626
627
628 template <int dim, int spacedim>
633 & neighbor_cell,
634 const unsigned int face_no,
635 const unsigned int subface_no)
636 {
637 Assert(hp_capability_enabled, ExcOnlyAvailableWithHP());
638
639 if (subface_no != numbers::invalid_unsigned_int)
640 {
641 if (!neighbor_hp_fe_subface_values)
642 neighbor_hp_fe_subface_values =
643 std::make_unique<hp::FESubfaceValues<dim, spacedim>>(
644 *mapping_collection,
645 *fe_collection,
646 face_quadrature_collection,
647 neighbor_face_update_flags);
648
649 if (neighbor_cell == cell)
650 {
651 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
652 face_no,
653 subface_no);
654 }
655 else
656 {
657 // When we want to ensure some agreement between the cell face and
658 // its neighbor on the quadrature order and mapping to use on this
659 // face, then we defer to the dominance of one FE over another. This
660 // should ensure that the optimal integration order and mapping
661 // order are selected for this situation.
662 const unsigned int dominated_fe_index =
663 fe_collection->find_dominated_fe(
664 {cell->active_fe_index(), neighbor_cell->active_fe_index()});
665
666 neighbor_hp_fe_subface_values->reinit(neighbor_cell,
667 face_no,
668 subface_no,
669 dominated_fe_index,
670 dominated_fe_index);
671 }
672
673 const auto &neighbor_fe_subface_values =
674 neighbor_hp_fe_subface_values->get_present_fe_values();
675 const auto &neighbor_fe =
676 (*fe_collection)[neighbor_cell->active_fe_index()];
677
678 neighbor_dof_indices.resize(neighbor_fe.n_dofs_per_cell());
679 neighbor_cell->get_dof_indices(neighbor_dof_indices);
680
681 current_neighbor_fe_values = &neighbor_fe_subface_values;
682 return neighbor_fe_subface_values;
683 }
684 else
685 {
686 return reinit_neighbor(cell, neighbor_cell, face_no);
687 }
688 }
689
690
691
692 template <int dim, int spacedim>
695 {
696 Assert(current_fe_values != nullptr,
697 ExcMessage("You have to initialize the cache using one of the "
698 "reinit functions first!"));
699 return *current_fe_values;
700 }
701
702
703
704 template <int dim, int spacedim>
707 {
708 Assert(interface_fe_values != nullptr,
709 ExcMessage("You have to initialize the cache using one of the "
710 "reinit functions first!"));
711 return *interface_fe_values;
712 }
713
714
715
716 template <int dim, int spacedim>
719 {
720 Assert(current_neighbor_fe_values != nullptr,
721 ExcMessage("You have to initialize the cache using one of the "
722 "reinit functions first!"));
723 return *current_neighbor_fe_values;
724 }
725
726
727
728 template <int dim, int spacedim>
729 const std::vector<Point<spacedim>> &
731 {
732 return get_current_fe_values().get_quadrature_points();
733 }
734
735
736
737 template <int dim, int spacedim>
738 const std::vector<double> &
740 {
741 return get_current_fe_values().get_JxW_values();
742 }
743
744
745
746 template <int dim, int spacedim>
747 const std::vector<double> &
749 {
750 return get_current_neighbor_fe_values().get_JxW_values();
751 }
752
753
754
755 template <int dim, int spacedim>
756 const std::vector<Tensor<1, spacedim>> &
758 {
759 return get_current_fe_values().get_normal_vectors();
760 }
761
762
763
764 template <int dim, int spacedim>
765 const std::vector<Tensor<1, spacedim>> &
767 {
768 return get_current_neighbor_fe_values().get_normal_vectors();
769 }
770
771
772
773 template <int dim, int spacedim>
774 const std::vector<types::global_dof_index> &
776 {
777 return local_dof_indices;
778 }
779
780
781
782 template <int dim, int spacedim>
783 unsigned int
785 {
786 return local_dof_indices.size();
787 }
788
789
790
791 template <int dim, int spacedim>
792 const std::vector<types::global_dof_index> &
794 {
795 return neighbor_dof_indices;
796 }
797
798
799
800 template <int dim, int spacedim>
801 unsigned int
803 {
804 return neighbor_dof_indices.size();
805 }
806
807
808
809 template <int dim, int spacedim>
812 {
813 return user_data_storage;
814 }
815
816
817
818 template <int dim, int spacedim>
819 const GeneralDataStorage &
821 {
822 return user_data_storage;
823 }
824
825
826
827 template <int dim, int spacedim>
830 {
831 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
832 Assert(mapping, ExcNotInitialized());
833 return *mapping;
834 }
835
836
837
838 template <int dim, int spacedim>
841 {
842 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
844 return *fe;
845 }
846
847
848
849 template <int dim, int spacedim>
850 const Quadrature<dim> &
852 {
853 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
854 return cell_quadrature;
855 }
856
857
858
859 template <int dim, int spacedim>
860 const Quadrature<dim - 1> &
862 {
863 Assert(hp_capability_enabled == false, ExcOnlyAvailableWithoutHP());
864 return face_quadrature;
865 }
866
867
868
869 template <int dim, int spacedim>
872 {
873 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
874 Assert(mapping_collection, ExcNotInitialized());
875 return *mapping_collection;
876 }
877
878
879
880 template <int dim, int spacedim>
883 {
884 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
885 Assert(fe_collection, ExcNotInitialized());
886 return *fe_collection;
887 }
888
889
890
891 template <int dim, int spacedim>
894 {
895 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
896 return cell_quadrature_collection;
897 }
898
899
900
901 template <int dim, int spacedim>
902 const hp::QCollection<dim - 1> &
904 {
905 Assert(hp_capability_enabled == true, ExcOnlyAvailableWithHP());
906 return face_quadrature_collection;
907 }
908
909
910
911 template <int dim, int spacedim>
912 bool
914 {
915 return hp_capability_enabled;
916 }
917
918
919
920 template <int dim, int spacedim>
923 {
924 return cell_update_flags;
925 }
926
927
928
929 template <int dim, int spacedim>
932 {
933 return neighbor_cell_update_flags;
934 }
935
936
937
938 template <int dim, int spacedim>
941 {
942 return face_update_flags;
943 }
944
945
946
947 template <int dim, int spacedim>
950 {
951 return neighbor_face_update_flags;
952 }
953
954} // namespace MeshWorker
956
957// Explicit instantiations
959namespace MeshWorker
960{
961#include "scratch_data.inst"
962}
void reinit(const Triangulation< dim, spacedim > &tria)
Abstract base class for mapping classes.
Definition: mapping.h:311
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)
Definition: scratch_data.cc:25
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:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition: mapping.cc:260
static const unsigned int invalid_unsigned_int
Definition: types.h:201