Reference documentation for deal.II version Git ede8f93e86 2020-12-03 14:59:20 -0700
\(\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\}}\)
fe.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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 
19 
21 
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/fe/mapping.h>
25 
26 #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 
39 
40 template <int dim, int spacedim>
42  : update_each(update_default)
43 {}
44 
45 
46 
47 template <int dim, int spacedim>
48 std::size_t
50 {
51  return sizeof(*this);
52 }
53 
54 
55 
56 template <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)
62  , adjust_line_dof_index_for_line_orientation_table(
63  dim == 3 ? this->n_dofs_per_line() : 0)
64  , system_to_base_table(this->n_dofs_per_cell())
65  , component_to_base_table(this->components,
66  std::make_pair(std::make_pair(0U, 0U), 0U))
67  ,
68 
69  // Special handling of vectors of length one: in this case, we
70  // assume that all entries were supposed to be equal
71  restriction_is_additive_flags(
72  r_i_a_f.size() == 1 ?
73  std::vector<bool>(fe_data.n_dofs_per_cell(), r_i_a_f[0]) :
74  r_i_a_f)
75  , nonzero_components(
76  nonzero_c.size() == 1 ?
77  std::vector<ComponentMask>(fe_data.n_dofs_per_cell(), nonzero_c[0]) :
78  nonzero_c)
79  , n_nonzero_components_table(compute_n_nonzero_components(nonzero_components))
80  , cached_primitivity(std::find_if(n_nonzero_components_table.begin(),
81  n_nonzero_components_table.end(),
82  [](const unsigned int n_components) {
83  return n_components != 1U;
84  }) == n_nonzero_components_table.end())
85 {
88  this->n_dofs_per_cell()));
90  for (unsigned int i = 0; i < nonzero_components.size(); ++i)
91  {
92  Assert(nonzero_components[i].size() == this->n_components(),
94  Assert(nonzero_components[i].n_selected_components() >= 1,
99  }
100 
101  // initialize some tables in the default way, i.e. if there is only one
102  // (vector-)component; if the element is not primitive, leave these tables
103  // empty.
104  if (this->is_primitive())
105  {
107  for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
108  system_to_component_table[j] = std::pair<unsigned, unsigned>(0, j);
109 
111  for (unsigned int f = 0; f < this->n_unique_faces(); ++f)
112  {
114  for (unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
116  std::pair<unsigned, unsigned>(0, j);
117  }
118  }
119 
120  for (unsigned int j = 0; j < this->n_dofs_per_cell(); ++j)
121  system_to_base_table[j] = std::make_pair(std::make_pair(0U, 0U), j);
122 
124  for (unsigned int f = 0; f < this->n_unique_faces(); ++f)
125  {
126  face_system_to_base_table[f].resize(this->n_dofs_per_face(f));
127  for (unsigned int j = 0; j < this->n_dofs_per_face(f); ++j)
129  std::make_pair(std::make_pair(0U, 0U), j);
130  }
131 
132  // Fill with default value; may be changed by constructor of derived class.
134 
135  // initialize the restriction and prolongation matrices. the default
136  // constructor of FullMatrix<dim> initializes them with size zero
139  for (unsigned int ref = RefinementCase<dim>::cut_x;
140  ref < RefinementCase<dim>::isotropic_refinement + 1;
141  ++ref)
142  {
144  RefinementCase<dim>(ref)),
147  RefinementCase<dim>(ref)),
149  }
150 
151 
152  if (dim == 3)
153  {
155  this->n_unique_quads());
156 
157  for (unsigned int f = 0; f < this->n_unique_quads(); ++f)
158  {
160  this->n_dofs_per_quad(f),
163  8 :
164  6);
166  }
167  }
168 
171 }
172 
173 
174 
175 template <int dim, int spacedim>
176 std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>
177 FiniteElement<dim, spacedim>::operator^(const unsigned int multiplicity) const
178 {
179  return {this->clone(), multiplicity};
180 }
181 
182 
183 
184 template <int dim, int spacedim>
185 double
187  const Point<dim> &) const
188 {
190  return 0.;
191 }
192 
193 
194 
195 template <int dim, int spacedim>
196 double
198  const Point<dim> &,
199  const unsigned int) const
200 {
202  return 0.;
203 }
204 
205 
206 
207 template <int dim, int spacedim>
210  const Point<dim> &) const
211 {
213  return Tensor<1, dim>();
214 }
215 
216 
217 
218 template <int dim, int spacedim>
221  const Point<dim> &,
222  const unsigned int) const
223 {
225  return Tensor<1, dim>();
226 }
227 
228 
229 
230 template <int dim, int spacedim>
233  const Point<dim> &) const
234 {
236  return Tensor<2, dim>();
237 }
238 
239 
240 
241 template <int dim, int spacedim>
244  const unsigned int,
245  const Point<dim> &,
246  const unsigned int) const
247 {
249  return Tensor<2, dim>();
250 }
251 
252 
253 
254 template <int dim, int spacedim>
257  const Point<dim> &) const
258 {
260  return Tensor<3, dim>();
261 }
262 
263 
264 
265 template <int dim, int spacedim>
268  const unsigned int,
269  const Point<dim> &,
270  const unsigned int) const
271 {
273  return Tensor<3, dim>();
274 }
275 
276 
277 
278 template <int dim, int spacedim>
281  const Point<dim> &) const
282 {
284  return Tensor<4, dim>();
285 }
286 
287 
288 
289 template <int dim, int spacedim>
292  const unsigned int,
293  const Point<dim> &,
294  const unsigned int) const
295 {
297  return Tensor<4, dim>();
298 }
299 
300 
301 template <int dim, int spacedim>
302 void
304  const bool isotropic_restriction_only,
305  const bool isotropic_prolongation_only)
306 {
307  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
308  ref_case <= RefinementCase<dim>::isotropic_refinement;
309  ++ref_case)
310  {
311  const unsigned int nc =
313 
314  for (unsigned int i = 0; i < nc; ++i)
315  {
316  if (this->restriction[ref_case - 1][i].m() !=
317  this->n_dofs_per_cell() &&
318  (!isotropic_restriction_only ||
320  this->restriction[ref_case - 1][i].reinit(this->n_dofs_per_cell(),
321  this->n_dofs_per_cell());
322  if (this->prolongation[ref_case - 1][i].m() !=
323  this->n_dofs_per_cell() &&
324  (!isotropic_prolongation_only ||
326  this->prolongation[ref_case - 1][i].reinit(this->n_dofs_per_cell(),
327  this->n_dofs_per_cell());
328  }
329  }
330 }
331 
332 
333 template <int dim, int spacedim>
334 const FullMatrix<double> &
336  const unsigned int child,
337  const RefinementCase<dim> &refinement_case) const
338 {
339  AssertIndexRange(refinement_case,
341  Assert(refinement_case != RefinementCase<dim>::no_refinement,
342  ExcMessage(
343  "Restriction matrices are only available for refined cells!"));
345  child, GeometryInfo<dim>::n_children(RefinementCase<dim>(refinement_case)));
346  // we use refinement_case-1 here. the -1 takes care of the origin of the
347  // vector, as for RefinementCase<dim>::no_refinement (=0) there is no data
348  // available and so the vector indices are shifted
349  Assert(restriction[refinement_case - 1][child].n() == this->n_dofs_per_cell(),
351  return restriction[refinement_case - 1][child];
352 }
353 
354 
355 
356 template <int dim, int spacedim>
357 const FullMatrix<double> &
359  const unsigned int child,
360  const RefinementCase<dim> &refinement_case) const
361 {
362  AssertIndexRange(refinement_case,
364  Assert(refinement_case != RefinementCase<dim>::no_refinement,
365  ExcMessage(
366  "Prolongation matrices are only available for refined cells!"));
368  child, GeometryInfo<dim>::n_children(RefinementCase<dim>(refinement_case)));
369  // we use refinement_case-1 here. the -1 takes care
370  // of the origin of the vector, as for
371  // RefinementCase::no_refinement (=0) there is no
372  // data available and so the vector indices
373  // are shifted
374  Assert(prolongation[refinement_case - 1][child].n() ==
375  this->n_dofs_per_cell(),
376  ExcEmbeddingVoid());
377  return prolongation[refinement_case - 1][child];
378 }
379 
380 
381 // TODO:[GK] This is probably not the most efficient way of doing this.
382 template <int dim, int spacedim>
383 unsigned int
385  const unsigned int index) const
386 {
387  AssertIndexRange(index, this->n_components());
388 
389  return first_block_of_base(component_to_base_table[index].first.first) +
390  component_to_base_table[index].second;
391 }
392 
393 
394 template <int dim, int spacedim>
397  const FEValuesExtractors::Scalar &scalar) const
398 {
399  AssertIndexRange(scalar.component, this->n_components());
400 
401  // TODO: it would be nice to verify that it is indeed possible
402  // to select this scalar component, i.e., that it is not part
403  // of a non-primitive element. unfortunately, there is no simple
404  // way to write such a condition...
405 
406  std::vector<bool> mask(this->n_components(), false);
407  mask[scalar.component] = true;
408  return mask;
409 }
410 
411 
412 template <int dim, int spacedim>
415  const FEValuesExtractors::Vector &vector) const
416 {
417  AssertIndexRange(vector.first_vector_component + dim - 1,
418  this->n_components());
419 
420  // TODO: it would be nice to verify that it is indeed possible
421  // to select these vector components, i.e., that they don't span
422  // beyond the beginning or end of anon-primitive element.
423  // unfortunately, there is no simple way to write such a condition...
424 
425  std::vector<bool> mask(this->n_components(), false);
426  for (unsigned int c = vector.first_vector_component;
427  c < vector.first_vector_component + dim;
428  ++c)
429  mask[c] = true;
430  return mask;
431 }
432 
433 
434 template <int dim, int spacedim>
437  const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const
438 {
441  this->n_components());
442 
443  // TODO: it would be nice to verify that it is indeed possible
444  // to select these vector components, i.e., that they don't span
445  // beyond the beginning or end of anon-primitive element.
446  // unfortunately, there is no simple way to write such a condition...
447 
448  std::vector<bool> mask(this->n_components(), false);
449  for (unsigned int c = sym_tensor.first_tensor_component;
450  c < sym_tensor.first_tensor_component +
452  ++c)
453  mask[c] = true;
454  return mask;
455 }
456 
457 
458 
459 template <int dim, int spacedim>
462 {
463  // if we get a block mask that represents all blocks, then
464  // do the same for the returned component mask
465  if (block_mask.represents_the_all_selected_mask())
466  return {};
467 
468  AssertDimension(block_mask.size(), this->n_blocks());
469 
470  std::vector<bool> component_mask(this->n_components(), false);
471  for (unsigned int c = 0; c < this->n_components(); ++c)
472  if (block_mask[component_to_block_index(c)] == true)
473  component_mask[c] = true;
474 
475  return component_mask;
476 }
477 
478 
479 
480 template <int dim, int spacedim>
481 BlockMask
483  const FEValuesExtractors::Scalar &scalar) const
484 {
485  // simply create the corresponding component mask (a simpler
486  // process) and then convert it to a block mask
487  return block_mask(component_mask(scalar));
488 }
489 
490 
491 template <int dim, int spacedim>
492 BlockMask
494  const FEValuesExtractors::Vector &vector) const
495 {
496  // simply create the corresponding component mask (a simpler
497  // process) and then convert it to a block mask
498  return block_mask(component_mask(vector));
499 }
500 
501 
502 template <int dim, int spacedim>
503 BlockMask
505  const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const
506 {
507  // simply create the corresponding component mask (a simpler
508  // process) and then convert it to a block mask
509  return block_mask(component_mask(sym_tensor));
510 }
511 
512 
513 
514 template <int dim, int spacedim>
515 BlockMask
517  const ComponentMask &component_mask) const
518 {
519  // if we get a component mask that represents all component, then
520  // do the same for the returned block mask
521  if (component_mask.represents_the_all_selected_mask())
522  return {};
523 
524  AssertDimension(component_mask.size(), this->n_components());
525 
526  // walk over all of the components
527  // of this finite element and see
528  // if we need to set the
529  // corresponding block. inside the
530  // block, walk over all the
531  // components that correspond to
532  // this block and make sure the
533  // component mask is set for all of
534  // them
535  std::vector<bool> block_mask(this->n_blocks(), false);
536  for (unsigned int c = 0; c < this->n_components();)
537  {
538  const unsigned int block = component_to_block_index(c);
539  if (component_mask[c] == true)
540  block_mask[block] = true;
541 
542  // now check all of the other
543  // components that correspond
544  // to this block
545  ++c;
546  while ((c < this->n_components()) &&
547  (component_to_block_index(c) == block))
548  {
549  Assert(component_mask[c] == block_mask[block],
550  ExcMessage(
551  "The component mask argument given to this function "
552  "is not a mask where the individual components belonging "
553  "to one block of the finite element are either all "
554  "selected or not selected. You can't call this function "
555  "with a component mask that splits blocks."));
556  ++c;
557  }
558  }
559 
560 
561  return block_mask;
562 }
563 
564 
565 
566 template <int dim, int spacedim>
567 unsigned int
569  const unsigned int face,
570  const bool face_orientation,
571  const bool face_flip,
572  const bool face_rotation) const
573 {
574  const auto &refence_cell =
576 
577  AssertIndexRange(face_index, this->n_dofs_per_face(face));
578  AssertIndexRange(face, refence_cell.n_faces());
579 
580  // TODO: we could presumably solve the 3d case below using the
581  // adjust_quad_dof_index_for_face_orientation_table field. for the
582  // 2d case, we can't use adjust_line_dof_index_for_line_orientation_table
583  // since that array is empty (presumably because we thought that
584  // there are no flipped edges in 2d, but these can happen in
585  // DoFTools::make_periodicity_constraints, for example). so we
586  // would need to either fill this field, or rely on derived classes
587  // implementing this function, as we currently do
588 
589  // see the function's documentation for an explanation of this
590  // assertion -- in essence, derived classes have to implement
591  // an overloaded version of this function if we are to use any
592  // other than standard orientation
593  if ((face_orientation != true) || (face_flip != false) ||
594  (face_rotation != false))
595  Assert((this->n_dofs_per_line() <= 1) && (this->n_dofs_per_quad(face) <= 1),
596  ExcMessage(
597  "The function in this base class can not handle this case. "
598  "Rather, the derived class you are using must provide "
599  "an overloaded version but apparently hasn't done so. See "
600  "the documentation of this function for more information."));
601 
602  // we need to distinguish between DoFs on vertices, lines and in 3d quads.
603  // do so in a sequence of if-else statements
604  if (face_index < this->get_first_face_line_index(face))
605  // DoF is on a vertex
606  {
607  // get the number of the vertex on the face that corresponds to this DoF,
608  // along with the number of the DoF on this vertex
609  const unsigned int face_vertex = face_index / this->n_dofs_per_vertex();
610  const unsigned int dof_index_on_vertex =
611  face_index % this->n_dofs_per_vertex();
612 
613  // then get the number of this vertex on the cell and translate
614  // this to a DoF number on the cell
615  return (refence_cell.face_to_cell_vertices(face,
616  face_vertex,
617  face_orientation +
618  2 * face_rotation +
619  4 * face_flip) *
620  this->n_dofs_per_vertex() +
621  dof_index_on_vertex);
622  }
623  else if (face_index < this->get_first_face_quad_index(face))
624  // DoF is on a face
625  {
626  // do the same kind of translation as before. we need to only consider
627  // DoFs on the lines, i.e., ignoring those on the vertices
628  const unsigned int index =
629  face_index - this->get_first_face_line_index(face);
630 
631  const unsigned int face_line = index / this->n_dofs_per_line();
632  const unsigned int dof_index_on_line = index % this->n_dofs_per_line();
633 
634  return (this->get_first_line_index() +
635  refence_cell.face_to_cell_lines(face,
636  face_line,
637  face_orientation +
638  2 * face_rotation +
639  4 * face_flip) *
640  this->n_dofs_per_line() +
641  dof_index_on_line);
642  }
643  else
644  // DoF is on a quad
645  {
646  Assert(dim >= 3, ExcInternalError());
647 
648  // ignore vertex and line dofs
649  const unsigned int index =
650  face_index - this->get_first_face_quad_index(face);
651 
652  return (this->get_first_quad_index(face) + index);
653  }
654 }
655 
656 
657 
658 template <int dim, int spacedim>
659 unsigned int
661  const unsigned int index,
662  const unsigned int face,
663  const bool face_orientation,
664  const bool face_flip,
665  const bool face_rotation) const
666 {
667  // general template for 1D and 2D: not
668  // implemented. in fact, the function
669  // shouldn't even be called unless we are
670  // in 3d, so throw an internal error
671  Assert(dim == 3, ExcInternalError());
672  if (dim < 3)
673  return index;
674 
675  // adjust dofs on 3d faces if the face is
676  // flipped. note that we query a table that
677  // derived elements need to have set up
678  // front. the exception are discontinuous
679  // elements for which there should be no
680  // face dofs anyway (i.e. dofs_per_quad==0
681  // in 3d), so we don't need the table, but
682  // the function should also not have been
683  // called
684  AssertIndexRange(index, this->n_dofs_per_quad(face));
686  [this->n_unique_quads() == 1 ? 0 : face]
687  .n_elements() ==
689  .face_reference_cell_type(face) == ReferenceCell::Type::Quad ?
690  8 :
691  6) *
692  this->n_dofs_per_quad(face),
693  ExcInternalError());
694  return index +
696  [this->n_unique_quads() == 1 ? 0 : face](
697  index, 4 * face_orientation + 2 * face_flip + face_rotation);
698 }
699 
700 
701 
702 template <int dim, int spacedim>
703 unsigned int
705  const unsigned int index,
706  const bool line_orientation) const
707 {
708  // general template for 1D and 2D: do
709  // nothing. Do not throw an Assertion,
710  // however, in order to allow to call this
711  // function in 2D as well
712  if (dim < 3)
713  return index;
714 
715  AssertIndexRange(index, this->n_dofs_per_line());
717  this->n_dofs_per_line(),
718  ExcInternalError());
719  if (line_orientation)
720  return index;
721  else
723 }
724 
725 
726 
727 template <int dim, int spacedim>
728 bool
730 {
731  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
732  ref_case < RefinementCase<dim>::isotropic_refinement + 1;
733  ++ref_case)
734  for (unsigned int c = 0;
735  c < GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
736  ++c)
737  {
738  // make sure also the lazily initialized matrices are created
740  Assert((prolongation[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
741  (prolongation[ref_case - 1][c].m() == 0),
742  ExcInternalError());
743  Assert((prolongation[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
744  (prolongation[ref_case - 1][c].n() == 0),
745  ExcInternalError());
746  if ((prolongation[ref_case - 1][c].m() == 0) ||
747  (prolongation[ref_case - 1][c].n() == 0))
748  return false;
749  }
750  return true;
751 }
752 
753 
754 
755 template <int dim, int spacedim>
756 bool
758 {
759  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
760  ref_case < RefinementCase<dim>::isotropic_refinement + 1;
761  ++ref_case)
762  for (unsigned int c = 0;
763  c < GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
764  ++c)
765  {
766  // make sure also the lazily initialized matrices are created
768  Assert((restriction[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
769  (restriction[ref_case - 1][c].m() == 0),
770  ExcInternalError());
771  Assert((restriction[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
772  (restriction[ref_case - 1][c].n() == 0),
773  ExcInternalError());
774  if ((restriction[ref_case - 1][c].m() == 0) ||
775  (restriction[ref_case - 1][c].n() == 0))
776  return false;
777  }
778  return true;
779 }
780 
781 
782 
783 template <int dim, int spacedim>
784 bool
786 {
787  const RefinementCase<dim> ref_case =
789 
790  for (unsigned int c = 0;
791  c < GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
792  ++c)
793  {
794  // make sure also the lazily initialized matrices are created
796  Assert((prolongation[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
797  (prolongation[ref_case - 1][c].m() == 0),
798  ExcInternalError());
799  Assert((prolongation[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
800  (prolongation[ref_case - 1][c].n() == 0),
801  ExcInternalError());
802  if ((prolongation[ref_case - 1][c].m() == 0) ||
803  (prolongation[ref_case - 1][c].n() == 0))
804  return false;
805  }
806  return true;
807 }
808 
809 
810 
811 template <int dim, int spacedim>
812 bool
814 {
815  const RefinementCase<dim> ref_case =
817 
818  for (unsigned int c = 0;
819  c < GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
820  ++c)
821  {
822  // make sure also the lazily initialized matrices are created
824  Assert((restriction[ref_case - 1][c].m() == this->n_dofs_per_cell()) ||
825  (restriction[ref_case - 1][c].m() == 0),
826  ExcInternalError());
827  Assert((restriction[ref_case - 1][c].n() == this->n_dofs_per_cell()) ||
828  (restriction[ref_case - 1][c].n() == 0),
829  ExcInternalError());
830  if ((restriction[ref_case - 1][c].m() == 0) ||
831  (restriction[ref_case - 1][c].n() == 0))
832  return false;
833  }
834  return true;
835 }
836 
837 
838 
839 template <int dim, int spacedim>
840 bool
842  const internal::SubfaceCase<dim> &subface_case) const
843 {
844  // TODO: the implementation makes the assumption that all faces have the
845  // same number of dofs
846  AssertDimension(this->n_unique_faces(), 1);
847  const unsigned int face_no = 0;
848 
849  if (subface_case == internal::SubfaceCase<dim>::case_isotropic)
850  return (this->n_dofs_per_face(face_no) == 0) ||
851  (interface_constraints.m() != 0);
852  else
853  return false;
854 }
855 
856 
857 
858 template <int dim, int spacedim>
859 bool
861 {
862  return false;
863 }
864 
865 
866 
867 template <int dim, int spacedim>
868 const FullMatrix<double> &
870  const internal::SubfaceCase<dim> &subface_case) const
871 {
872  // TODO: the implementation makes the assumption that all faces have the
873  // same number of dofs
874  AssertDimension(this->n_unique_faces(), 1);
875  const unsigned int face_no = 0;
876  (void)face_no;
877 
878  (void)subface_case;
880  ExcMessage("Constraints for this element are only implemented "
881  "for the case that faces are refined isotropically "
882  "(which is always the case in 2d, and in 3d requires "
883  "that the neighboring cell of a coarse cell presents "
884  "exactly four children on the common face)."));
885  Assert((this->n_dofs_per_face(face_no) == 0) ||
886  (interface_constraints.m() != 0),
887  ExcMessage("The finite element for which you try to obtain "
888  "hanging node constraints does not appear to "
889  "implement them."));
890 
891  if (dim == 1)
895 
896  return interface_constraints;
897 }
898 
899 
900 
901 template <int dim, int spacedim>
904 {
905  // TODO: the implementation makes the assumption that all faces have the
906  // same number of dofs
907  AssertDimension(this->n_unique_faces(), 1);
908  const unsigned int face_no = 0;
909 
910  switch (dim)
911  {
912  case 1:
913  return {0U, 0U};
914  case 2:
915  return {this->n_dofs_per_vertex() + 2 * this->n_dofs_per_line(),
916  this->n_dofs_per_face(face_no)};
917  case 3:
918  return {5 * this->n_dofs_per_vertex() + 12 * this->n_dofs_per_line() +
919  4 * this->n_dofs_per_quad(face_no),
920  this->n_dofs_per_face(face_no)};
921  default:
922  Assert(false, ExcNotImplemented());
923  }
925 }
926 
927 
928 
929 template <int dim, int spacedim>
930 void
933  FullMatrix<double> &) const
934 {
935  // by default, no interpolation
936  // implemented. so throw exception,
937  // as documentation says
938  AssertThrow(
939  false,
941 }
942 
943 
944 
945 template <int dim, int spacedim>
946 void
950  const unsigned int) const
951 {
952  // by default, no interpolation
953  // implemented. so throw exception,
954  // as documentation says
955  AssertThrow(
956  false,
958 }
959 
960 
961 
962 template <int dim, int spacedim>
963 void
966  const unsigned int,
968  const unsigned int) const
969 {
970  // by default, no interpolation
971  // implemented. so throw exception,
972  // as documentation says
973  AssertThrow(
974  false,
976 }
977 
978 
979 
980 template <int dim, int spacedim>
981 std::vector<std::pair<unsigned int, unsigned int>>
983  const FiniteElement<dim, spacedim> &) const
984 {
985  Assert(false, ExcNotImplemented());
986  return std::vector<std::pair<unsigned int, unsigned int>>();
987 }
988 
989 
990 
991 template <int dim, int spacedim>
992 std::vector<std::pair<unsigned int, unsigned int>>
994  const FiniteElement<dim, spacedim> &) const
995 {
996  Assert(false, ExcNotImplemented());
997  return std::vector<std::pair<unsigned int, unsigned int>>();
998 }
999 
1000 
1001 
1002 template <int dim, int spacedim>
1003 std::vector<std::pair<unsigned int, unsigned int>>
1006  const unsigned int) const
1007 {
1008  Assert(false, ExcNotImplemented());
1009  return std::vector<std::pair<unsigned int, unsigned int>>();
1010 }
1011 
1012 
1013 
1014 template <int dim, int spacedim>
1018  const unsigned int) const
1019 {
1020  Assert(false, ExcNotImplemented());
1022 }
1023 
1024 
1025 
1026 template <int dim, int spacedim>
1027 bool
1030 {
1031  // Compare fields in roughly increasing order of how expensive the
1032  // comparison is
1033  return ((typeid(*this) == typeid(f)) && (this->get_name() == f.get_name()) &&
1034  (static_cast<const FiniteElementData<dim> &>(*this) ==
1035  static_cast<const FiniteElementData<dim> &>(f)) &&
1037 }
1038 
1039 
1040 
1041 template <int dim, int spacedim>
1042 bool
1045 {
1046  return !(*this == f);
1047 }
1048 
1049 
1050 
1051 template <int dim, int spacedim>
1052 const std::vector<Point<dim>> &
1054 {
1055  // a finite element may define
1056  // support points, but only if
1057  // there are as many as there are
1058  // degrees of freedom
1059  Assert((unit_support_points.size() == 0) ||
1060  (unit_support_points.size() == this->n_dofs_per_cell()),
1061  ExcInternalError());
1062  return unit_support_points;
1063 }
1064 
1065 
1066 
1067 template <int dim, int spacedim>
1068 bool
1070 {
1071  return (unit_support_points.size() != 0);
1072 }
1073 
1074 
1075 
1076 template <int dim, int spacedim>
1077 const std::vector<Point<dim>> &
1079 {
1080  // If the finite element implements generalized support points, return
1081  // those. Otherwise fall back to unit support points.
1082  return ((generalized_support_points.size() == 0) ?
1085 }
1086 
1087 
1088 
1089 template <int dim, int spacedim>
1090 bool
1092 {
1093  return (get_generalized_support_points().size() != 0);
1094 }
1095 
1096 
1097 
1098 template <int dim, int spacedim>
1099 Point<dim>
1101 {
1102  AssertIndexRange(index, this->n_dofs_per_cell());
1103  Assert(unit_support_points.size() == this->n_dofs_per_cell(),
1105  return unit_support_points[index];
1106 }
1107 
1108 
1109 
1110 template <int dim, int spacedim>
1111 const std::vector<Point<dim - 1>> &
1113  const unsigned int face_no) const
1114 {
1115  // a finite element may define
1116  // support points, but only if
1117  // there are as many as there are
1118  // degrees of freedom on a face
1119  Assert((unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1120  .size() == 0) ||
1121  (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1122  .size() == this->n_dofs_per_face(face_no)),
1123  ExcInternalError());
1124  return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no];
1125 }
1126 
1127 
1128 
1129 template <int dim, int spacedim>
1130 bool
1132  const unsigned int face_no) const
1133 {
1134  return (unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1135  .size() != 0);
1136 }
1137 
1138 
1139 
1140 template <int dim, int spacedim>
1141 Point<dim - 1>
1143  const unsigned int index,
1144  const unsigned int face_no) const
1145 {
1146  AssertIndexRange(index, this->n_dofs_per_face(face_no));
1147  Assert(unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1148  .size() == this->n_dofs_per_face(face_no),
1150  return unit_face_support_points[this->n_unique_faces() == 1 ? 0 : face_no]
1151  [index];
1152 }
1153 
1154 
1155 
1156 template <int dim, int spacedim>
1157 bool
1159  const unsigned int) const
1160 {
1161  return true;
1162 }
1163 
1164 
1165 
1166 template <int dim, int spacedim>
1169 {
1170  // Translate the ComponentMask into first_selected and n_components after
1171  // some error checking:
1172  const unsigned int n_total_components = this->n_components();
1173  Assert((n_total_components == mask.size()) || (mask.size() == 0),
1174  ExcMessage("The given ComponentMask has the wrong size."));
1175 
1176  const unsigned int n_selected =
1177  mask.n_selected_components(n_total_components);
1178  Assert(n_selected > 0,
1179  ExcMessage("You need at least one selected component."));
1180 
1181  const unsigned int first_selected =
1182  mask.first_selected_component(n_total_components);
1183 
1184 #ifdef DEBUG
1185  // check that it is contiguous:
1186  for (unsigned int c = 0; c < n_total_components; ++c)
1187  Assert((c < first_selected && (!mask[c])) ||
1188  (c >= first_selected && c < first_selected + n_selected &&
1189  mask[c]) ||
1190  (c >= first_selected + n_selected && !mask[c]),
1191  ExcMessage("Error: the given ComponentMask is not contiguous!"));
1192 #endif
1193 
1194  return get_sub_fe(first_selected, n_selected);
1195 }
1196 
1197 
1198 
1199 template <int dim, int spacedim>
1202  const unsigned int first_component,
1203  const unsigned int n_selected_components) const
1204 {
1205  // No complicated logic is needed here, because it is overridden in
1206  // FESystem<dim,spacedim>. Just make sure that what the user chose is valid:
1207  Assert(first_component == 0 && n_selected_components == this->n_components(),
1208  ExcMessage(
1209  "You can only select a whole FiniteElement, not a part of one."));
1210 
1211  (void)first_component;
1212  (void)n_selected_components;
1213  return *this;
1214 }
1215 
1216 
1217 
1218 template <int dim, int spacedim>
1219 std::pair<Table<2, bool>, std::vector<unsigned int>>
1221 {
1222  Assert(false, ExcNotImplemented());
1223  return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1224  Table<2, bool>(this->n_components(), this->n_dofs_per_cell()),
1225  std::vector<unsigned int>(this->n_components()));
1226 }
1227 
1228 
1229 
1230 template <int dim, int spacedim>
1231 void
1234  const std::vector<Vector<double>> &,
1235  std::vector<double> &) const
1236 {
1238  ExcMessage("The element for which you are calling the current "
1239  "function does not have generalized support points (see "
1240  "the glossary for a definition of generalized support "
1241  "points). Consequently, the current function can not "
1242  "be defined and is not implemented by the element."));
1243  Assert(false, ExcNotImplemented());
1244 }
1245 
1246 
1247 
1248 template <int dim, int spacedim>
1249 std::size_t
1251 {
1252  return (
1253  sizeof(FiniteElementData<dim>) +
1265 }
1266 
1267 
1268 
1269 template <int dim, int spacedim>
1270 std::vector<unsigned int>
1272  const std::vector<ComponentMask> &nonzero_components)
1273 {
1274  std::vector<unsigned int> retval(nonzero_components.size());
1275  for (unsigned int i = 0; i < nonzero_components.size(); ++i)
1276  retval[i] = nonzero_components[i].n_selected_components();
1277  return retval;
1278 }
1279 
1280 
1281 
1282 /*------------------------------- FiniteElement ----------------------*/
1283 
1284 template <int dim, int spacedim>
1285 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1287  const UpdateFlags flags,
1288  const Mapping<dim, spacedim> &mapping,
1289  const Quadrature<dim - 1> & quadrature,
1291  spacedim>
1292  &output_data) const
1293 {
1294  return get_data(flags,
1295  mapping,
1297  this->reference_cell_type(), quadrature),
1298  output_data);
1299 }
1300 
1301 
1302 
1303 template <int dim, int spacedim>
1304 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1306  const UpdateFlags flags,
1307  const Mapping<dim, spacedim> &mapping,
1308  const Quadrature<dim - 1> & quadrature,
1310  spacedim>
1311  &output_data) const
1312 {
1313  return get_data(flags,
1314  mapping,
1316  this->reference_cell_type(), quadrature),
1317  output_data);
1318 }
1319 
1320 
1321 
1322 template <int dim, int spacedim>
1324 FiniteElement<dim, spacedim>::base_element(const unsigned int index) const
1325 {
1326  (void)index;
1327  AssertIndexRange(index, 1);
1328  // This function should not be
1329  // called for a system element
1331  return *this;
1332 }
1333 
1334 
1335 
1336 /*------------------------------- Explicit Instantiations -------------*/
1337 #include "fe.inst"
1338 
1339 
static ::ExceptionBase & ExcFEHasNoSupportPoints()
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > face_system_to_component_table
Definition: fe.h:2510
bool represents_the_all_selected_mask() const
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:196
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:482
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:243
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2404
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation) const
Definition: fe.cc:660
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1623
virtual std::size_t memory_consumption() const
Definition: fe.cc:49
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
Definition: fe.cc:1029
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2455
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:256
unsigned int get_first_line_index() const
FullMatrix< double > interface_constraints
Definition: fe.h:2430
bool operator!=(const FiniteElement< dim, spacedim > &) const
Definition: fe.cc:1044
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:57
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
Definition: fe.cc:1168
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:197
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:931
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
virtual Point< dim > unit_support_point(const unsigned int index) const
Definition: fe.cc:1100
bool restriction_is_implemented() const
Definition: fe.cc:757
const std::vector< Point< dim > > & get_generalized_support_points() const
Definition: fe.cc:1078
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:209
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:267
const ReferenceCell::internal::Info::Base & get_cell(const ReferenceCell::Type &type)
bool has_face_support_points(const unsigned int face_no=0) const
Definition: fe.cc:1131
#define AssertIndexRange(index, range)
Definition: exceptions.h:1691
bool isotropic_prolongation_is_implemented() const
Definition: fe.cc:785
unsigned int n_dofs_per_vertex() const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition: fe.h:2565
bool has_generalized_support_points() const
Definition: fe.cc:1091
bool represents_the_all_selected_mask() const
unsigned int n_unique_quads() const
STL namespace.
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:396
static const char U
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const =0
#define AssertThrow(cond, exc)
Definition: exceptions.h:1576
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:1271
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:358
bool is_primitive() const
Definition: fe.h:3275
const std::vector< unsigned int > n_nonzero_components_table
Definition: fe.h:2590
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
unsigned int n_blocks() const
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2442
ReferenceCell::Type reference_cell_type() const
bool has_support_points() const
Definition: fe.cc:1069
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
Definition: fe.cc:1305
unsigned int component_to_block_index(const unsigned int component) const
Definition: fe.cc:384
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index, const unsigned int face_no=0) const
Definition: fe.cc:1142
static ::ExceptionBase & ExcMessage(std::string arg1)
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
Definition: fe.cc:869
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:220
virtual std::size_t memory_consumption() const
Definition: fe.cc:1250
size_type n() const
unsigned int n_dofs_per_line() const
No update.
virtual ReferenceCell::Type face_reference_cell_type(const unsigned int face_no) const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2418
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
Definition: fe.cc:947
#define Assert(cond, exc)
Definition: exceptions.h:1466
UpdateFlags
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:335
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: mapping.h:301
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
Definition: fe.cc:1112
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
Definition: fe.h:2449
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:372
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
Definition: fe.cc:841
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
VectorType::value_type * end(VectorType &V)
unsigned int size() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1324
virtual std::string get_name() const =0
static ::ExceptionBase & ExcProjectionVoid()
std::vector< std::vector< Point< dim - 1 > > > generalized_face_support_points
Definition: fe.h:2461
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:982
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe.cc:568
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_unique_faces() 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
Definition: fe.cc:1233
Point< 2 > first
Definition: grid_out.cc:4340
types::global_dof_index first_block_of_base(const unsigned int b) const
Definition: fe.h:3195
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const bool line_orientation) const
Definition: fe.cc:704
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:1053
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
Definition: fe.cc:964
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe.cc:1286
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:993
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
Definition: fe.cc:1016
static ::ExceptionBase & ExcEmbeddingVoid()
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:186
unsigned int n_components() const
unsigned int n_dofs_per_cell() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:371
VectorType::value_type * begin(VectorType &V)
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
Definition: fe.h:2498
virtual bool hp_constraints_are_implemented() const
Definition: fe.cc:860
std::vector< std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > > face_system_to_base_table
Definition: fe.h:2536
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
bool prolongation_is_implemented() const
Definition: fe.cc:729
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:291
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:232
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2581
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:303
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition: fe.h:2529
static ::ExceptionBase & ExcNotImplemented()
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
Definition: fe.cc:1004
bool isotropic_restriction_is_implemented() const
Definition: fe.cc:813
BlockIndices base_to_block_indices
Definition: fe.h:2542
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:280
unsigned int size() const
unsigned int size() const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition: fe.h:2493
TableIndices< 2 > interface_constraints_size() const
Definition: fe.cc:903
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe.cc:1220
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2572
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
Definition: fe.cc:177
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe.cc:1158
virtual std::unique_ptr< InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
std::vector< Table< 2, int > > adjust_quad_dof_index_for_face_orientation_table
Definition: fe.h:2478
static ::ExceptionBase & ExcInternalError()