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.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
18
20
21#include <deal.II/fe/fe.h>
23#include <deal.II/fe/mapping.h>
24
25#include <deal.II/grid/tria.h>
28
29#include <algorithm>
30#include <functional>
31#include <numeric>
32#include <typeinfo>
33
35
36
37/*------------------------------- FiniteElement ----------------------*/
38#ifndef DOXYGEN
39
40template <int dim, int spacedim>
42 : update_each(update_default)
43{}
44
45
46
47template <int dim, int spacedim>
48std::size_t
50{
51 return sizeof(*this);
52}
53
54
55
56template <int dim, int spacedim>
58 const FiniteElementData<dim> &fe_data,
59 const std::vector<bool> &r_i_a_f,
60 const std::vector<ComponentMask> &nonzero_c)
61 : FiniteElementData<dim>(fe_data)
65 std::make_pair(std::make_pair(0U, 0U), 0U))
66 ,
67
68 // Special handling of vectors of length one: in this case, we
69 // assume that all entries were supposed to be equal
71 r_i_a_f.size() == 1 ?
72 std::vector<bool>(fe_data.n_dofs_per_cell(), r_i_a_f[0]) :
73 r_i_a_f)
75 nonzero_c.size() == 1 ?
76 std::vector<ComponentMask>(fe_data.n_dofs_per_cell(), nonzero_c[0]) :
77 nonzero_c)
81 [](const unsigned int n_components) {
82 return n_components != 1U;
83 }) == n_nonzero_components_table.end())
84{
85 Assert(restriction_is_additive_flags.size() == this->n_dofs_per_cell(),
86 ExcDimensionMismatch(restriction_is_additive_flags.size(),
87 this->n_dofs_per_cell()));
88 AssertDimension(nonzero_components.size(), this->n_dofs_per_cell());
89 for (unsigned int i = 0; i < nonzero_components.size(); ++i)
90 {
91 Assert(nonzero_components[i].size() == this->n_components(),
93 Assert(nonzero_components[i].n_selected_components() >= 1,
95 Assert(n_nonzero_components_table[i] >= 1, ExcInternalError());
96 Assert(n_nonzero_components_table[i] <= this->n_components(),
98 }
99
100 // initialize some tables in the default way, i.e. if there is only one
101 // (vector-)component; if the element is not primitive, leave these tables
102 // empty.
103 if (this->is_primitive())
104 {
105 system_to_component_table.resize(this->n_dofs_per_cell());
106 for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
107 system_to_component_table[j] = std::pair<unsigned, unsigned>(0, j);
108
109 face_system_to_component_table.resize(this->n_unique_faces());
110 for (unsigned int f = 0; f < this->n_unique_faces(); ++f)
111 {
112 face_system_to_component_table[f].resize(this->n_dofs_per_face(f));
113 for (unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
114 face_system_to_component_table[f][j] =
115 std::pair<unsigned, unsigned>(0, j);
116 }
117 }
118
119 for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
120 system_to_base_table[j] = std::make_pair(std::make_pair(0U, 0U), j);
121
122 face_system_to_base_table.resize(this->n_unique_faces());
123 for (unsigned int f = 0; f < this->n_unique_faces(); ++f)
124 {
125 face_system_to_base_table[f].resize(this->n_dofs_per_face(f));
126 for (unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
127 face_system_to_base_table[f][j] =
128 std::make_pair(std::make_pair(0U, 0U), j);
129 }
130
131 // Fill with default value; may be changed by constructor of derived class.
132 base_to_block_indices.reinit(1, 1);
133
134 // initialize the restriction and prolongation matrices. the default
135 // constructor of FullMatrix<dim> initializes them with size zero
138 for (const unsigned int ref_case :
139 RefinementCase<dim>::all_refinement_cases())
140 if (ref_case != RefinementCase<dim>::no_refinement)
141 {
142 prolongation[ref_case - 1].resize(this->reference_cell().n_children(
143 RefinementCase<dim>(ref_case)),
145 restriction[ref_case - 1].resize(this->reference_cell().n_children(
146 RefinementCase<dim>(ref_case)),
148 }
149
150
151 if (dim == 3)
152 {
153 adjust_quad_dof_index_for_face_orientation_table.resize(
154 this->n_unique_2d_subobjects());
155
156 for (unsigned int f = 0; f < this->n_unique_2d_subobjects(); ++f)
157 {
158 adjust_quad_dof_index_for_face_orientation_table[f] =
159 Table<2, int>(this->n_dofs_per_quad(f),
160 this->reference_cell().n_face_orientations(f));
161 adjust_quad_dof_index_for_face_orientation_table[f].fill(0);
162 }
163 }
164
165 unit_face_support_points.resize(this->n_unique_faces());
166 generalized_face_support_points.resize(this->n_unique_faces());
167}
168
169
170
171template <int dim, int spacedim>
172std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>
173FiniteElement<dim, spacedim>::operator^(const unsigned int multiplicity) const
174{
175 return {this->clone(), multiplicity};
176}
177
178
179
180template <int dim, int spacedim>
181double
183 const Point<dim> &) const
184{
185 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
186 return 0.;
187}
188
189
190
191template <int dim, int spacedim>
192double
194 const Point<dim> &,
195 const unsigned int) const
196{
197 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
198 return 0.;
199}
200
201
202
203template <int dim, int spacedim>
206 const Point<dim> &) const
207{
208 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
209 return Tensor<1, dim>();
210}
211
212
213
214template <int dim, int spacedim>
217 const Point<dim> &,
218 const unsigned int) const
219{
220 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
221 return Tensor<1, dim>();
222}
223
224
225
226template <int dim, int spacedim>
229 const Point<dim> &) const
230{
231 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
232 return Tensor<2, dim>();
233}
234
235
236
237template <int dim, int spacedim>
240 const unsigned int,
241 const Point<dim> &,
242 const unsigned int) const
243{
244 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
245 return Tensor<2, dim>();
246}
247
248
249
250template <int dim, int spacedim>
253 const Point<dim> &) const
254{
255 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
256 return Tensor<3, dim>();
257}
258
259
260
261template <int dim, int spacedim>
264 const unsigned int,
265 const Point<dim> &,
266 const unsigned int) const
267{
268 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
269 return Tensor<3, dim>();
270}
271
272
273
274template <int dim, int spacedim>
277 const Point<dim> &) const
278{
279 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
280 return Tensor<4, dim>();
281}
282
283
284
285template <int dim, int spacedim>
288 const unsigned int,
289 const Point<dim> &,
290 const unsigned int) const
291{
292 AssertThrow(false, ExcUnitShapeValuesDoNotExist());
293 return Tensor<4, dim>();
294}
295
296
297template <int dim, int spacedim>
298void
300 const bool isotropic_restriction_only,
301 const bool isotropic_prolongation_only)
302{
303 for (const unsigned int ref_case :
304 RefinementCase<dim>::all_refinement_cases())
305 if (ref_case != RefinementCase<dim>::no_refinement)
306 {
307 const unsigned int nc =
308 this->reference_cell().n_children(RefinementCase<dim>(ref_case));
309
310 for (unsigned int i = 0; i < nc; ++i)
311 {
312 if (this->restriction[ref_case - 1][i].m() !=
313 this->n_dofs_per_cell() &&
314 (!isotropic_restriction_only ||
316 this->restriction[ref_case - 1][i].reinit(
317 this->n_dofs_per_cell(), this->n_dofs_per_cell());
318 if (this->prolongation[ref_case - 1][i].m() !=
319 this->n_dofs_per_cell() &&
320 (!isotropic_prolongation_only ||
322 this->prolongation[ref_case - 1][i].reinit(
323 this->n_dofs_per_cell(), this->n_dofs_per_cell());
324 }
325 }
326}
327
328
329template <int dim, int spacedim>
330const FullMatrix<double> &
332 const unsigned int child,
333 const RefinementCase<dim> &refinement_case) const
334{
335 AssertIndexRange(refinement_case,
339 "Restriction matrices are only available for refined cells!"));
340 AssertIndexRange(child,
341 this->reference_cell().n_children(
342 RefinementCase<dim>(refinement_case)));
343 // we use refinement_case-1 here. the -1 takes care of the origin of the
344 // vector, as for RefinementCase<dim>::no_refinement (=0) there is no data
345 // available and so the vector indices are shifted
346 Assert(restriction[refinement_case - 1][child].n() == this->n_dofs_per_cell(),
347 ExcProjectionVoid());
348 return restriction[refinement_case - 1][child];
349}
350
351
352
353template <int dim, int spacedim>
354const FullMatrix<double> &
356 const unsigned int child,
357 const RefinementCase<dim> &refinement_case) const
358{
359 AssertIndexRange(refinement_case,
363 "Prolongation matrices are only available for refined cells!"));
364 AssertIndexRange(child,
365 this->reference_cell().n_children(
366 RefinementCase<dim>(refinement_case)));
367 // we use refinement_case-1 here. the -1 takes care
368 // of the origin of the vector, as for
369 // RefinementCase::no_refinement (=0) there is no
370 // data available and so the vector indices
371 // are shifted
372 Assert(prolongation[refinement_case - 1][child].n() ==
373 this->n_dofs_per_cell(),
374 ExcEmbeddingVoid());
375 return prolongation[refinement_case - 1][child];
376}
377
378
379// TODO:[GK] This is probably not the most efficient way of doing this.
380template <int dim, int spacedim>
381unsigned int
383 const unsigned int index) const
384{
385 AssertIndexRange(index, this->n_components());
386
387 return first_block_of_base(component_to_base_table[index].first.first) +
388 component_to_base_table[index].second;
389}
390
391
392template <int dim, int spacedim>
395 const FEValuesExtractors::Scalar &scalar) const
396{
397 AssertIndexRange(scalar.component, this->n_components());
398
399 // TODO: it would be nice to verify that it is indeed possible
400 // to select this scalar component, i.e., that it is not part
401 // of a non-primitive element. unfortunately, there is no simple
402 // way to write such a condition...
403
404 std::vector<bool> mask(this->n_components(), false);
405 mask[scalar.component] = true;
406 return ComponentMask(mask);
407}
408
409
410template <int dim, int spacedim>
413 const FEValuesExtractors::Vector &vector) const
414{
416 this->n_components());
417
418 // TODO: it would be nice to verify that it is indeed possible
419 // to select these vector components, i.e., that they don't span
420 // beyond the beginning or end of anon-primitive element.
421 // unfortunately, there is no simple way to write such a condition...
422
423 std::vector<bool> mask(this->n_components(), false);
424 for (unsigned int c = vector.first_vector_component;
425 c < vector.first_vector_component + dim;
426 ++c)
427 mask[c] = true;
428 return ComponentMask(mask);
429}
430
431
432template <int dim, int spacedim>
435 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const
436{
439 this->n_components());
440
441 // TODO: it would be nice to verify that it is indeed possible
442 // to select these vector components, i.e., that they don't span
443 // beyond the beginning or end of anon-primitive element.
444 // unfortunately, there is no simple way to write such a condition...
445
446 std::vector<bool> mask(this->n_components(), false);
447 for (unsigned int c = sym_tensor.first_tensor_component;
448 c < sym_tensor.first_tensor_component +
450 ++c)
451 mask[c] = true;
452 return ComponentMask(mask);
453}
454
455
456
457template <int dim, int spacedim>
460{
461 // if we get a block mask that represents all blocks, then
462 // do the same for the returned component mask
463 if (block_mask.represents_the_all_selected_mask())
464 return {};
465
466 AssertDimension(block_mask.size(), this->n_blocks());
467
468 std::vector<bool> component_mask(this->n_components(), false);
469 for (unsigned int c = 0; c < this->n_components(); ++c)
470 if (block_mask[component_to_block_index(c)] == true)
471 component_mask[c] = true;
472
473 return ComponentMask(component_mask);
474}
475
476
477
478template <int dim, int spacedim>
481 const FEValuesExtractors::Scalar &scalar) const
482{
483 // simply create the corresponding component mask (a simpler
484 // process) and then convert it to a block mask
485 return block_mask(component_mask(scalar));
486}
487
488
489template <int dim, int spacedim>
492 const FEValuesExtractors::Vector &vector) const
493{
494 // simply create the corresponding component mask (a simpler
495 // process) and then convert it to a block mask
496 return block_mask(component_mask(vector));
497}
498
499
500template <int dim, int spacedim>
503 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const
504{
505 // simply create the corresponding component mask (a simpler
506 // process) and then convert it to a block mask
507 return block_mask(component_mask(sym_tensor));
508}
509
510
511
512template <int dim, int spacedim>
515 const ComponentMask &component_mask) const
516{
517 // if we get a component mask that represents all component, then
518 // do the same for the returned block mask
519 if (component_mask.represents_the_all_selected_mask())
520 return {};
521
522 AssertDimension(component_mask.size(), this->n_components());
523
524 // walk over all of the components
525 // of this finite element and see
526 // if we need to set the
527 // corresponding block. inside the
528 // block, walk over all the
529 // components that correspond to
530 // this block and make sure the
531 // component mask is set for all of
532 // them
533 std::vector<bool> block_mask(this->n_blocks(), false);
534 for (unsigned int c = 0; c < this->n_components();)
535 {
536 const unsigned int block = component_to_block_index(c);
537 if (component_mask[c] == true)
538 block_mask[block] = true;
539
540 // now check all of the other
541 // components that correspond
542 // to this block
543 ++c;
544 while ((c < this->n_components()) &&
545 (component_to_block_index(c) == block))
546 {
547 Assert(component_mask[c] == block_mask[block],
549 "The component mask argument given to this function "
550 "is not a mask where the individual components belonging "
551 "to one block of the finite element are either all "
552 "selected or not selected. You can't call this function "
553 "with a component mask that splits blocks."));
554 ++c;
555 }
556 }
557
558
559 return BlockMask(block_mask);
560}
561
562
563
564template <int dim, int spacedim>
565unsigned int
567 const unsigned int face_index,
568 const unsigned int face,
569 const unsigned char combined_orientation) const
570{
571 AssertIndexRange(face_index, this->n_dofs_per_face(face));
572 AssertIndexRange(face, this->reference_cell().n_faces());
573
574 // see the function's documentation for an explanation of this
575 // assertion -- in essence, derived classes have to implement
576 // an overloaded version of this function if we are to use any
577 // other than default (standard) orientation
578 if (combined_orientation !=
580 Assert((this->n_dofs_per_line() <= 1) && (this->n_dofs_per_quad(face) <= 1),
582 "The function in this base class can not handle this case. "
583 "Rather, the derived class you are using must provide "
584 "an overloaded version but apparently hasn't done so. See "
585 "the documentation of this function for more information."));
586
587 // we need to distinguish between DoFs on vertices, lines and (in 3d) quads.
588 // do so in a sequence of if-else statements
589 if (face_index < this->get_first_face_line_index(face))
590 // DoF is on a vertex
591 {
592 // get the number of the vertex on the face that corresponds to this DoF,
593 // along with the number of the DoF on this vertex
594 const unsigned int face_vertex = face_index / this->n_dofs_per_vertex();
595 const unsigned int dof_index_on_vertex =
596 face_index % this->n_dofs_per_vertex();
597
598 // then get the number of this vertex on the cell and translate
599 // this to a DoF number on the cell
600 return (this->reference_cell().face_to_cell_vertices(
601 face, face_vertex, combined_orientation) *
602 this->n_dofs_per_vertex() +
603 dof_index_on_vertex);
604 }
605 else if (face_index < this->get_first_face_quad_index(face))
606 // DoF is on a face
607 {
608 // do the same kind of translation as before. we need to only consider
609 // DoFs on the lines, i.e., ignoring those on the vertices
610 const unsigned int index =
611 face_index - this->get_first_face_line_index(face);
612
613 const unsigned int face_line = index / this->n_dofs_per_line();
614 const unsigned int dof_index_on_line = index % this->n_dofs_per_line();
615
616 return (this->get_first_line_index() +
617 this->reference_cell().face_to_cell_lines(face,
618 face_line,
619 combined_orientation) *
620 this->n_dofs_per_line() +
621 dof_index_on_line);
622 }
623 else
624 // DoF is on a quad
625 {
626 Assert(dim >= 3, ExcInternalError());
627
628 // ignore vertex and line dofs
629 const unsigned int index =
630 face_index - this->get_first_face_quad_index(face);
631
632 return (this->get_first_quad_index(face) + index);
633 }
634}
635
636
637
638template <int dim, int spacedim>
639unsigned int
641 const unsigned int index,
642 const unsigned int face,
643 const unsigned char combined_orientation) const
644{
645 // general template for 1d and 2d: not
646 // implemented. in fact, the function
647 // shouldn't even be called unless we are
648 // in 3d, so throw an internal error
649 Assert(dim == 3, ExcInternalError());
650 if (dim < 3)
651 return index;
652
653 // adjust dofs on 3d faces if the face is
654 // flipped. note that we query a table that
655 // derived elements need to have set up
656 // front. the exception are discontinuous
657 // elements for which there should be no
658 // face dofs anyway (i.e. dofs_per_quad==0
659 // in 3d), so we don't need the table, but
660 // the function should also not have been
661 // called
662 AssertIndexRange(index, this->n_dofs_per_quad(face));
663 const auto table_n = this->n_unique_2d_subobjects() == 1 ? 0 : face;
664 Assert(
665 adjust_quad_dof_index_for_face_orientation_table[table_n].n_elements() ==
666 (this->reference_cell().n_face_orientations(face)) *
667 this->n_dofs_per_quad(face),
669 return index + adjust_quad_dof_index_for_face_orientation_table[table_n](
670 index, combined_orientation);
671}
672
673
674
675template <int dim, int spacedim>
676unsigned int
678 const unsigned int index,
679 const unsigned char combined_orientation) const
680{
681 Assert(combined_orientation ==
683 combined_orientation ==
686
687 AssertIndexRange(index, this->n_dofs_per_line());
688 Assert(adjust_line_dof_index_for_line_orientation_table.size() ==
689 this->n_dofs_per_line(),
691 if (combined_orientation ==
693 return index;
694 else
695 return index + adjust_line_dof_index_for_line_orientation_table[index];
696}
697
698
699
700template <int dim, int spacedim>
701bool
703{
704 for (const unsigned int ref_case :
705 RefinementCase<dim>::all_refinement_cases())
706 if (ref_case != RefinementCase<dim>::no_refinement)
707 for (unsigned int c = 0;
708 c < this->reference_cell().n_children(RefinementCase<dim>(ref_case));
709 ++c)
710 {
711 // make sure also the lazily initialized matrices are created
712 get_prolongation_matrix(c, RefinementCase<dim>(ref_case));
713 Assert((prolongation[ref_case - 1][c].m() ==
714 this->n_dofs_per_cell()) ||
715 (prolongation[ref_case - 1][c].m() == 0),
717 Assert((prolongation[ref_case - 1][c].n() ==
718 this->n_dofs_per_cell()) ||
719 (prolongation[ref_case - 1][c].n() == 0),
721 if ((prolongation[ref_case - 1][c].m() == 0) ||
722 (prolongation[ref_case - 1][c].n() == 0))
723 return false;
724 }
725 return true;
726}
727
728
729
730template <int dim, int spacedim>
731bool
733{
734 for (const unsigned int ref_case :
735 RefinementCase<dim>::all_refinement_cases())
736 if (ref_case != RefinementCase<dim>::no_refinement)
737 for (unsigned int c = 0;
738 c < this->reference_cell().n_children(RefinementCase<dim>(ref_case));
739 ++c)
740 {
741 // make sure also the lazily initialized matrices are created
742 get_restriction_matrix(c, RefinementCase<dim>(ref_case));
743 Assert((restriction[ref_case - 1][c].m() ==
744 this->n_dofs_per_cell()) ||
745 (restriction[ref_case - 1][c].m() == 0),
747 Assert((restriction[ref_case - 1][c].n() ==
748 this->n_dofs_per_cell()) ||
749 (restriction[ref_case - 1][c].n() == 0),
751 if ((restriction[ref_case - 1][c].m() == 0) ||
752 (restriction[ref_case - 1][c].n() == 0))
753 return false;
754 }
755 return true;
756}
757
758
759
760template <int dim, int spacedim>
761bool
763{
764 const RefinementCase<dim> ref_case =
766
767 for (unsigned int c = 0;
768 c < this->reference_cell().n_children(RefinementCase<dim>(ref_case));
769 ++c)
770 {
771 // make sure also the lazily initialized matrices are created
772 get_prolongation_matrix(c, RefinementCase<dim>(ref_case));
773 Assert((prolongation[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
774 (prolongation[ref_case - 1][c].m() == 0),
776 Assert((prolongation[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
777 (prolongation[ref_case - 1][c].n() == 0),
779 if ((prolongation[ref_case - 1][c].m() == 0) ||
780 (prolongation[ref_case - 1][c].n() == 0))
781 return false;
782 }
783 return true;
784}
785
786
787
788template <int dim, int spacedim>
789bool
791{
792 const RefinementCase<dim> ref_case =
794
795 for (unsigned int c = 0;
796 c < this->reference_cell().n_children(RefinementCase<dim>(ref_case));
797 ++c)
798 {
799 // make sure also the lazily initialized matrices are created
800 get_restriction_matrix(c, RefinementCase<dim>(ref_case));
801 Assert((restriction[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
802 (restriction[ref_case - 1][c].m() == 0),
804 Assert((restriction[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
805 (restriction[ref_case - 1][c].n() == 0),
807 if ((restriction[ref_case - 1][c].m() == 0) ||
808 (restriction[ref_case - 1][c].n() == 0))
809 return false;
810 }
811 return true;
812}
813
814
815
816template <int dim, int spacedim>
817bool
819 const internal::SubfaceCase<dim> &subface_case) const
820{
822 {
823 unsigned int n_dofs_on_faces = 0;
824
825 for (const auto face_no : this->reference_cell().face_indices())
826 n_dofs_on_faces += this->n_dofs_per_face(face_no);
827
828 return (n_dofs_on_faces == 0) || (interface_constraints.m() != 0);
829 }
830 else
831 return false;
832}
833
834
835
836template <int dim, int spacedim>
837bool
839{
840 return false;
841}
842
843
844
845template <int dim, int spacedim>
846const FullMatrix<double> &
848 const internal::SubfaceCase<dim> &subface_case) const
849{
850 // TODO: the implementation makes the assumption that all faces have the
851 // same number of dofs
852 AssertDimension(this->n_unique_faces(), 1);
853 const unsigned int face_no = 0;
854 (void)face_no;
855
856 (void)subface_case;
858 ExcMessage("Constraints for this element are only implemented "
859 "for the case that faces are refined isotropically "
860 "(which is always the case in 2d, and in 3d requires "
861 "that the neighboring cell of a coarse cell presents "
862 "exactly four children on the common face)."));
863 Assert((this->n_dofs_per_face(face_no) == 0) ||
864 (interface_constraints.m() != 0),
865 ExcMessage("The finite element for which you try to obtain "
866 "hanging node constraints does not appear to "
867 "implement them."));
868
869 if (dim == 1)
870 Assert((interface_constraints.m() == 0) && (interface_constraints.n() == 0),
871 ExcWrongInterfaceMatrixSize(interface_constraints.m(),
872 interface_constraints.n()));
873
874 return interface_constraints;
875}
876
877
878
879template <int dim, int spacedim>
882{
883 // TODO: the implementation makes the assumption that all faces have the
884 // same number of dofs
885 AssertDimension(this->n_unique_faces(), 1);
886 const unsigned int face_no = 0;
887
888 switch (dim)
889 {
890 case 1:
891 return {0U, 0U};
892 case 2:
893 return {this->n_dofs_per_vertex() + 2 * this->n_dofs_per_line(),
894 this->n_dofs_per_face(face_no)};
895 case 3:
896 return {5 * this->n_dofs_per_vertex() + 12 * this->n_dofs_per_line() +
897 4 * this->n_dofs_per_quad(face_no),
898 this->n_dofs_per_face(face_no)};
899 default:
901 }
903}
904
905
906
907template <int dim, int spacedim>
908void
911 FullMatrix<double> &) const
912{
913 // by default, no interpolation
914 // implemented. so throw exception,
915 // as documentation says
917 false,
919}
920
921
922
923template <int dim, int spacedim>
924void
928 const unsigned int) const
929{
930 // by default, no interpolation
931 // implemented. so throw exception,
932 // as documentation says
934 false,
936}
937
938
939
940template <int dim, int spacedim>
941void
944 const unsigned int,
946 const unsigned int) const
947{
948 // by default, no interpolation
949 // implemented. so throw exception,
950 // as documentation says
952 false,
954}
955
956
957
958template <int dim, int spacedim>
959std::vector<std::pair<unsigned int, unsigned int>>
961 const FiniteElement<dim, spacedim> &) const
962{
964 return std::vector<std::pair<unsigned int, unsigned int>>();
965}
966
967
968
969template <int dim, int spacedim>
970std::vector<std::pair<unsigned int, unsigned int>>
972 const FiniteElement<dim, spacedim> &) const
973{
975 return std::vector<std::pair<unsigned int, unsigned int>>();
976}
977
978
979
980template <int dim, int spacedim>
981std::vector<std::pair<unsigned int, unsigned int>>
984 const unsigned int) const
985{
987 return std::vector<std::pair<unsigned int, unsigned int>>();
988}
989
990
991
992template <int dim, int spacedim>
996 const unsigned int) const
997{
1000}
1001
1002
1003
1004template <int dim, int spacedim>
1005bool
1007 const FiniteElement<dim, spacedim> &f) const
1008{
1009 // Compare fields in roughly increasing order of how expensive the
1010 // comparison is
1011 return ((typeid(*this) == typeid(f)) && (this->get_name() == f.get_name()) &&
1012 (static_cast<const FiniteElementData<dim> &>(*this) ==
1013 static_cast<const FiniteElementData<dim> &>(f)) &&
1014 (interface_constraints == f.interface_constraints));
1015}
1016
1017
1018
1019template <int dim, int spacedim>
1020bool
1022 const FiniteElement<dim, spacedim> &f) const
1023{
1024 return !(*this == f);
1025}
1026
1027
1028
1029template <int dim, int spacedim>
1030const std::vector<Point<dim>> &
1032{
1033 // a finite element may define
1034 // support points, but only if
1035 // there are as many as there are
1036 // degrees of freedom
1037 Assert((unit_support_points.empty()) ||
1038 (unit_support_points.size() == this->n_dofs_per_cell()),
1040 return unit_support_points;
1041}
1042
1043
1044
1045template <int dim, int spacedim>
1046bool
1048{
1049 return (unit_support_points.size() != 0);
1050}
1051
1052
1053
1054template <int dim, int spacedim>
1055const std::vector<Point<dim>> &
1057{
1058 // If the finite element implements generalized support points, return
1059 // those. Otherwise fall back to unit support points.
1060 return ((generalized_support_points.empty()) ? unit_support_points :
1061 generalized_support_points);
1062}
1063
1064
1065
1066template <int dim, int spacedim>
1067bool
1069{
1070 return (get_generalized_support_points().size() != 0);
1071}
1072
1073
1074
1075template <int dim, int spacedim>
1077FiniteElement<dim, spacedim>::unit_support_point(const unsigned int index) const
1078{
1079 AssertIndexRange(index, this->n_dofs_per_cell());
1080 Assert(unit_support_points.size() == this->n_dofs_per_cell(),
1081 ExcFEHasNoSupportPoints());
1082 return unit_support_points[index];
1083}
1084
1085
1086
1087template <int dim, int spacedim>
1088const std::vector<Point<dim - 1>> &
1090 const unsigned int face_no) const
1091{
1092 // a finite element may define
1093 // support points, but only if
1094 // there are as many as there are
1095 // degrees of freedom on a face
1096 Assert((unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1097 .empty()) ||
1098 (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1099 .size() == this->n_dofs_per_face(face_no)),
1101 return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no];
1102}
1103
1104
1105
1106template <int dim, int spacedim>
1107bool
1109 const unsigned int face_no) const
1110{
1111 return (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1112 .size() != 0);
1113}
1114
1115
1116
1117template <int dim, int spacedim>
1118Point<dim - 1>
1120 const unsigned int index,
1121 const unsigned int face_no) const
1122{
1123 AssertIndexRange(index, this->n_dofs_per_face(face_no));
1124 Assert(unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1125 .size() == this->n_dofs_per_face(face_no),
1126 ExcFEHasNoSupportPoints());
1127 return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1128 [index];
1129}
1130
1131
1132
1133template <int dim, int spacedim>
1134bool
1136 const unsigned int) const
1137{
1138 return true;
1139}
1140
1141
1142
1143template <int dim, int spacedim>
1146{
1147 // Translate the ComponentMask into first_selected and n_components after
1148 // some error checking:
1149 const unsigned int n_total_components = this->n_components();
1150 Assert((n_total_components == mask.size()) || (mask.size() == 0),
1151 ExcMessage("The given ComponentMask has the wrong size."));
1152
1153 const unsigned int n_selected =
1154 mask.n_selected_components(n_total_components);
1155 Assert(n_selected > 0,
1156 ExcMessage("You need at least one selected component."));
1157
1158 const unsigned int first_selected =
1159 mask.first_selected_component(n_total_components);
1160
1161# ifdef DEBUG
1162 // check that it is contiguous:
1163 for (unsigned int c = 0; c < n_total_components; ++c)
1164 Assert((c < first_selected && (!mask[c])) ||
1165 (c >= first_selected && c < first_selected + n_selected &&
1166 mask[c]) ||
1167 (c >= first_selected + n_selected && !mask[c]),
1168 ExcMessage("Error: the given ComponentMask is not contiguous!"));
1169# endif
1170
1171 return get_sub_fe(first_selected, n_selected);
1172}
1173
1174
1175
1176template <int dim, int spacedim>
1179 const unsigned int first_component,
1180 const unsigned int n_selected_components) const
1181{
1182 // No complicated logic is needed here, because it is overridden in
1183 // FESystem<dim,spacedim>. Just make sure that what the user chose is valid:
1184 Assert(first_component == 0 && n_selected_components == this->n_components(),
1185 ExcMessage(
1186 "You can only select a whole FiniteElement, not a part of one."));
1187
1188 (void)first_component;
1189 (void)n_selected_components;
1190 return *this;
1191}
1192
1193
1194
1195template <int dim, int spacedim>
1196std::pair<Table<2, bool>, std::vector<unsigned int>>
1198{
1200 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1201 Table<2, bool>(this->n_components(), this->n_dofs_per_cell()),
1202 std::vector<unsigned int>(this->n_components()));
1203}
1204
1205
1206
1207template <int dim, int spacedim>
1208void
1211 const std::vector<Vector<double>> &,
1212 std::vector<double> &) const
1213{
1214 Assert(has_generalized_support_points(),
1215 ExcMessage("The element for which you are calling the current "
1216 "function does not have generalized support points (see "
1217 "the glossary for a definition of generalized support "
1218 "points). Consequently, the current function can not "
1219 "be defined and is not implemented by the element."));
1221}
1222
1223
1224
1225template <int dim, int spacedim>
1226std::size_t
1228{
1229 return (
1230 sizeof(FiniteElementData<dim>) +
1233 MemoryConsumption::memory_consumption(interface_constraints) +
1234 MemoryConsumption::memory_consumption(system_to_component_table) +
1235 MemoryConsumption::memory_consumption(face_system_to_component_table) +
1236 MemoryConsumption::memory_consumption(system_to_base_table) +
1237 MemoryConsumption::memory_consumption(face_system_to_base_table) +
1238 MemoryConsumption::memory_consumption(component_to_base_table) +
1239 MemoryConsumption::memory_consumption(restriction_is_additive_flags) +
1240 MemoryConsumption::memory_consumption(nonzero_components) +
1241 MemoryConsumption::memory_consumption(n_nonzero_components_table));
1242}
1243
1244
1245
1246template <int dim, int spacedim>
1247std::vector<unsigned int>
1249 const std::vector<ComponentMask> &nonzero_components)
1250{
1251 std::vector<unsigned int> retval(nonzero_components.size());
1252 for (unsigned int i = 0; i < nonzero_components.size(); ++i)
1253 retval[i] = nonzero_components[i].n_selected_components();
1254 return retval;
1255}
1256
1257
1258
1259/*------------------------------- FiniteElement ----------------------*/
1260
1261# ifndef DOXYGEN
1262template <int dim, int spacedim>
1263std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1265 const UpdateFlags flags,
1266 const Mapping<dim, spacedim> &mapping,
1267 const hp::QCollection<dim - 1> &quadrature,
1269 spacedim>
1270 &output_data) const
1271{
1272 return get_data(flags,
1273 mapping,
1275 quadrature),
1276 output_data);
1277}
1278
1279
1280
1281template <int dim, int spacedim>
1282std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1284 const UpdateFlags flags,
1285 const Mapping<dim, spacedim> &mapping,
1286 const Quadrature<dim - 1> &quadrature,
1288 spacedim>
1289 &output_data) const
1290{
1291 return get_data(flags,
1292 mapping,
1294 quadrature),
1295 output_data);
1296}
1297
1298
1299
1300template <int dim, int spacedim>
1301inline void
1304 const unsigned int face_no,
1305 const hp::QCollection<dim - 1> &quadrature,
1306 const Mapping<dim, spacedim> &mapping,
1307 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
1309 &mapping_data,
1310 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
1312 spacedim>
1313 &output_data) const
1314{
1315 // base class version, implement overridden function in derived classes
1316 AssertDimension(quadrature.size(), 1);
1317 fill_fe_face_values(cell,
1318 face_no,
1319 quadrature[0],
1320 mapping,
1321 mapping_internal,
1322 mapping_data,
1323 fe_internal,
1324 output_data);
1325}
1326
1327
1328
1329template <int dim, int spacedim>
1330inline void
1333 const unsigned int face_no,
1334 const Quadrature<dim - 1> &quadrature,
1335 const Mapping<dim, spacedim> &mapping,
1336 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
1338 &mapping_data,
1339 const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
1341 spacedim>
1342 &output_data) const
1343{
1344 Assert(false,
1345 ExcMessage("Use of a deprecated interface, please implement "
1346 "fill_fe_face_values taking a hp::QCollection argument"));
1347 (void)cell;
1348 (void)face_no;
1349 (void)quadrature;
1350 (void)mapping;
1351 (void)mapping_internal;
1352 (void)mapping_data;
1353 (void)fe_internal;
1354 (void)output_data;
1355}
1356# endif
1357
1358
1359
1360template <int dim, int spacedim>
1361std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1363 const UpdateFlags flags,
1364 const Mapping<dim, spacedim> &mapping,
1365 const Quadrature<dim - 1> &quadrature,
1367 spacedim>
1368 &output_data) const
1369{
1370 return get_data(flags,
1371 mapping,
1373 this->reference_cell(), quadrature),
1374 output_data);
1375}
1376
1377
1378
1379template <int dim, int spacedim>
1381FiniteElement<dim, spacedim>::base_element(const unsigned int index) const
1382{
1383 (void)index;
1384 AssertIndexRange(index, 1);
1385 // This function should not be
1386 // called for a system element
1387 Assert(base_to_block_indices.size() == 1, ExcInternalError());
1388 return *this;
1389}
1390
1391
1392#endif
1393/*------------------------------- Explicit Instantiations -------------*/
1394#include "fe.inst"
1395
1396
bool represents_the_all_selected_mask() const
unsigned int size() const
bool represents_the_all_selected_mask() const
unsigned int size() const
const unsigned int components
Definition fe_data.h:446
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int n_components() const
virtual std::size_t memory_consumption() const
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
bool isotropic_prolongation_is_implemented() const
virtual std::string get_name() const =0
virtual Point< dim > unit_support_point(const unsigned int index) const
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
const std::vector< unsigned int > n_nonzero_components_table
Definition fe.h:2597
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool prolongation_is_implemented() const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index, const unsigned int face_no=0) const
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
virtual std::unique_ptr< InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
bool has_face_support_points(const unsigned int face_no=0) const
const bool cached_primitivity
Definition fe.h:2604
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
bool operator!=(const FiniteElement< dim, spacedim > &) const
bool has_support_points() const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const unsigned char combined_orientation=ReferenceCell::default_combined_face_orientation()) const
bool has_generalized_support_points() const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
const std::vector< Point< dim > > & get_unit_support_points() const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
bool isotropic_restriction_is_implemented() const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition fe.h:2500
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const
unsigned int component_to_block_index(const unsigned int component) const
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition fe.h:2536
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const unsigned int face_no, const unsigned char combined_orientation) const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
const std::vector< bool > restriction_is_additive_flags
Definition fe.h:2579
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const unsigned char combined_orientation) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition fe.h:2572
TableIndices< 2 > interface_constraints_size() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
bool restriction_is_implemented() const
FullMatrix< double > interface_constraints
Definition fe.h:2439
const std::vector< Point< dim > > & get_generalized_support_points() const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
virtual std::size_t memory_consumption() const
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
const std::vector< ComponentMask > nonzero_components
Definition fe.h:2588
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual bool hp_constraints_are_implemented() const
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
static constexpr unsigned char default_combined_face_orientation()
static constexpr unsigned char reversed_combined_line_orientation()
unsigned int size() const
Definition collection.h:308
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 2 > first
Definition grid_out.cc:4623
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
UpdateFlags
@ update_default
No update.
#define DEAL_II_NOT_IMPLEMENTED()
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Definition operators.h:49
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::vector< Point< dim > > unit_support_points(const std::vector< Point< 1 > > &line_support_points, const std::vector< unsigned int > &renumbering)
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.