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