Reference documentation for deal.II version 9.0.0
fe.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/memory_consumption.h>
17 #include <deal.II/fe/mapping.h>
18 #include <deal.II/fe/fe.h>
19 #include <deal.II/fe/fe_values.h>
20 #include <deal.II/base/quadrature.h>
21 #include <deal.II/base/qprojector.h>
22 #include <deal.II/grid/tria.h>
23 #include <deal.II/grid/tria_iterator.h>
24 #include <deal.II/dofs/dof_accessor.h>
25 #include <deal.II/grid/tria_boundary.h>
26 
27 #include <algorithm>
28 #include <functional>
29 #include <numeric>
30 #include <typeinfo>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
35 /*------------------------------- FiniteElement ----------------------*/
36 
37 
38 template <int dim, int spacedim>
40  update_each(update_default)
41 {}
42 
43 
44 
45 template <int dim, int spacedim>
46 std::size_t
48 {
49  return sizeof(*this);
50 }
51 
52 
53 
54 template <int dim, int spacedim>
57  const std::vector<bool> &r_i_a_f,
58  const std::vector<ComponentMask> &nonzero_c)
59  :
60  FiniteElementData<dim> (fe_data),
62  this->dofs_per_quad : 0,
63  dim==3 ? 8 : 0),
65  this->dofs_per_line : 0),
69  std::make_pair(std::make_pair(0U, 0U), 0U)),
70 
71  // Special handling of vectors of length one: in this case, we
72  // assume that all entries were supposed to be equal
73  restriction_is_additive_flags(r_i_a_f.size() == 1
74  ?
75  std::vector<bool> (fe_data.dofs_per_cell, r_i_a_f[0])
76  :
77  r_i_a_f),
78  nonzero_components (nonzero_c.size() == 1
79  ?
80  std::vector<ComponentMask> (fe_data.dofs_per_cell, nonzero_c[0])
81  :
82  nonzero_c),
86  std::bind (std::not_equal_to<unsigned int>(),
87  std::placeholders::_1,
88  1U))
90 {
93  this->dofs_per_cell));
95  for (unsigned int i=0; i<nonzero_components.size(); ++i)
96  {
97  Assert (nonzero_components[i].size() == this->n_components(),
99  Assert (nonzero_components[i].n_selected_components ()
100  >= 1,
101  ExcInternalError());
103  ExcInternalError());
105  ExcInternalError());
106  }
107 
108  // initialize some tables in the default way, i.e. if there is only one
109  // (vector-)component; if the element is not primitive, leave these tables
110  // empty.
111  if (this->is_primitive())
112  {
115  for (unsigned int j=0 ; j<this->dofs_per_cell ; ++j)
116  system_to_component_table[j] = std::pair<unsigned,unsigned>(0,j);
117  for (unsigned int j=0 ; j<this->dofs_per_face ; ++j)
118  face_system_to_component_table[j] = std::pair<unsigned,unsigned>(0,j);
119  }
120 
121  for (unsigned int j=0 ; j<this->dofs_per_cell ; ++j)
122  system_to_base_table[j] = std::make_pair(std::make_pair(0U,0U),j);
123  for (unsigned int j=0 ; j<this->dofs_per_face ; ++j)
124  face_system_to_base_table[j] = std::make_pair(std::make_pair(0U,0U),j);
125 
126  // Fill with default value; may be changed by constructor of derived class.
128 
129  // initialize the restriction and prolongation matrices. the default
130  // constructor of FullMatrix<dim> initializes them with size zero
133  for (unsigned int ref=RefinementCase<dim>::cut_x;
134  ref<RefinementCase<dim>::isotropic_refinement+1; ++ref)
135  {
136  prolongation[ref-1].resize (GeometryInfo<dim>::
137  n_children(RefinementCase<dim>(ref)),
139  restriction[ref-1].resize (GeometryInfo<dim>::
140  n_children(RefinementCase<dim>(ref)),
142  }
143 
145 }
146 
147 
148 
149 template <int dim, int spacedim>
150 std::pair<std::unique_ptr<FiniteElement<dim, spacedim> >, unsigned int>
151 FiniteElement<dim, spacedim>::operator^ (const unsigned int multiplicity) const
152 {
153  return {this->clone(), multiplicity};
154 }
155 
156 
157 
158 template <int dim, int spacedim>
159 double
161  const Point<dim> &) const
162 {
163  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
164  return 0.;
165 }
166 
167 
168 
169 template <int dim, int spacedim>
170 double
172  const Point<dim> &,
173  const unsigned int) const
174 {
175  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
176  return 0.;
177 }
178 
179 
180 
181 template <int dim, int spacedim>
184  const Point<dim> &) const
185 {
186  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
187  return Tensor<1,dim> ();
188 }
189 
190 
191 
192 template <int dim, int spacedim>
195  const Point<dim> &,
196  const unsigned int) const
197 {
198  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
199  return Tensor<1,dim> ();
200 }
201 
202 
203 
204 template <int dim, int spacedim>
207  const Point<dim> &) const
208 {
209  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
210  return Tensor<2,dim> ();
211 }
212 
213 
214 
215 template <int dim, int spacedim>
218  const Point<dim> &,
219  const unsigned int) const
220 {
221  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
222  return Tensor<2,dim> ();
223 }
224 
225 
226 
227 template <int dim, int spacedim>
230  const Point<dim> &) const
231 {
232  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
233  return Tensor<3,dim> ();
234 }
235 
236 
237 
238 template <int dim, int spacedim>
241  const Point<dim> &,
242  const unsigned int) const
243 {
244  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
245  return Tensor<3,dim> ();
246 }
247 
248 
249 
250 template <int dim, int spacedim>
253  const Point<dim> &) const
254 {
255  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
256  return Tensor<4,dim> ();
257 }
258 
259 
260 
261 template <int dim, int spacedim>
264  const Point<dim> &,
265  const unsigned int) const
266 {
267  AssertThrow(false, ExcUnitShapeValuesDoNotExist());
268  return Tensor<4,dim> ();
269 }
270 
271 
272 template <int dim, int spacedim>
273 void
275  const bool isotropic_restriction_only,
276  const bool isotropic_prolongation_only)
277 {
278  for (unsigned int ref_case=RefinementCase<dim>::cut_x;
279  ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
280  {
281  const unsigned int nc = GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
282 
283  for (unsigned int i=0; i<nc; ++i)
284  {
285  if (this->restriction[ref_case-1][i].m() != this->dofs_per_cell
286  &&
287  (!isotropic_restriction_only || ref_case==RefinementCase<dim>::isotropic_refinement))
288  this->restriction[ref_case-1][i].reinit (this->dofs_per_cell,
289  this->dofs_per_cell);
290  if (this->prolongation[ref_case-1][i].m() != this->dofs_per_cell
291  &&
292  (!isotropic_prolongation_only || ref_case==RefinementCase<dim>::isotropic_refinement))
293  this->prolongation[ref_case-1][i].reinit (this->dofs_per_cell,
294  this->dofs_per_cell);
295  }
296  }
297 }
298 
299 
300 template <int dim, int spacedim>
301 const FullMatrix<double> &
303  const RefinementCase<dim> &refinement_case) const
304 {
307  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
308  ExcMessage("Restriction matrices are only available for refined cells!"));
311  // we use refinement_case-1 here. the -1 takes care of the origin of the
312  // vector, as for RefinementCase<dim>::no_refinement (=0) there is no data
313  // available and so the vector indices are shifted
314  Assert (restriction[refinement_case-1][child].n() == this->dofs_per_cell, ExcProjectionVoid());
315  return restriction[refinement_case-1][child];
316 }
317 
318 
319 
320 template <int dim, int spacedim>
321 const FullMatrix<double> &
323  const RefinementCase<dim> &refinement_case) const
324 {
327  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
328  ExcMessage("Prolongation matrices are only available for refined cells!"));
331  // we use refinement_case-1 here. the -1 takes care
332  // of the origin of the vector, as for
333  // RefinementCase::no_refinement (=0) there is no
334  // data available and so the vector indices
335  // are shifted
336  Assert (prolongation[refinement_case-1][child].n() == this->dofs_per_cell, ExcEmbeddingVoid());
337  return prolongation[refinement_case-1][child];
338 }
339 
340 
341 //TODO:[GK] This is probably not the most efficient way of doing this.
342 template <int dim, int spacedim>
343 unsigned int
345 {
346  Assert (index < this->n_components(),
347  ExcIndexRange(index, 0, this->n_components()));
348 
349  return first_block_of_base(component_to_base_table[index].first.first)
350  + component_to_base_table[index].second;
351 }
352 
353 
354 template <int dim, int spacedim>
358 {
359  AssertIndexRange(scalar.component, this->n_components());
360 
361 //TODO: it would be nice to verify that it is indeed possible
362 // to select this scalar component, i.e., that it is not part
363 // of a non-primitive element. unfortunately, there is no simple
364 // way to write such a condition...
365 
366  std::vector<bool> mask (this->n_components(), false);
367  mask[scalar.component] = true;
368  return mask;
369 }
370 
371 
372 template <int dim, int spacedim>
376 {
377  AssertIndexRange(vector.first_vector_component+dim-1, this->n_components());
378 
379  //TODO: it would be nice to verify that it is indeed possible
380  // to select these vector components, i.e., that they don't span
381  // beyond the beginning or end of anon-primitive element.
382  // unfortunately, there is no simple way to write such a condition...
383 
384  std::vector<bool> mask (this->n_components(), false);
385  for (unsigned int c=vector.first_vector_component; c<vector.first_vector_component+dim; ++c)
386  mask[c] = true;
387  return mask;
388 }
389 
390 
391 template <int dim, int spacedim>
395 {
398  this->n_components());
399 
400  //TODO: it would be nice to verify that it is indeed possible
401  // to select these vector components, i.e., that they don't span
402  // beyond the beginning or end of anon-primitive element.
403  // unfortunately, there is no simple way to write such a condition...
404 
405  std::vector<bool> mask (this->n_components(), false);
406  for (unsigned int c=sym_tensor.first_tensor_component;
408  mask[c] = true;
409  return mask;
410 }
411 
412 
413 
414 template <int dim, int spacedim>
417 component_mask (const BlockMask &block_mask) const
418 {
419  // if we get a block mask that represents all blocks, then
420  // do the same for the returned component mask
421  if (block_mask.represents_the_all_selected_mask())
422  return ComponentMask();
423 
424  AssertDimension(block_mask.size(), this->n_blocks());
425 
426  std::vector<bool> component_mask (this->n_components(), false);
427  for (unsigned int c=0; c<this->n_components(); ++c)
428  if (block_mask[component_to_block_index(c)] == true)
429  component_mask[c] = true;
430 
431  return component_mask;
432 }
433 
434 
435 
436 template <int dim, int spacedim>
437 BlockMask
440 {
441  // simply create the corresponding component mask (a simpler
442  // process) and then convert it to a block mask
443  return block_mask(component_mask(scalar));
444 }
445 
446 
447 template <int dim, int spacedim>
448 BlockMask
451 {
452  // simply create the corresponding component mask (a simpler
453  // process) and then convert it to a block mask
454  return block_mask(component_mask(vector));
455 }
456 
457 
458 template <int dim, int spacedim>
459 BlockMask
462 {
463  // simply create the corresponding component mask (a simpler
464  // process) and then convert it to a block mask
465  return block_mask(component_mask(sym_tensor));
466 }
467 
468 
469 
470 template <int dim, int spacedim>
471 BlockMask
473 block_mask (const ComponentMask &component_mask) const
474 {
475  // if we get a component mask that represents all component, then
476  // do the same for the returned block mask
477  if (component_mask.represents_the_all_selected_mask())
478  return BlockMask();
479 
480  AssertDimension(component_mask.size(), this->n_components());
481 
482  // walk over all of the components
483  // of this finite element and see
484  // if we need to set the
485  // corresponding block. inside the
486  // block, walk over all the
487  // components that correspond to
488  // this block and make sure the
489  // component mask is set for all of
490  // them
491  std::vector<bool> block_mask (this->n_blocks(), false);
492  for (unsigned int c=0; c<this->n_components();)
493  {
494  const unsigned int block = component_to_block_index(c);
495  if (component_mask[c] == true)
496  block_mask[block] = true;
497 
498  // now check all of the other
499  // components that correspond
500  // to this block
501  ++c;
502  while ((c<this->n_components())
503  &&
504  (component_to_block_index(c) == block))
505  {
506  Assert (component_mask[c] == block_mask[block],
507  ExcMessage ("The component mask argument given to this function "
508  "is not a mask where the individual components belonging "
509  "to one block of the finite element are either all "
510  "selected or not selected. You can't call this function "
511  "with a component mask that splits blocks."));
512  ++c;
513  }
514  }
515 
516 
517  return block_mask;
518 }
519 
520 
521 
522 template <int dim, int spacedim>
523 unsigned int
525 face_to_cell_index (const unsigned int face_index,
526  const unsigned int face,
527  const bool face_orientation,
528  const bool face_flip,
529  const bool face_rotation) const
530 {
531  Assert (face_index < this->dofs_per_face,
532  ExcIndexRange(face_index, 0, this->dofs_per_face));
535 
536 //TODO: we could presumably solve the 3d case below using the
537 // adjust_quad_dof_index_for_face_orientation_table field. for the
538 // 2d case, we can't use adjust_line_dof_index_for_line_orientation_table
539 // since that array is empty (presumably because we thought that
540 // there are no flipped edges in 2d, but these can happen in
541 // DoFTools::make_periodicity_constraints, for example). so we
542 // would need to either fill this field, or rely on derived classes
543 // implementing this function, as we currently do
544 
545  // see the function's documentation for an explanation of this
546  // assertion -- in essence, derived classes have to implement
547  // an overloaded version of this function if we are to use any
548  // other than standard orientation
549  if ((face_orientation != true) || (face_flip != false) || (face_rotation != false))
550  Assert ((this->dofs_per_line <= 1) && (this->dofs_per_quad <= 1),
551  ExcMessage ("The function in this base class can not handle this case. "
552  "Rather, the derived class you are using must provide "
553  "an overloaded version but apparently hasn't done so. See "
554  "the documentation of this function for more information."));
555 
556  // we need to distinguish between DoFs on vertices, lines and in 3d quads.
557  // do so in a sequence of if-else statements
558  if (face_index < this->first_face_line_index)
559  // DoF is on a vertex
560  {
561  // get the number of the vertex on the face that corresponds to this DoF,
562  // along with the number of the DoF on this vertex
563  const unsigned int face_vertex = face_index / this->dofs_per_vertex;
564  const unsigned int dof_index_on_vertex = face_index % this->dofs_per_vertex;
565 
566  // then get the number of this vertex on the cell and translate
567  // this to a DoF number on the cell
568  return (GeometryInfo<dim>::face_to_cell_vertices(face, face_vertex,
569  face_orientation,
570  face_flip,
571  face_rotation)
572  * this->dofs_per_vertex
573  +
574  dof_index_on_vertex);
575  }
576  else if (face_index < this->first_face_quad_index)
577  // DoF is on a face
578  {
579  // do the same kind of translation as before. we need to only consider
580  // DoFs on the lines, i.e., ignoring those on the vertices
581  const unsigned int index = face_index - this->first_face_line_index;
582 
583  const unsigned int face_line = index / this->dofs_per_line;
584  const unsigned int dof_index_on_line = index % this->dofs_per_line;
585 
586  return (this->first_line_index
587  + GeometryInfo<dim>::face_to_cell_lines(face, face_line,
588  face_orientation,
589  face_flip,
590  face_rotation)
591  * this->dofs_per_line
592  +
593  dof_index_on_line);
594  }
595  else
596  // DoF is on a quad
597  {
598  Assert (dim >= 3, ExcInternalError());
599 
600  // ignore vertex and line dofs
601  const unsigned int index = face_index - this->first_face_quad_index;
602 
603  return (this->first_quad_index
604  + face * this->dofs_per_quad
605  + index);
606  }
607 }
608 
609 
610 
611 
612 template <int dim, int spacedim>
613 unsigned int
615  const bool face_orientation,
616  const bool face_flip,
617  const bool face_rotation) const
618 {
619  // general template for 1D and 2D: not
620  // implemented. in fact, the function
621  // shouldn't even be called unless we are
622  // in 3d, so throw an internal error
623  Assert (dim==3, ExcInternalError());
624  if (dim < 3)
625  return index;
626 
627  // adjust dofs on 3d faces if the face is
628  // flipped. note that we query a table that
629  // derived elements need to have set up
630  // front. the exception are discontinuous
631  // elements for which there should be no
632  // face dofs anyway (i.e. dofs_per_quad==0
633  // in 3d), so we don't need the table, but
634  // the function should also not have been
635  // called
636  Assert (index<this->dofs_per_quad, ExcIndexRange(index,0,this->dofs_per_quad));
637  Assert (adjust_quad_dof_index_for_face_orientation_table.n_elements()==8*this->dofs_per_quad,
638  ExcInternalError());
639  return index+adjust_quad_dof_index_for_face_orientation_table(index,4*face_orientation+2*face_flip+face_rotation);
640 }
641 
642 
643 
644 template <int dim, int spacedim>
645 unsigned int
647  const bool line_orientation) const
648 {
649  // general template for 1D and 2D: do
650  // nothing. Do not throw an Assertion,
651  // however, in order to allow to call this
652  // function in 2D as well
653  if (dim<3)
654  return index;
655 
656  Assert (index<this->dofs_per_line, ExcIndexRange(index,0,this->dofs_per_line));
657  Assert (adjust_line_dof_index_for_line_orientation_table.size()==this->dofs_per_line,
658  ExcInternalError());
659  if (line_orientation)
660  return index;
661  else
662  return index+adjust_line_dof_index_for_line_orientation_table[index];
663 }
664 
665 
666 
667 template <int dim, int spacedim>
668 bool
670 {
671  for (unsigned int ref_case=RefinementCase<dim>::cut_x;
672  ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
673  for (unsigned int c=0;
674  c<GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case)); ++c)
675  {
676  // make sure also the lazily initialized matrices are created
677  get_prolongation_matrix(c, RefinementCase<dim>(ref_case));
678  Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
679  (prolongation[ref_case-1][c].m() == 0),
680  ExcInternalError());
681  Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
682  (prolongation[ref_case-1][c].n() == 0),
683  ExcInternalError());
684  if ((prolongation[ref_case-1][c].m() == 0) ||
685  (prolongation[ref_case-1][c].n() == 0))
686  return false;
687  }
688  return true;
689 }
690 
691 
692 
693 template <int dim, int spacedim>
694 bool
696 {
697  for (unsigned int ref_case=RefinementCase<dim>::cut_x;
698  ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
699  for (unsigned int c=0;
700  c<GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case)); ++c)
701  {
702  // make sure also the lazily initialized matrices are created
703  get_restriction_matrix(c, RefinementCase<dim>(ref_case));
704  Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
705  (restriction[ref_case-1][c].m() == 0),
706  ExcInternalError());
707  Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
708  (restriction[ref_case-1][c].n() == 0),
709  ExcInternalError());
710  if ((restriction[ref_case-1][c].m() == 0) ||
711  (restriction[ref_case-1][c].n() == 0))
712  return false;
713  }
714  return true;
715 }
716 
717 
718 
719 template <int dim, int spacedim>
720 bool
722 {
724 
725  for (unsigned int c=0;
726  c<GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case)); ++c)
727  {
728  // make sure also the lazily initialized matrices are created
729  get_prolongation_matrix(c, RefinementCase<dim>(ref_case));
730  Assert ((prolongation[ref_case-1][c].m() == this->dofs_per_cell) ||
731  (prolongation[ref_case-1][c].m() == 0),
732  ExcInternalError());
733  Assert ((prolongation[ref_case-1][c].n() == this->dofs_per_cell) ||
734  (prolongation[ref_case-1][c].n() == 0),
735  ExcInternalError());
736  if ((prolongation[ref_case-1][c].m() == 0) ||
737  (prolongation[ref_case-1][c].n() == 0))
738  return false;
739  }
740  return true;
741 }
742 
743 
744 
745 template <int dim, int spacedim>
746 bool
748 {
750 
751  for (unsigned int c=0;
752  c<GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case)); ++c)
753  {
754  // make sure also the lazily initialized matrices are created
755  get_restriction_matrix(c, RefinementCase<dim>(ref_case));
756  Assert ((restriction[ref_case-1][c].m() == this->dofs_per_cell) ||
757  (restriction[ref_case-1][c].m() == 0),
758  ExcInternalError());
759  Assert ((restriction[ref_case-1][c].n() == this->dofs_per_cell) ||
760  (restriction[ref_case-1][c].n() == 0),
761  ExcInternalError());
762  if ((restriction[ref_case-1][c].m() == 0) ||
763  (restriction[ref_case-1][c].n() == 0))
764  return false;
765  }
766  return true;
767 }
768 
769 
770 
771 template <int dim, int spacedim>
772 bool
774 {
776  return (this->dofs_per_face == 0) || (interface_constraints.m() != 0);
777  else
778  return false;
779 }
780 
781 
782 
783 template <int dim, int spacedim>
784 bool
786 {
787  return false;
788 }
789 
790 
791 
792 template <int dim, int spacedim>
793 const FullMatrix<double> &
795 {
796  (void)subface_case;
798  ExcMessage("Constraints for this element are only implemented "
799  "for the case that faces are refined isotropically "
800  "(which is always the case in 2d, and in 3d requires "
801  "that the neighboring cell of a coarse cell presents "
802  "exactly four children on the common face)."));
803  Assert ((this->dofs_per_face == 0) || (interface_constraints.m() != 0),
804  ExcMessage ("The finite element for which you try to obtain "
805  "hanging node constraints does not appear to "
806  "implement them."));
807 
808  if (dim==1)
809  Assert ((interface_constraints.m()==0) && (interface_constraints.n()==0),
810  ExcWrongInterfaceMatrixSize(interface_constraints.m(),
811  interface_constraints.n()));
812 
813  return interface_constraints;
814 }
815 
816 
817 
818 template <int dim, int spacedim>
821 {
822  switch (dim)
823  {
824  case 1:
825  return TableIndices<2> (0U, 0U);
826  case 2:
827  return TableIndices<2> (this->dofs_per_vertex +
828  2*this->dofs_per_line,
829  this->dofs_per_face);
830  case 3:
831  return TableIndices<2> (5*this->dofs_per_vertex +
832  12*this->dofs_per_line +
833  4*this->dofs_per_quad,
834  this->dofs_per_face);
835  default:
836  Assert (false, ExcNotImplemented());
837  };
840 }
841 
842 
843 
844 template <int dim, int spacedim>
845 void
848  FullMatrix<double> &) const
849 {
850  // by default, no interpolation
851  // implemented. so throw exception,
852  // as documentation says
853  AssertThrow (false,
854  (typename FiniteElement<dim,spacedim>::
855  ExcInterpolationNotImplemented()));
856 }
857 
858 
859 
860 template <int dim, int spacedim>
861 void
864  FullMatrix<double> &) const
865 {
866  // by default, no interpolation
867  // implemented. so throw exception,
868  // as documentation says
869  AssertThrow (false,
870  (typename FiniteElement<dim,spacedim>::
871  ExcInterpolationNotImplemented()));
872 }
873 
874 
875 
876 template <int dim, int spacedim>
877 void
880  const unsigned int,
881  FullMatrix<double> &) const
882 {
883  // by default, no interpolation
884  // implemented. so throw exception,
885  // as documentation says
886  AssertThrow (false,
887  (typename FiniteElement<dim,spacedim>::
888  ExcInterpolationNotImplemented()));
889 }
890 
891 
892 
893 template <int dim, int spacedim>
894 std::vector<std::pair<unsigned int, unsigned int> >
897 {
898  Assert (false, ExcNotImplemented());
899  return std::vector<std::pair<unsigned int, unsigned int> > ();
900 }
901 
902 
903 
904 template <int dim, int spacedim>
905 std::vector<std::pair<unsigned int, unsigned int> >
908 {
909  Assert (false, ExcNotImplemented());
910  return std::vector<std::pair<unsigned int, unsigned int> > ();
911 }
912 
913 
914 
915 template <int dim, int spacedim>
916 std::vector<std::pair<unsigned int, unsigned int> >
919 {
920  Assert (false, ExcNotImplemented());
921  return std::vector<std::pair<unsigned int, unsigned int> > ();
922 }
923 
924 
925 
926 template <int dim, int spacedim>
930 {
931  Assert (false, ExcNotImplemented());
933 }
934 
935 
936 
937 template <int dim, int spacedim>
938 bool
940 {
941  // Compare fields in roughly increasing order of how expensive the
942  // comparison is
943  return ((typeid(*this) == typeid(f))
944  &&
945  (this->get_name() == f.get_name())
946  &&
947  (static_cast<const FiniteElementData<dim>&>(*this) ==
948  static_cast<const FiniteElementData<dim>&>(f))
949  &&
950  (interface_constraints == f.interface_constraints));
951 }
952 
953 
954 
955 template <int dim, int spacedim>
956 bool
958 {
959  return !(*this == f);
960 }
961 
962 
963 
964 template <int dim, int spacedim>
965 const std::vector<Point<dim> > &
967 {
968  // a finite element may define
969  // support points, but only if
970  // there are as many as there are
971  // degrees of freedom
972  Assert ((unit_support_points.size() == 0) ||
973  (unit_support_points.size() == this->dofs_per_cell),
974  ExcInternalError());
975  return unit_support_points;
976 }
977 
978 
979 
980 template <int dim, int spacedim>
981 bool
983 {
984  return (unit_support_points.size() != 0);
985 }
986 
987 
988 
989 template <int dim, int spacedim>
990 const std::vector<Point<dim> > &
992 {
993  // If the finite element implements generalized support points, return
994  // those. Otherwise fall back to unit support points.
995  return ((generalized_support_points.size() == 0)
996  ? unit_support_points
997  : generalized_support_points);
998 }
999 
1000 
1001 
1002 template <int dim, int spacedim>
1003 bool
1005 {
1006  return (get_generalized_support_points().size() != 0);
1007 }
1008 
1009 
1010 
1011 template <int dim, int spacedim>
1012 Point<dim>
1013 FiniteElement<dim,spacedim>::unit_support_point (const unsigned int index) const
1014 {
1015  Assert (index < this->dofs_per_cell,
1016  ExcIndexRange (index, 0, this->dofs_per_cell));
1017  Assert (unit_support_points.size() == this->dofs_per_cell,
1018  ExcFEHasNoSupportPoints ());
1019  return unit_support_points[index];
1020 }
1021 
1022 
1023 
1024 template <int dim, int spacedim>
1025 const std::vector<Point<dim-1> > &
1027 {
1028  // a finite element may define
1029  // support points, but only if
1030  // there are as many as there are
1031  // degrees of freedom on a face
1032  Assert ((unit_face_support_points.size() == 0) ||
1033  (unit_face_support_points.size() == this->dofs_per_face),
1034  ExcInternalError());
1035  return unit_face_support_points;
1036 }
1037 
1038 
1039 
1040 template <int dim, int spacedim>
1041 bool
1043 {
1044  return (unit_face_support_points.size() != 0);
1045 }
1046 
1047 
1048 
1049 template <int dim, int spacedim>
1050 const std::vector<Point<dim-1> > &
1052 {
1053  // a finite element may define
1054  // support points, but only if
1055  // there are as many as there are
1056  // degrees of freedom on a face
1057  return ((generalized_face_support_points.size() == 0)
1058  ? unit_face_support_points
1059  : generalized_face_support_points);
1060 }
1061 
1062 
1063 
1064 template <int dim, int spacedim>
1065 bool
1067 {
1068  return (generalized_face_support_points.size() != 0);
1069 }
1070 
1071 
1072 
1073 template <int dim, int spacedim>
1074 Point<dim-1>
1076 {
1077  Assert (index < this->dofs_per_face,
1078  ExcIndexRange (index, 0, this->dofs_per_face));
1079  Assert (unit_face_support_points.size() == this->dofs_per_face,
1080  ExcFEHasNoSupportPoints ());
1081  return unit_face_support_points[index];
1082 }
1083 
1084 
1085 
1086 template <int dim, int spacedim>
1087 bool
1089  const unsigned int,
1090  const unsigned int) const
1091 {
1092  return true;
1093 }
1094 
1095 
1096 
1097 template <int dim, int spacedim>
1100 get_sub_fe (const ComponentMask &mask) const
1101 {
1102  // Translate the ComponentMask into first_selected and n_components after
1103  // some error checking:
1104  const unsigned int n_total_components = this->n_components();
1105  Assert((n_total_components == mask.size()) || (mask.size()==0),
1106  ExcMessage("The given ComponentMask has the wrong size."));
1107 
1108  const unsigned int n_selected = mask.n_selected_components(n_total_components);
1109  Assert(n_selected>0,
1110  ExcMessage("You need at least one selected component."));
1111 
1112  const unsigned int first_selected = mask.first_selected_component(n_total_components);
1113 
1114 #ifdef DEBUG
1115  // check that it is contiguous:
1116  for (unsigned int c=0; c<n_total_components; ++c)
1117  Assert((c<first_selected && (!mask[c]))
1118  ||
1119  (c>=first_selected && c<first_selected+n_selected && mask[c])
1120  ||
1121  (c>=first_selected+n_selected && !mask[c]),
1122  ExcMessage("Error: the given ComponentMask is not contiguous!"));
1123 #endif
1124 
1125  return get_sub_fe(first_selected, n_selected);
1126 }
1127 
1128 
1129 
1130 template <int dim, int spacedim>
1133 get_sub_fe (const unsigned int first_component,
1134  const unsigned int n_selected_components) const
1135 {
1136  // No complicated logic is needed here, because it is overridden in
1137  // FESystem<dim,spacedim>. Just make sure that what the user chose is valid:
1138  Assert(first_component == 0 && n_selected_components == this->n_components(),
1139  ExcMessage("You can only select a whole FiniteElement, not a part of one."));
1140 
1141  (void)first_component;
1142  (void)n_selected_components;
1143  return *this;
1144 }
1145 
1146 
1147 
1148 template <int dim, int spacedim>
1149 std::pair<Table<2,bool>, std::vector<unsigned int> >
1151 {
1152  Assert (false, ExcNotImplemented());
1153  return std::pair<Table<2,bool>, std::vector<unsigned int> >
1154  (Table<2,bool>(this->n_components(), this->dofs_per_cell),
1155  std::vector<unsigned int>(this->n_components()));
1156 }
1157 
1158 
1159 
1160 template <int dim, int spacedim>
1161 void
1164  std::vector<double> &) const
1165 {
1166  Assert (has_generalized_support_points(),
1167  ExcMessage ("The element for which you are calling the current "
1168  "function does not have generalized support points (see "
1169  "the glossary for a definition of generalized support "
1170  "points). Consequently, the current function can not "
1171  "be defined and is not implemented by the element."));
1172  Assert (false, ExcNotImplemented());
1173 }
1174 
1175 
1176 
1177 template <int dim, int spacedim>
1178 std::size_t
1180 {
1181  return (sizeof(FiniteElementData<dim>) +
1183  MemoryConsumption::memory_consumption (prolongation) +
1184  MemoryConsumption::memory_consumption (interface_constraints) +
1185  MemoryConsumption::memory_consumption (system_to_component_table) +
1186  MemoryConsumption::memory_consumption (face_system_to_component_table) +
1187  MemoryConsumption::memory_consumption (system_to_base_table) +
1188  MemoryConsumption::memory_consumption (face_system_to_base_table) +
1189  MemoryConsumption::memory_consumption (component_to_base_table) +
1190  MemoryConsumption::memory_consumption (restriction_is_additive_flags) +
1191  MemoryConsumption::memory_consumption (nonzero_components) +
1192  MemoryConsumption::memory_consumption (n_nonzero_components_table));
1193 }
1194 
1195 
1196 
1197 template <int dim, int spacedim>
1198 std::vector<unsigned int>
1200  const std::vector<ComponentMask> &nonzero_components)
1201 {
1202  std::vector<unsigned int> retval (nonzero_components.size());
1203  for (unsigned int i=0; i<nonzero_components.size(); ++i)
1204  retval[i] = nonzero_components[i].n_selected_components();
1205  return retval;
1206 }
1207 
1208 
1209 
1210 /*------------------------------- FiniteElement ----------------------*/
1211 
1212 template <int dim, int spacedim>
1213 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
1215  const Mapping<dim,spacedim> &mapping,
1216  const Quadrature<dim-1> &quadrature,
1218 {
1219  return get_data (flags, mapping,
1221  output_data);
1222 }
1223 
1224 
1225 
1226 template <int dim, int spacedim>
1227 std::unique_ptr<typename FiniteElement<dim,spacedim>::InternalDataBase>
1229  const Mapping<dim,spacedim> &mapping,
1230  const Quadrature<dim-1> &quadrature,
1232 {
1233  return get_data (flags, mapping,
1235  output_data);
1236 }
1237 
1238 
1239 
1240 template <int dim, int spacedim>
1242 FiniteElement<dim,spacedim>::base_element(const unsigned int index) const
1243 {
1244  (void)index;
1245  Assert (index==0, ExcIndexRange(index,0,1));
1246  // This function should not be
1247  // called for a system element
1248  Assert (base_to_block_indices.size() == 1, ExcInternalError());
1249  return *this;
1250 }
1251 
1252 
1253 
1254 /*------------------------------- Explicit Instantiations -------------*/
1255 #include "fe.inst"
1256 
1257 
1258 DEAL_II_NAMESPACE_CLOSE
bool represents_the_all_selected_mask() const
Definition: block_mask.h:335
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:1214
static const unsigned int invalid_unsigned_int
Definition: types.h:173
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:439
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:217
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
const unsigned int components
Definition: fe_base.h:305
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
virtual std::size_t memory_consumption() const
Definition: fe.cc:47
const std::vector< Point< dim-1 > > & get_unit_face_support_points() const
Definition: fe.cc:1026
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:863
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
Definition: fe.cc:939
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:229
FullMatrix< double > interface_constraints
Definition: fe.h:2401
bool operator!=(const FiniteElement< dim, spacedim > &) const
Definition: fe.cc:957
const std::vector< Point< dim-1 > > & get_generalized_face_support_points() const
Definition: fe.cc:1051
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:56
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
Definition: fe.cc:1100
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:171
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:847
virtual Point< dim > unit_support_point(const unsigned int index) const
Definition: fe.cc:1013
bool restriction_is_implemented() const
Definition: fe.cc:695
const std::vector< Point< dim > > & get_generalized_support_points() const
Definition: fe.cc:991
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:183
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:240
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
const unsigned int dofs_per_quad
Definition: fe_base.h:250
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
Definition: fe.h:2449
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:918
bool isotropic_prolongation_is_implemented() const
Definition: fe.cc:721
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition: fe.h:2534
bool has_generalized_support_points() const
Definition: fe.cc:1004
bool represents_the_all_selected_mask() const
STL namespace.
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:357
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:1199
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
Definition: fe.h:2505
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:322
bool is_primitive() const
Definition: fe.h:3228
const std::vector< unsigned int > n_nonzero_components_table
Definition: fe.h:2559
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
const unsigned int dofs_per_line
Definition: fe_base.h:244
bool has_face_support_points() const
Definition: fe.cc:1042
bool has_support_points() const
Definition: fe.cc:982
unsigned int component_to_block_index(const unsigned int component) const
Definition: fe.cc:344
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:794
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:194
virtual std::size_t memory_consumption() const
Definition: fe.cc:1179
No update.
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
#define Assert(cond, exc)
Definition: exceptions.h:1142
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:302
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe.cc:879
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
Definition: fe.cc:773
unsigned int size() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1242
virtual std::string get_name() const =0
bool has_generalized_face_support_points() const
Definition: fe.cc:1066
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:896
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:525
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
Definition: fe.cc:1075
const unsigned int dofs_per_cell
Definition: fe_base.h:295
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:646
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:966
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:907
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:160
unsigned int n_components() const
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
Definition: fe.h:2469
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:1228
virtual bool hp_constraints_are_implemented() const
Definition: fe.cc:785
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
const unsigned int dofs_per_face
Definition: fe_base.h:288
bool prolongation_is_implemented() const
Definition: fe.cc:669
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:263
const bool cached_primitivity
Definition: fe.h:2566
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:206
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2550
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
static ::ExceptionBase & ExcNotImplemented()
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
Definition: fe.h:2480
bool isotropic_restriction_is_implemented() const
Definition: fe.cc:747
BlockIndices base_to_block_indices
Definition: fe.h:2511
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:252
unsigned int size() const
Definition: block_mask.h:250
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
Definition: fe.cc:614
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition: fe.h:2464
TableIndices< 2 > interface_constraints_size() const
Definition: fe.cc:820
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe.cc:1150
void fill(InputIterator entries, const bool C_style_indexing=true)
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2541
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
Definition: fe.cc:151
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:1163
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:929
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe.cc:1088
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition: fe.h:2499
static ::ExceptionBase & ExcInternalError()