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
fe_values.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 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
22
24
26
27#include <deal.II/fe/fe.h>
29#include <deal.II/fe/mapping.h>
30
33
34#include <deal.II/lac/vector.h>
35
36#include <boost/container/small_vector.hpp>
37
38#include <iomanip>
39#include <memory>
40#include <type_traits>
41
43
44
45namespace internal
46{
47 template <int dim, int spacedim>
48 inline std::vector<unsigned int>
50 {
51 std::vector<unsigned int> shape_function_to_row_table(
53 unsigned int row = 0;
54 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
55 {
56 // loop over all components that are nonzero for this particular
57 // shape function. if a component is zero then we leave the
58 // value in the table unchanged (at the invalid value)
59 // otherwise it is mapped to the next free entry
60 unsigned int nth_nonzero_component = 0;
61 for (unsigned int c = 0; c < fe.n_components(); ++c)
62 if (fe.get_nonzero_components(i)[c] == true)
63 {
64 shape_function_to_row_table[i * fe.n_components() + c] =
65 row + nth_nonzero_component;
66 ++nth_nonzero_component;
67 }
68 row += fe.n_nonzero_components(i);
69 }
70
71 return shape_function_to_row_table;
72 }
73} // namespace internal
74
75
76
77namespace internal
78{
79 namespace FEValuesImplementation
80 {
81 template <int dim, int spacedim>
82 void
84 const unsigned int n_quadrature_points,
86 const UpdateFlags flags)
87 {
88 // initialize the table mapping from shape function number to
89 // the rows in the tables storing the data by shape function and
90 // nonzero component
91 this->shape_function_to_row_table =
93
94 // count the total number of non-zero components accumulated
95 // over all shape functions
96 unsigned int n_nonzero_shape_components = 0;
97 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
98 n_nonzero_shape_components += fe.n_nonzero_components(i);
99 Assert(n_nonzero_shape_components >= fe.n_dofs_per_cell(),
101
102 // with the number of rows now known, initialize those fields
103 // that we will need to their correct size
104 if (flags & update_values)
105 {
106 this->shape_values.reinit(n_nonzero_shape_components,
107 n_quadrature_points);
108 this->shape_values.fill(numbers::signaling_nan<double>());
109 }
110
111 if (flags & update_gradients)
112 {
113 this->shape_gradients.reinit(n_nonzero_shape_components,
114 n_quadrature_points);
115 this->shape_gradients.fill(
117 }
118
119 if (flags & update_hessians)
120 {
121 this->shape_hessians.reinit(n_nonzero_shape_components,
122 n_quadrature_points);
123 this->shape_hessians.fill(
125 }
126
127 if (flags & update_3rd_derivatives)
128 {
129 this->shape_3rd_derivatives.reinit(n_nonzero_shape_components,
130 n_quadrature_points);
131 this->shape_3rd_derivatives.fill(
133 }
134 }
135
136
137
138 template <int dim, int spacedim>
139 std::size_t
141 {
142 return (
146 MemoryConsumption::memory_consumption(shape_3rd_derivatives) +
147 MemoryConsumption::memory_consumption(shape_function_to_row_table));
148 }
149 } // namespace FEValuesImplementation
150} // namespace internal
151
152/*------------------------------- FEValues -------------------------------*/
153#ifndef DOXYGEN
154
155template <int dim, int spacedim>
157
158
159
160template <int dim, int spacedim>
163 const Quadrature<dim> &q,
164 const UpdateFlags update_flags)
165 : FEValuesBase<dim, spacedim>(q.size(),
166 fe.n_dofs_per_cell(),
168 mapping,
169 fe)
170 , quadrature(q)
171{
172 initialize(update_flags);
173}
174
175
176
177template <int dim, int spacedim>
180 const hp::QCollection<dim> &q,
181 const UpdateFlags update_flags)
182 : FEValues(mapping, fe, q[0], update_flags)
183{
184 AssertDimension(q.size(), 1);
185}
186
187
188
189template <int dim, int spacedim>
191 const Quadrature<dim> &q,
192 const UpdateFlags update_flags)
193 : FEValuesBase<dim, spacedim>(
194 q.size(),
195 fe.n_dofs_per_cell(),
197 fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
198 fe)
199 , quadrature(q)
200{
201 initialize(update_flags);
202}
203
204
205
206template <int dim, int spacedim>
208 const hp::QCollection<dim> &q,
209 const UpdateFlags update_flags)
210 : FEValues(fe, q[0], update_flags)
211{
212 AssertDimension(q.size(), 1);
213}
214
215
216
217template <int dim, int spacedim>
218void
220{
221 // You can compute normal vectors to the cells only in the
222 // codimension one case.
223 if (dim != spacedim - 1)
224 Assert((update_flags & update_normal_vectors) == false,
225 ExcMessage("You can only pass the 'update_normal_vectors' "
226 "flag to FEFaceValues or FESubfaceValues objects, "
227 "but not to an FEValues object unless the "
228 "triangulation it refers to is embedded in a higher "
229 "dimensional space."));
230
231 const UpdateFlags flags = this->compute_update_flags(update_flags);
232
233 // initialize the base classes
234 if (flags & update_mapping)
235 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
236 this->finite_element_output.initialize(this->max_n_quadrature_points,
237 *this->fe,
238 flags);
239
240 // then get objects into which the FE and the Mapping can store
241 // intermediate data used across calls to reinit. we can do this in parallel
243 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
244 fe_get_data = Threads::new_task([&]() {
245 return this->fe->get_data(flags,
246 *this->mapping,
247 quadrature,
248 this->finite_element_output);
249 });
250
252 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
253 mapping_get_data;
254 if (flags & update_mapping)
255 mapping_get_data = Threads::new_task(
256 [&]() { return this->mapping->get_data(flags, quadrature); });
257
258 this->update_flags = flags;
259
260 // then collect answers from the two task above
261 this->fe_data = std::move(fe_get_data.return_value());
262 if (flags & update_mapping)
263 this->mapping_data = std::move(mapping_get_data.return_value());
264 else
265 this->mapping_data =
266 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
267}
268
269
270
271template <int dim, int spacedim>
272void
275{
276 // Check that mapping and reference cell type are compatible:
277 Assert(this->get_mapping().is_compatible_with(cell->reference_cell()),
279 "You are trying to call FEValues::reinit() with a cell of type " +
280 cell->reference_cell().to_string() +
281 " with a Mapping that is not compatible with it."));
282
283 // no FE in this cell, so no assertion
284 // necessary here
285 this->maybe_invalidate_previous_present_cell(cell);
286 this->check_cell_similarity(cell);
287
288 this->present_cell = {cell};
289
290 // this was the part of the work that is dependent on the actual
291 // data type of the iterator. now pass on to the function doing
292 // the real work.
293 do_reinit();
294}
295
296
297
298template <int dim, int spacedim>
299template <bool lda>
300void
303{
304 // assert that the finite elements passed to the constructor and
305 // used by the DoFHandler used by this cell, are the same
306 Assert(static_cast<const FiniteElementData<dim> &>(*this->fe) ==
307 static_cast<const FiniteElementData<dim> &>(cell->get_fe()),
309
310 // Check that mapping and reference cell type are compatible:
311 Assert(this->get_mapping().is_compatible_with(cell->reference_cell()),
313 "You are trying to call FEValues::reinit() with a cell of type " +
314 cell->reference_cell().to_string() +
315 " with a Mapping that is not compatible with it."));
316
317 this->maybe_invalidate_previous_present_cell(cell);
318 this->check_cell_similarity(cell);
319
320 this->present_cell = {cell};
321
322 // this was the part of the work that is dependent on the actual
323 // data type of the iterator. now pass on to the function doing
324 // the real work.
325 do_reinit();
326}
327
328
329
330template <int dim, int spacedim>
331void
333{
334 // first call the mapping and let it generate the data
335 // specific to the mapping. also let it inspect the
336 // cell similarity flag and, if necessary, update
337 // it
338 if (this->update_flags & update_mapping)
339 {
340 this->cell_similarity =
341 this->get_mapping().fill_fe_values(this->present_cell,
342 this->cell_similarity,
343 quadrature,
344 *this->mapping_data,
345 this->mapping_output);
346 }
347
348 // then call the finite element and, with the data
349 // already filled by the mapping, let it compute the
350 // data for the mapped shape function values, gradients,
351 // etc.
352 this->get_fe().fill_fe_values(this->present_cell,
353 this->cell_similarity,
354 this->quadrature,
355 this->get_mapping(),
356 *this->mapping_data,
357 this->mapping_output,
358 *this->fe_data,
359 this->finite_element_output);
360}
361
362
363
364template <int dim, int spacedim>
365std::size_t
367{
370}
371
372#endif
373/*------------------------------- FEFaceValuesBase --------------------------*/
374#ifndef DOXYGEN
375
376template <int dim, int spacedim>
378 const unsigned int dofs_per_cell,
379 const UpdateFlags flags,
380 const Mapping<dim, spacedim> &mapping,
382 const Quadrature<dim - 1> &quadrature)
383 : FEFaceValuesBase<dim, spacedim>(dofs_per_cell,
384 flags,
385 mapping,
386 fe,
387 hp::QCollection<dim - 1>(quadrature))
388{}
389
390
391
392template <int dim, int spacedim>
394 const unsigned int dofs_per_cell,
395 const UpdateFlags,
396 const Mapping<dim, spacedim> &mapping,
398 const hp::QCollection<dim - 1> &quadrature)
399 : FEValuesBase<dim, spacedim>(quadrature.max_n_quadrature_points(),
400 dofs_per_cell,
402 mapping,
403 fe)
404 , present_face_index(numbers::invalid_unsigned_int)
405 , quadrature(quadrature)
406{
407 Assert(quadrature.size() == 1 ||
408 quadrature.size() == fe.reference_cell().n_faces(),
410}
411
412
413
414template <int dim, int spacedim>
415const std::vector<Tensor<1, spacedim>> &
418 Assert(this->update_flags & update_boundary_forms,
420 "update_boundary_forms")));
421 return this->mapping_output.boundary_forms;
422}
423
424
425
426template <int dim, int spacedim>
427std::size_t
429{
432}
433#endif
434
435/*------------------------------- FEFaceValues -------------------------------*/
436
437#ifndef DOXYGEN
438
439template <int dim, int spacedim>
441
442
443
444template <int dim, int spacedim>
446
447
448
449template <int dim, int spacedim>
451 const Mapping<dim, spacedim> &mapping,
453 const Quadrature<dim - 1> &quadrature,
454 const UpdateFlags update_flags)
455 : FEFaceValues<dim, spacedim>(mapping,
456 fe,
457 hp::QCollection<dim - 1>(quadrature),
458 update_flags)
459{}
460
461
462
463template <int dim, int spacedim>
465 const Mapping<dim, spacedim> &mapping,
467 const hp::QCollection<dim - 1> &quadrature,
468 const UpdateFlags update_flags)
469 : FEFaceValuesBase<dim, spacedim>(fe.n_dofs_per_cell(),
470 update_flags,
471 mapping,
472 fe,
473 quadrature)
474{
475 initialize(update_flags);
476}
477
478
479
480template <int dim, int spacedim>
483 const Quadrature<dim - 1> &quadrature,
484 const UpdateFlags update_flags)
485 : FEFaceValues<dim, spacedim>(fe,
486 hp::QCollection<dim - 1>(quadrature),
487 update_flags)
488{}
489
490
491
492template <int dim, int spacedim>
495 const hp::QCollection<dim - 1> &quadrature,
496 const UpdateFlags update_flags)
497 : FEFaceValuesBase<dim, spacedim>(
498 fe.n_dofs_per_cell(),
499 update_flags,
500 fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
501 fe,
502 quadrature)
503{
504 initialize(update_flags);
505}
506
507
508
509template <int dim, int spacedim>
510void
512{
513 const UpdateFlags flags = this->compute_update_flags(update_flags);
514
515 // initialize the base classes
516 if (flags & update_mapping)
517 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
518 this->finite_element_output.initialize(this->max_n_quadrature_points,
519 *this->fe,
520 flags);
521
522 // then get objects into which the FE and the Mapping can store
523 // intermediate data used across calls to reinit. this can be done in parallel
524
525 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase> (
526 FiniteElement<dim, spacedim>::*finite_element_get_face_data)(
527 const UpdateFlags,
531 spacedim>
533
534 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> (
535 Mapping<dim, spacedim>::*mapping_get_face_data)(
536 const UpdateFlags, const hp::QCollection<dim - 1> &) const =
538
539
541 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
542 fe_get_data = Threads::new_task(finite_element_get_face_data,
543 *this->fe,
544 flags,
545 *this->mapping,
546 this->quadrature,
547 this->finite_element_output);
549 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
550 mapping_get_data;
551 if (flags & update_mapping)
552 mapping_get_data = Threads::new_task(mapping_get_face_data,
553 *this->mapping,
554 flags,
555 this->quadrature);
556
557 this->update_flags = flags;
558
559 // then collect answers from the two task above
560 this->fe_data = std::move(fe_get_data.return_value());
561 if (flags & update_mapping)
562 this->mapping_data = std::move(mapping_get_data.return_value());
563 else
564 this->mapping_data =
565 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
566}
567
568
569
570template <int dim, int spacedim>
571template <bool lda>
572void
575 const unsigned int face_no)
576{
577 // assert that the finite elements passed to the constructor and
578 // used by the DoFHandler used by this cell, are the same
579 Assert(static_cast<const FiniteElementData<dim> &>(*this->fe) ==
580 static_cast<const FiniteElementData<dim> &>(
581 cell->get_dof_handler().get_fe(cell->active_fe_index())),
583
585
586 this->maybe_invalidate_previous_present_cell(cell);
587 this->present_cell = {cell};
588
589 // this was the part of the work that is dependent on the actual
590 // data type of the iterator. now pass on to the function doing
591 // the real work.
592 do_reinit(face_no);
593}
594
595
596
597template <int dim, int spacedim>
598template <bool lda>
599void
603{
604 const auto face_n = cell->face_iterator_to_index(face);
605 reinit(cell, face_n);
606}
607
608
609
610template <int dim, int spacedim>
611void
614 const unsigned int face_no)
615{
617
618 this->maybe_invalidate_previous_present_cell(cell);
619 this->present_cell = {cell};
620
621 // this was the part of the work that is dependent on the actual
622 // data type of the iterator. now pass on to the function doing
623 // the real work.
624 do_reinit(face_no);
625}
626
627
628
629template <int dim, int spacedim>
630void
634{
635 const auto face_n = cell->face_iterator_to_index(face);
636 reinit(cell, face_n);
637}
638
639
640
641template <int dim, int spacedim>
642void
643FEFaceValues<dim, spacedim>::do_reinit(const unsigned int face_no)
644{
645 this->present_face_no = face_no;
646
647 // first of all, set the present_face_index (if available)
649 this->present_cell;
650 this->present_face_index = cell->face_index(face_no);
651
652 if (this->update_flags & update_mapping)
653 {
654 this->get_mapping().fill_fe_face_values(this->present_cell,
655 face_no,
656 this->quadrature,
657 *this->mapping_data,
658 this->mapping_output);
659 }
660
661 this->get_fe().fill_fe_face_values(this->present_cell,
662 face_no,
663 this->quadrature,
664 this->get_mapping(),
665 *this->mapping_data,
666 this->mapping_output,
667 *this->fe_data,
668 this->finite_element_output);
669
670 const_cast<unsigned int &>(this->n_quadrature_points) =
671 this->quadrature[this->quadrature.size() == 1 ? 0 : face_no].size();
672}
673
674
675/* ---------------------------- FESubFaceValues ---------------------------- */
676
677
678template <int dim, int spacedim>
680
681
682
683template <int dim, int spacedim>
685
686
687
688template <int dim, int spacedim>
690 const Mapping<dim, spacedim> &mapping,
692 const Quadrature<dim - 1> &quadrature,
693 const UpdateFlags update_flags)
694 : FEFaceValuesBase<dim, spacedim>(fe.n_dofs_per_cell(),
695 update_flags,
696 mapping,
697 fe,
698 quadrature)
699{
700 initialize(update_flags);
701}
702
703
704
705template <int dim, int spacedim>
707 const Mapping<dim, spacedim> &mapping,
709 const hp::QCollection<dim - 1> &quadrature,
710 const UpdateFlags update_flags)
711 : FESubfaceValues(mapping, fe, quadrature[0], update_flags)
712{
713 AssertDimension(quadrature.size(), 1);
714}
715
716
717
718template <int dim, int spacedim>
721 const Quadrature<dim - 1> &quadrature,
722 const UpdateFlags update_flags)
723 : FEFaceValuesBase<dim, spacedim>(
724 fe.n_dofs_per_cell(),
725 update_flags,
726 fe.reference_cell().template get_default_linear_mapping<dim, spacedim>(),
727 fe,
728 quadrature)
729{
730 initialize(update_flags);
731}
732
733
734
735template <int dim, int spacedim>
738 const hp::QCollection<dim - 1> &quadrature,
739 const UpdateFlags update_flags)
740 : FESubfaceValues(fe, quadrature[0], update_flags)
741{
742 AssertDimension(quadrature.size(), 1);
743}
744
745
746
747template <int dim, int spacedim>
748void
750{
751 const UpdateFlags flags = this->compute_update_flags(update_flags);
752
753 // initialize the base classes
754 if (flags & update_mapping)
755 this->mapping_output.initialize(this->max_n_quadrature_points, flags);
756 this->finite_element_output.initialize(this->max_n_quadrature_points,
757 *this->fe,
758 flags);
759
760 // then get objects into which the FE and the Mapping can store
761 // intermediate data used across calls to reinit. this can be done
762 // in parallel
764 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>>
765 fe_get_data =
767 *this->fe,
768 flags,
769 *this->mapping,
770 this->quadrature[0],
771 this->finite_element_output);
773 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
774 mapping_get_data;
775 if (flags & update_mapping)
776 mapping_get_data =
778 *this->mapping,
779 flags,
780 this->quadrature[0]);
781
782 this->update_flags = flags;
783
784 // then collect answers from the two task above
785 this->fe_data = std::move(fe_get_data.return_value());
786 if (flags & update_mapping)
787 this->mapping_data = std::move(mapping_get_data.return_value());
788 else
789 this->mapping_data =
790 std::make_unique<typename Mapping<dim, spacedim>::InternalDataBase>();
791}
792
793
794
795template <int dim, int spacedim>
796template <bool lda>
797void
800 const unsigned int face_no,
801 const unsigned int subface_no)
802{
803 // assert that the finite elements passed to the constructor and
804 // used by the DoFHandler used by this cell, are the same
805 Assert(static_cast<const FiniteElementData<dim> &>(*this->fe) ==
806 static_cast<const FiniteElementData<dim> &>(
807 cell->get_dof_handler().get_fe(cell->active_fe_index())),
810 // We would like to check for subface_no < cell->face(face_no)->n_children(),
811 // but unfortunately the current function is also called for
812 // faces without children (see tests/fe/mapping.cc). Therefore,
813 // we must use following workaround of two separate assertions
814 Assert(cell->face(face_no)->has_children() ||
816 ExcIndexRange(subface_no,
817 0,
819 Assert(!cell->face(face_no)->has_children() ||
820 subface_no < cell->face(face_no)->n_active_descendants(),
821 ExcIndexRange(subface_no,
822 0,
823 cell->face(face_no)->n_active_descendants()));
824 Assert(cell->has_children() == false,
825 ExcMessage("You can't use subface data for cells that are "
826 "already refined. Iterate over their children "
827 "instead in these cases."));
828
829 this->maybe_invalidate_previous_present_cell(cell);
830 this->present_cell = {cell};
831
832 // this was the part of the work that is dependent on the actual
833 // data type of the iterator. now pass on to the function doing
834 // the real work.
835 do_reinit(face_no, subface_no);
836}
837
838
839
840template <int dim, int spacedim>
841template <bool lda>
842void
846 const typename Triangulation<dim, spacedim>::face_iterator &subface)
847{
848 reinit(cell,
849 cell->face_iterator_to_index(face),
850 face->child_iterator_to_index(subface));
851}
852
853
854
855template <int dim, int spacedim>
856void
859 const unsigned int face_no,
860 const unsigned int subface_no)
861{
863 // We would like to check for subface_no < cell->face(face_no)->n_children(),
864 // but unfortunately the current function is also called for
865 // faces without children for periodic faces, which have hanging nodes on
866 // the other side (see include/deal.II/matrix_free/mapping_info.templates.h).
867 AssertIndexRange(subface_no,
868 (cell->has_periodic_neighbor(face_no) ?
869 cell->periodic_neighbor(face_no)
870 ->face(cell->periodic_neighbor_face_no(face_no))
871 ->n_children() :
872 cell->face(face_no)->n_children()));
873
874 this->maybe_invalidate_previous_present_cell(cell);
875 this->present_cell = {cell};
876
877 // this was the part of the work that is dependent on the actual
878 // data type of the iterator. now pass on to the function doing
879 // the real work.
880 do_reinit(face_no, subface_no);
881}
882
883
884
885template <int dim, int spacedim>
886void
890 const typename Triangulation<dim, spacedim>::face_iterator &subface)
891{
892 reinit(cell,
893 cell->face_iterator_to_index(face),
894 face->child_iterator_to_index(subface));
895}
896
897
898
899template <int dim, int spacedim>
900void
901FESubfaceValues<dim, spacedim>::do_reinit(const unsigned int face_no,
902 const unsigned int subface_no)
903{
904 this->present_face_no = face_no;
905
906 // first of all, set the present_face_index (if available)
908 this->present_cell;
909
910 if (!cell->face(face_no)->has_children())
911 // no subfaces at all, so set present_face_index to this face rather
912 // than any subface
913 this->present_face_index = cell->face_index(face_no);
914 else if (dim != 3)
915 this->present_face_index = cell->face(face_no)->child_index(subface_no);
916 else
917 {
918 // this is the same logic we use in cell->neighbor_child_on_subface(). See
919 // there for an explanation of the different cases
920 unsigned int subface_index = numbers::invalid_unsigned_int;
921 switch (cell->subface_case(face_no))
922 {
926 subface_index = cell->face(face_no)->child_index(subface_no);
927 break;
930 subface_index = cell->face(face_no)
931 ->child(subface_no / 2)
932 ->child_index(subface_no % 2);
933 break;
936 switch (subface_no)
937 {
938 case 0:
939 case 1:
940 subface_index =
941 cell->face(face_no)->child(0)->child_index(subface_no);
942 break;
943 case 2:
944 subface_index = cell->face(face_no)->child_index(1);
945 break;
946 default:
948 }
949 break;
952 switch (subface_no)
953 {
954 case 0:
955 subface_index = cell->face(face_no)->child_index(0);
956 break;
957 case 1:
958 case 2:
959 subface_index =
960 cell->face(face_no)->child(1)->child_index(subface_no - 1);
961 break;
962 default:
964 }
965 break;
966 default:
968 break;
969 }
970 Assert(subface_index != numbers::invalid_unsigned_int,
972 this->present_face_index = subface_index;
973 }
974
975 // now ask the mapping and the finite element to do the actual work
976 if (this->update_flags & update_mapping)
977 {
978 this->get_mapping().fill_fe_subface_values(this->present_cell,
979 face_no,
980 subface_no,
981 this->quadrature[0],
982 *this->mapping_data,
983 this->mapping_output);
984 }
985
986 this->get_fe().fill_fe_subface_values(this->present_cell,
987 face_no,
988 subface_no,
989 this->quadrature[0],
990 this->get_mapping(),
991 *this->mapping_data,
992 this->mapping_output,
993 *this->fe_data,
994 this->finite_element_output);
995}
996#endif
997
998/*------------------------- Explicit Instantiations --------------------------*/
999
1000#include "fe_values.inst"
1001
FEFaceValuesBase(const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature)
std::size_t memory_consumption() const
const std::vector< Tensor< 1, spacedim > > & get_boundary_forms() const
void initialize(const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no)
void do_reinit(const unsigned int face_no)
FEFaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
FESubfaceValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &face_quadrature, const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell, const unsigned int face_no, const unsigned int subface_no)
void initialize(const UpdateFlags update_flags)
void do_reinit(const unsigned int face_no, const unsigned int subface_no)
FEValues(const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
void do_reinit()
void initialize(const UpdateFlags update_flags)
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
std::size_t memory_consumption() const
unsigned int n_dofs_per_cell() const
unsigned int n_components() const
ReferenceCell reference_cell() const
const ComponentMask & get_nonzero_components(const unsigned int i) const
unsigned int n_nonzero_components(const unsigned int i) const
Abstract base class for mapping classes.
Definition mapping.h:318
unsigned int n_faces() const
unsigned int size() const
Definition collection.h:308
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition fe_values.cc:83
#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)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
UpdateFlags
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_normal_vectors
Normal vectors.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_mapping
@ update_gradients
Shape function gradients.
@ update_default
No update.
@ update_boundary_forms
Outer normal vector, not normalized.
Task< RT > new_task(const std::function< RT()> &function)
#define DEAL_II_ASSERT_UNREACHABLE()
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
Definition mapping.cc:294
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Definition hp.h:117
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
std::vector< unsigned int > make_shape_function_to_row_table(const FiniteElement< dim, spacedim > &fe)
Definition fe_values.cc:49
static const unsigned int invalid_unsigned_int
Definition types.h:220
T signaling_nan()