16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/quadrature.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/grid/tria.h> 20 #include <deal.II/grid/tria_iterator.h> 21 #include <deal.II/dofs/dof_accessor.h> 22 #include <deal.II/fe/mapping.h> 23 #include <deal.II/fe/fe_system.h> 24 #include <deal.II/fe/fe_values.h> 29 DEAL_II_NAMESPACE_OPEN
33 bool IsNonZero (
unsigned int i)
38 unsigned int count_nonzeros(
const std::vector<unsigned int> &vec)
40 return std::count_if(vec.begin(), vec.end(), IsNonZero);
47 template <
int dim,
int spacedim>
50 base_fe_datas(n_base_elements),
51 base_fe_output_objects(n_base_elements)
56 template <
int dim,
int spacedim>
70 template <
int dim,
int spacedim>
82 template <
int dim,
int spacedim>
95 template <
int dim,
int spacedim>
110 template <
int dim,
int spacedim>
121 template <
int dim,
int spacedim>
124 const std::vector<unsigned int> &multiplicities)
128 unsigned int multiplied_dofs_per_vertex = 0;
129 unsigned int multiplied_dofs_per_line = 0;
130 unsigned int multiplied_dofs_per_quad = 0;
131 unsigned int multiplied_dofs_per_hex = 0;
133 unsigned int multiplied_n_components = 0;
137 for (
unsigned int i=0; i<fes.size(); i++)
138 if (multiplicities[i]>0)
140 multiplied_dofs_per_vertex += fes[i]->dofs_per_vertex * multiplicities[i];
141 multiplied_dofs_per_line += fes[i]->dofs_per_line * multiplicities[i];
142 multiplied_dofs_per_quad += fes[i]->dofs_per_quad * multiplicities[i];
143 multiplied_dofs_per_hex += fes[i]->dofs_per_hex * multiplicities[i];
145 multiplied_n_components+=fes[i]->n_components() * multiplicities[i];
156 unsigned int index = 0;
157 for (index=0; index<fes.size(); ++index)
158 if (multiplicities[index]>0)
164 for (; index<fes.size(); ++index)
165 if (multiplicities[index]>0)
172 std::vector<unsigned int> dpo;
173 dpo.push_back(multiplied_dofs_per_vertex);
174 dpo.push_back(multiplied_dofs_per_line);
175 if (dim>1) dpo.push_back(multiplied_dofs_per_quad);
176 if (dim>2) dpo.push_back(multiplied_dofs_per_hex);
180 for (
unsigned int base=0; base < fes.size(); ++base)
181 for (
unsigned int m = 0; m < multiplicities[base]; ++m)
182 block_indices.
push_back(fes[base]->dofs_per_cell);
185 multiplied_n_components,
194 template <
int dim,
int spacedim>
197 const unsigned int N1,
199 const unsigned int N2=0,
201 const unsigned int N3=0,
203 const unsigned int N4=0,
205 const unsigned int N5=0)
207 std::vector<const FiniteElement<dim,spacedim>*> fes;
214 std::vector<unsigned int> mult;
220 return multiply_dof_numbers(fes, mult);
230 template <
int dim,
int spacedim>
233 const std::vector<unsigned int> &multiplicities)
239 unsigned int n_shape_functions = 0;
240 for (
unsigned int i=0; i<fes.size(); ++i)
241 if (multiplicities[i]>0)
242 n_shape_functions += fes[i]->dofs_per_cell * multiplicities[i];
245 std::vector<bool> retval (n_shape_functions,
false);
254 unsigned int total_index = 0;
255 for (
unsigned int vertex_number=0;
256 vertex_number<GeometryInfo<dim>::vertices_per_cell;
259 for (
unsigned int base=0; base<fes.size(); ++base)
260 for (
unsigned int m=0; m<multiplicities[base]; ++m)
261 for (
unsigned int local_index = 0;
262 local_index < fes[base]->dofs_per_vertex;
263 ++local_index, ++total_index)
265 const unsigned int index_in_base
266 = (fes[base]->dofs_per_vertex*vertex_number +
271 retval[total_index] = fes[base]->restriction_is_additive(index_in_base);
277 for (
unsigned int line_number= 0;
281 for (
unsigned int base=0; base<fes.size(); ++base)
282 for (
unsigned int m=0; m<multiplicities[base]; ++m)
283 for (
unsigned int local_index = 0;
284 local_index < fes[base]->dofs_per_line;
285 ++local_index, ++total_index)
287 const unsigned int index_in_base
288 = (fes[base]->dofs_per_line*line_number +
290 fes[base]->first_line_index);
294 retval[total_index] = fes[base]->restriction_is_additive(index_in_base);
300 for (
unsigned int quad_number= 0;
304 for (
unsigned int base=0; base<fes.size(); ++base)
305 for (
unsigned int m=0; m<multiplicities[base]; ++m)
306 for (
unsigned int local_index = 0;
307 local_index < fes[base]->dofs_per_quad;
308 ++local_index, ++total_index)
310 const unsigned int index_in_base
311 = (fes[base]->dofs_per_quad*quad_number +
313 fes[base]->first_quad_index);
317 retval[total_index] = fes[base]->restriction_is_additive(index_in_base);
323 for (
unsigned int hex_number= 0;
327 for (
unsigned int base=0; base<fes.size(); ++base)
328 for (
unsigned int m=0; m<multiplicities[base]; ++m)
329 for (
unsigned int local_index = 0;
330 local_index < fes[base]->dofs_per_hex;
331 ++local_index, ++total_index)
333 const unsigned int index_in_base
334 = (fes[base]->dofs_per_hex*hex_number +
336 fes[base]->first_hex_index);
340 retval[total_index] = fes[base]->restriction_is_additive(index_in_base);
344 Assert (total_index == n_shape_functions, ExcInternalError());
357 template <
int dim,
int spacedim>
360 const unsigned int N1,
362 const unsigned int N2=0,
364 const unsigned int N3=0,
366 const unsigned int N4=0,
368 const unsigned int N5=0)
370 std::vector<const FiniteElement<dim,spacedim>*> fe_list;
371 std::vector<unsigned int> multiplicities;
373 fe_list.push_back (fe1);
374 multiplicities.push_back (N1);
376 fe_list.push_back (fe2);
377 multiplicities.push_back (N2);
379 fe_list.push_back (fe3);
380 multiplicities.push_back (N3);
382 fe_list.push_back (fe4);
383 multiplicities.push_back (N4);
385 fe_list.push_back (fe5);
386 multiplicities.push_back (N5);
387 return compute_restriction_is_additive_flags (fe_list, multiplicities);
396 template <
int dim,
int spacedim>
397 std::vector<ComponentMask>
399 const std::vector<unsigned int> &multiplicities)
405 unsigned int n_shape_functions = 0;
406 for (
unsigned int i=0; i<fes.size(); ++i)
407 if (multiplicities[i]>0)
408 n_shape_functions += fes[i]->dofs_per_cell * multiplicities[i];
411 for (
unsigned int i=0; i<fes.size(); ++i)
412 if (multiplicities[i]>0)
413 n_components += fes[i]->n_components() * multiplicities[i];
416 std::vector<std::vector<bool> >
417 retval (n_shape_functions, std::vector<bool> (n_components,
false));
427 unsigned int total_index = 0;
428 for (
unsigned int vertex_number=0;
429 vertex_number<GeometryInfo<dim>::vertices_per_cell;
432 unsigned int comp_start = 0;
433 for (
unsigned int base=0; base<fes.size(); ++base)
434 for (
unsigned int m=0; m<multiplicities[base];
435 ++m, comp_start+=fes[base]->n_components())
436 for (
unsigned int local_index = 0;
437 local_index < fes[base]->dofs_per_vertex;
438 ++local_index, ++total_index)
440 const unsigned int index_in_base
441 = (fes[base]->dofs_per_vertex*vertex_number +
445 retval[total_index].size(),
447 for (
unsigned int c=0; c<fes[base]->n_components(); ++c)
451 retval[total_index][comp_start+c]
452 = fes[base]->get_nonzero_components(index_in_base)[c];
459 for (
unsigned int line_number= 0;
463 unsigned int comp_start = 0;
464 for (
unsigned int base=0; base<fes.size(); ++base)
465 for (
unsigned int m=0; m<multiplicities[base];
466 ++m, comp_start+=fes[base]->n_components())
467 for (
unsigned int local_index = 0;
468 local_index < fes[base]->dofs_per_line;
469 ++local_index, ++total_index)
471 const unsigned int index_in_base
472 = (fes[base]->dofs_per_line*line_number +
474 fes[base]->first_line_index);
477 retval[total_index].size(),
479 for (
unsigned int c=0; c<fes[base]->n_components(); ++c)
483 retval[total_index][comp_start+c]
484 = fes[base]->get_nonzero_components(index_in_base)[c];
491 for (
unsigned int quad_number= 0;
495 unsigned int comp_start = 0;
496 for (
unsigned int base=0; base<fes.size(); ++base)
497 for (
unsigned int m=0; m<multiplicities[base];
498 ++m, comp_start+=fes[base]->n_components())
499 for (
unsigned int local_index = 0;
500 local_index < fes[base]->dofs_per_quad;
501 ++local_index, ++total_index)
503 const unsigned int index_in_base
504 = (fes[base]->dofs_per_quad*quad_number +
506 fes[base]->first_quad_index);
509 retval[total_index].size(),
511 for (
unsigned int c=0; c<fes[base]->n_components(); ++c)
515 retval[total_index][comp_start+c]
516 = fes[base]->get_nonzero_components(index_in_base)[c];
523 for (
unsigned int hex_number= 0;
527 unsigned int comp_start = 0;
528 for (
unsigned int base=0; base<fes.size(); ++base)
529 for (
unsigned int m=0; m<multiplicities[base];
530 ++m, comp_start+=fes[base]->n_components())
531 for (
unsigned int local_index = 0;
532 local_index < fes[base]->dofs_per_hex;
533 ++local_index, ++total_index)
535 const unsigned int index_in_base
536 = (fes[base]->dofs_per_hex*hex_number +
538 fes[base]->first_hex_index);
541 retval[total_index].size(),
543 for (
unsigned int c=0; c<fes[base]->n_components(); ++c)
547 retval[total_index][comp_start+c]
548 = fes[base]->get_nonzero_components(index_in_base)[c];
553 Assert (total_index == n_shape_functions, ExcInternalError());
559 std::vector<ComponentMask> xretval (retval.size());
560 for (
unsigned int i=0; i<retval.size(); ++i)
569 template <
int dim,
int spacedim>
570 std::vector<ComponentMask>
572 const unsigned int N1,
574 const unsigned int N2=0,
576 const unsigned int N3=0,
578 const unsigned int N4=0,
580 const unsigned int N5=0)
582 std::vector<const FiniteElement<dim,spacedim>*> fe_list;
583 std::vector<unsigned int> multiplicities;
585 fe_list.push_back (fe1);
586 multiplicities.push_back (N1);
588 fe_list.push_back (fe2);
589 multiplicities.push_back (N2);
591 fe_list.push_back (fe3);
592 multiplicities.push_back (N3);
594 fe_list.push_back (fe4);
595 multiplicities.push_back (N4);
597 fe_list.push_back (fe5);
598 multiplicities.push_back (N5);
600 return compute_nonzero_components (fe_list, multiplicities);
606 template <
int dim,
int spacedim>
608 const unsigned int n_elements) :
609 FiniteElement<dim,spacedim> (multiply_dof_numbers(&fe, n_elements),
610 compute_restriction_is_additive_flags (&fe, n_elements),
611 compute_nonzero_components(&fe, n_elements)),
614 std::vector<const FiniteElement<dim,spacedim>*> fes;
616 std::vector<unsigned int> multiplicities;
617 multiplicities.push_back(n_elements);
623 template <
int dim,
int spacedim>
625 const unsigned int n1,
627 const unsigned int n2) :
628 FiniteElement<dim,spacedim> (multiply_dof_numbers(&fe1, n1, &fe2, n2),
629 compute_restriction_is_additive_flags (&fe1, n1,
631 compute_nonzero_components(&fe1, n1,
635 std::vector<const FiniteElement<dim,spacedim>*> fes;
638 std::vector<unsigned int> multiplicities;
639 multiplicities.push_back(n1);
640 multiplicities.push_back(n2);
646 template <
int dim,
int spacedim>
648 const unsigned int n1,
650 const unsigned int n2,
652 const unsigned int n3) :
656 compute_restriction_is_additive_flags (&fe1, n1,
659 compute_nonzero_components(&fe1, n1,
664 std::vector<const FiniteElement<dim,spacedim>*> fes;
668 std::vector<unsigned int> multiplicities;
669 multiplicities.push_back(n1);
670 multiplicities.push_back(n2);
671 multiplicities.push_back(n3);
677 template <
int dim,
int spacedim>
679 const unsigned int n1,
681 const unsigned int n2,
683 const unsigned int n3,
685 const unsigned int n4) :
690 compute_restriction_is_additive_flags (&fe1, n1,
694 compute_nonzero_components(&fe1, n1,
700 std::vector<const FiniteElement<dim,spacedim>*> fes;
705 std::vector<unsigned int> multiplicities;
706 multiplicities.push_back(n1);
707 multiplicities.push_back(n2);
708 multiplicities.push_back(n3);
709 multiplicities.push_back(n4);
715 template <
int dim,
int spacedim>
717 const unsigned int n1,
719 const unsigned int n2,
721 const unsigned int n3,
723 const unsigned int n4,
725 const unsigned int n5) :
731 compute_restriction_is_additive_flags (&fe1, n1,
736 compute_nonzero_components(&fe1, n1,
743 std::vector<const FiniteElement<dim,spacedim>*> fes;
749 std::vector<unsigned int> multiplicities;
750 multiplicities.push_back(n1);
751 multiplicities.push_back(n2);
752 multiplicities.push_back(n3);
753 multiplicities.push_back(n4);
754 multiplicities.push_back(n5);
760 template <
int dim,
int spacedim>
763 const std::vector<unsigned int> &multiplicities)
765 FiniteElement<dim,spacedim> (multiply_dof_numbers(fes, multiplicities),
766 compute_restriction_is_additive_flags (fes, multiplicities),
767 compute_nonzero_components(fes, multiplicities)),
775 template <
int dim,
int spacedim>
781 template <
int dim,
int spacedim>
792 std::ostringstream namebuf;
794 namebuf <<
"FESystem<" 807 return namebuf.str();
812 template <
int dim,
int spacedim>
816 std::vector< const FiniteElement<dim,spacedim>* > fes;
817 std::vector<unsigned int> multiplicities;
829 template <
int dim,
int spacedim>
839 .shape_value(this->system_to_base_table[i].second, p));
844 template <
int dim,
int spacedim>
848 const unsigned int component)
const 881 template <
int dim,
int spacedim>
891 .shape_grad(this->system_to_base_table[i].second, p));
896 template <
int dim,
int spacedim>
900 const unsigned int component)
const 926 template <
int dim,
int spacedim>
936 .shape_grad_grad(this->system_to_base_table[i].second, p));
941 template <
int dim,
int spacedim>
945 const unsigned int component)
const 971 template <
int dim,
int spacedim>
981 .shape_3rd_derivative(this->system_to_base_table[i].second, p));
986 template <
int dim,
int spacedim>
990 const unsigned int component)
const 1011 component_in_base));
1016 template <
int dim,
int spacedim>
1026 .shape_4th_derivative(this->system_to_base_table[i].second, p));
1031 template <
int dim,
int spacedim>
1035 const unsigned int component)
const 1056 component_in_base));
1061 template <
int dim,
int spacedim>
1074 ExcDimensionMismatch (interpolation_matrix.
m(),
1079 ExcDimensionMismatch (interpolation_matrix.
m(),
1092 ExcInterpolationNotImplemented());
1101 ExcInterpolationNotImplemented());
1106 source_fe.element_multiplicity(i),
1108 ExcInterpolationNotImplemented());
1117 std::vector<FullMatrix<double> > base_matrices (this->
n_base_elements());
1121 source_fe.base_element(i).dofs_per_cell);
1122 base_element(i).get_interpolation_matrix (source_fe.base_element(i),
1129 interpolation_matrix = 0;
1131 for (
unsigned int j=0; j<source_fe.dofs_per_cell; ++j)
1133 source_fe.system_to_base_table[j].first)
1134 interpolation_matrix(i,j)
1137 source_fe.system_to_base_table[j].second));
1142 template <
int dim,
int spacedim>
1151 ExcMessage(
"Restriction matrices are only available for refined cells!"));
1156 if (this->
restriction[refinement_case-1][child].n() == 0)
1161 if (this->
restriction[refinement_case-1][child].n() ==
1163 return this->
restriction[refinement_case-1][child];
1166 bool do_restriction =
true;
1169 std::vector<const FullMatrix<double> *>
1175 &
base_element(i).get_restriction_matrix(child, refinement_case);
1177 do_restriction =
false;
1218 restriction(i,j) = (*base_matrices[base])(base_index_i,base_index_j);
1222 (this->restriction[refinement_case-1][child]));
1226 return this->
restriction[refinement_case-1][child];
1231 template <
int dim,
int spacedim>
1240 ExcMessage(
"Restriction matrices are only available for refined cells!"));
1246 if (this->
prolongation[refinement_case-1][child].n() == 0)
1250 if (this->
prolongation[refinement_case-1][child].n() ==
1254 bool do_prolongation =
true;
1255 std::vector<const FullMatrix<double> *>
1260 &
base_element(i).get_prolongation_matrix(child, refinement_case);
1262 do_prolongation =
false;
1267 if (do_prolongation)
1284 prolongate(i,j) = (*base_matrices[base])(base_index_i,base_index_j);
1295 template <
int dim,
int spacedim>
1299 const unsigned int face,
1300 const bool face_orientation,
1301 const bool face_flip,
1302 const bool face_rotation)
const 1307 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1311 base_face_to_cell_index
1312 = this->
base_element(face_base_index.first.first).face_to_cell_index (face_base_index.second,
1323 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1324 target = std::make_pair (face_base_index.first, base_face_to_cell_index);
1329 Assert (
false, ExcInternalError());
1341 template <
int dim,
int spacedim>
1348 for (
unsigned int base_no=0; base_no<this->
n_base_elements(); ++base_no)
1349 out |=
base_element(base_no).requires_update_flags (flags);
1355 template <
int dim,
int spacedim>
1385 for (
unsigned int base_no=0; base_no<this->
n_base_elements(); ++base_no)
1390 flags |
base_element(base_no).requires_update_flags(flags));
1399 base_element(base_no).get_data (flags, mapping, quadrature,
1400 base_fe_output_object);
1411 template <
int dim,
int spacedim>
1441 for (
unsigned int base_no=0; base_no<this->
n_base_elements(); ++base_no)
1446 flags |
base_element(base_no).requires_update_flags(flags));
1455 base_element(base_no).get_face_data (flags, mapping, quadrature,
1456 base_fe_output_object);
1469 template <
int dim,
int spacedim>
1499 for (
unsigned int base_no=0; base_no<this->
n_base_elements(); ++base_no)
1504 flags |
base_element(base_no).requires_update_flags(flags));
1513 base_element(base_no).get_subface_data (flags, mapping, quadrature,
1514 base_fe_output_object);
1524 template <
int dim,
int spacedim>
1528 const CellSimilarity::Similarity cell_similarity,
1532 const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1537 quadrature, cell_similarity, mapping_internal, fe_internal,
1538 mapping_data, output_data);
1543 template <
int dim,
int spacedim>
1547 const unsigned int face_no,
1551 const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1556 CellSimilarity::none, mapping_internal, fe_internal,
1557 mapping_data, output_data);
1563 template <
int dim,
int spacedim>
1567 const unsigned int face_no,
1568 const unsigned int sub_no,
1572 const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1576 compute_fill (mapping, cell, face_no, sub_no, quadrature,
1577 CellSimilarity::none, mapping_internal, fe_internal,
1578 mapping_data, output_data);
1583 template <
int dim,
int spacedim>
1584 template <
int dim_1>
1589 const unsigned int face_no,
1590 const unsigned int sub_no,
1592 const CellSimilarity::Similarity cell_similarity,
1602 Assert (dynamic_cast<const InternalData *> (&fe_internal) != 0, ExcInternalError());
1608 Assert(dim_1==dim || dim_1==dim-1, ExcInternalError());
1624 for (
unsigned int base_no=0; base_no<this->
n_base_elements(); ++base_no)
1636 const Quadrature<dim-1> *face_quadrature = 0;
1637 const unsigned int n_q_points = quadrature.
size();
1641 const Subscriptor *quadrature_base_pointer = &quadrature;
1645 Assert(dim_1==dim, ExcDimensionMismatch(dim_1,dim));
1647 ExcInternalError());
1654 Assert(dim_1==dim-1, ExcDimensionMismatch(dim_1,dim-1));
1655 Assert (
dynamic_cast<const Quadrature<dim-1
> *>(quadrature_base_pointer) != 0,
1656 ExcInternalError());
1659 =
static_cast<const Quadrature<dim-1
> *>(quadrature_base_pointer);
1673 mapping, mapping_internal, mapping_data,
1674 base_fe_data, base_data);
1678 mapping, mapping_internal, mapping_data,
1679 base_fe_data, base_data);
1683 mapping, mapping_internal, mapping_data,
1684 base_fe_data, base_data);
1703 if (cell_similarity != CellSimilarity::translation)
1704 for (
unsigned int system_index=0; system_index<this->
dofs_per_cell;
1718 unsigned int out_index = 0;
1719 for (
unsigned int i=0; i<system_index; ++i)
1721 unsigned int in_index = 0;
1722 for (
unsigned int i=0; i<base_index; ++i)
1728 ExcInternalError());
1732 for (
unsigned int q=0; q<n_q_points; ++q)
1738 for (
unsigned int q=0; q<n_q_points; ++q)
1744 for (
unsigned int q=0; q<n_q_points; ++q)
1750 for (
unsigned int q=0; q<n_q_points; ++q)
1760 template <
int dim,
int spacedim>
1769 unsigned int total_index = 0;
1774 for (
unsigned int k=0; k<
base_element(base).n_components(); ++k)
1776 = std::make_pair(std::make_pair(base,k), m);
1779 ExcInternalError());
1784 const std::pair<unsigned int, unsigned int>
1792 for (
unsigned int vertex_number=0;
1793 vertex_number<GeometryInfo<dim>::vertices_per_cell;
1796 unsigned int comp_start = 0;
1800 for (
unsigned int local_index = 0;
1802 ++local_index, ++total_index)
1804 const unsigned int index_in_base
1809 = std::make_pair (std::make_pair(base, m), index_in_base);
1813 const unsigned int comp_in_base
1814 =
base_element(base).system_to_component_index(index_in_base).first;
1815 const unsigned int comp
1816 = comp_start + comp_in_base;
1817 const unsigned int index_in_comp
1818 =
base_element(base).system_to_component_index(index_in_base).second;
1820 = std::make_pair (comp, index_in_comp);
1829 for (
unsigned int line_number= 0;
1833 unsigned int comp_start = 0;
1837 for (
unsigned int local_index = 0;
1839 ++local_index, ++total_index)
1841 const unsigned int index_in_base
1847 = std::make_pair (std::make_pair(base,m), index_in_base);
1851 const unsigned int comp_in_base
1852 =
base_element(base).system_to_component_index(index_in_base).first;
1853 const unsigned int comp
1854 = comp_start + comp_in_base;
1855 const unsigned int index_in_comp
1856 =
base_element(base).system_to_component_index(index_in_base).second;
1858 = std::make_pair (comp, index_in_comp);
1867 for (
unsigned int quad_number= 0;
1871 unsigned int comp_start = 0;
1875 for (
unsigned int local_index = 0;
1877 ++local_index, ++total_index)
1879 const unsigned int index_in_base
1885 = std::make_pair (std::make_pair(base,m), index_in_base);
1889 const unsigned int comp_in_base
1890 =
base_element(base).system_to_component_index(index_in_base).first;
1891 const unsigned int comp
1892 = comp_start + comp_in_base;
1893 const unsigned int index_in_comp
1894 =
base_element(base).system_to_component_index(index_in_base).second;
1896 = std::make_pair (comp, index_in_comp);
1905 for (
unsigned int hex_number= 0;
1909 unsigned int comp_start = 0;
1913 for (
unsigned int local_index = 0;
1915 ++local_index, ++total_index)
1917 const unsigned int index_in_base
1923 = std::make_pair (std::make_pair(base,m), index_in_base);
1927 const unsigned int comp_in_base
1928 =
base_element(base).system_to_component_index(index_in_base).first;
1929 const unsigned int comp
1930 = comp_start + comp_in_base;
1931 const unsigned int index_in_comp
1932 =
base_element(base).system_to_component_index(index_in_base).second;
1934 = std::make_pair (comp, index_in_comp);
1944 template <
int dim,
int spacedim>
1952 const std::pair<unsigned int, unsigned int>
1957 unsigned int total_index = 0;
1958 for (
unsigned int vertex_number=0;
1959 vertex_number<GeometryInfo<dim>::vertices_per_face;
1962 unsigned int comp_start = 0;
1966 for (
unsigned int local_index = 0;
1968 ++local_index, ++total_index)
1978 const unsigned int index_in_base
1982 const unsigned int face_index_in_base
1987 = std::make_pair (std::make_pair(base,m), face_index_in_base);
1991 const unsigned int comp_in_base
1992 =
base_element(base).face_system_to_component_index(face_index_in_base).first;
1993 const unsigned int comp
1994 = comp_start + comp_in_base;
1995 const unsigned int face_index_in_comp
1996 =
base_element(base).face_system_to_component_index(face_index_in_base).second;
1998 = std::make_pair (comp, face_index_in_comp);
2007 for (
unsigned int line_number= 0;
2011 unsigned int comp_start = 0;
2015 for (
unsigned int local_index = 0;
2017 ++local_index, ++total_index)
2020 const unsigned int index_in_base
2025 const unsigned int face_index_in_base
2031 = std::make_pair (std::make_pair(base,m), face_index_in_base);
2035 const unsigned int comp_in_base
2036 =
base_element(base).face_system_to_component_index(face_index_in_base).first;
2037 const unsigned int comp
2038 = comp_start + comp_in_base;
2039 const unsigned int face_index_in_comp
2040 =
base_element(base).face_system_to_component_index(face_index_in_base).second;
2042 = std::make_pair (comp, face_index_in_comp);
2051 for (
unsigned int quad_number= 0;
2055 unsigned int comp_start = 0;
2059 for (
unsigned int local_index = 0;
2061 ++local_index, ++total_index)
2064 const unsigned int index_in_base
2069 const unsigned int face_index_in_base
2075 = std::make_pair (std::make_pair(base,m), face_index_in_base);
2079 const unsigned int comp_in_base
2080 =
base_element(base).face_system_to_component_index(face_index_in_base).first;
2081 const unsigned int comp
2082 = comp_start + comp_in_base;
2083 const unsigned int face_index_in_comp
2084 =
base_element(base).face_system_to_component_index(face_index_in_base).second;
2086 = std::make_pair (comp, face_index_in_comp);
2094 ExcInternalError());
2096 ExcInternalError());
2101 template <
int dim,
int spacedim>
2109 if (
base_element(base).constraints_are_implemented() ==
false)
2133 const std::pair<std::pair<unsigned int,unsigned int>,
unsigned int> n_index
2138 std::pair<std::pair<unsigned int,unsigned int>,
unsigned int> m_index;
2145 Assert (
false, ExcInternalError());
2160 const unsigned int index_in_line
2162 const unsigned int sub_line
2164 Assert (sub_line < 2, ExcInternalError());
2170 const unsigned int tmp1 = 2*this->dofs_per_vertex+index_in_line;
2183 ExcInternalError());
2187 ExcInternalError());
2188 m_index.second =
base_element(m_index.first.first).dofs_per_vertex +
2189 base_element(m_index.first.first).dofs_per_line*sub_line +
2209 const unsigned int index_in_line
2211 const unsigned int sub_line
2213 Assert (sub_line < 12, ExcInternalError());
2215 const unsigned int tmp1 = 4*this->dofs_per_vertex+index_in_line;
2220 ExcInternalError());
2224 ExcInternalError());
2225 m_index.second = 5*
base_element(m_index.first.first).dofs_per_vertex +
2226 base_element(m_index.first.first).dofs_per_line*sub_line +
2233 const unsigned int index_in_quad
2237 ExcInternalError());
2238 const unsigned int sub_quad
2240 this->dofs_per_quad);
2241 Assert (sub_quad < 4, ExcInternalError());
2243 const unsigned int tmp1 = 4*this->dofs_per_vertex +
2247 ExcInternalError());
2253 ExcInternalError());
2258 ExcInternalError());
2259 m_index.second = 5*
base_element(m_index.first.first).dofs_per_vertex +
2261 base_element(m_index.first.first).dofs_per_quad*sub_quad +
2269 Assert (
false, ExcNotImplemented());
2275 if (n_index.first == m_index.first)
2277 = (
base_element(n_index.first.first).constraints()(m_index.second,
2284 template <
int dim,
int spacedim>
2286 const std::vector<unsigned int> &multiplicities)
2288 Assert (fes.size() == multiplicities.size(),
2289 ExcDimensionMismatch (fes.size(), multiplicities.size()) );
2291 ExcMessage (
"Need to pass at least one finite element."));
2292 Assert (count_nonzeros(multiplicities) > 0,
2293 ExcMessage(
"You only passed FiniteElements with multiplicity 0."));
2299 for (
unsigned int i=0; i<fes.size(); i++)
2300 if (multiplicities[i]>0)
2303 std::vector<Threads::Task<FiniteElement<dim,spacedim>*> > clone_base_elements;
2305 for (
unsigned int i=0; i<fes.size(); i++)
2306 if (multiplicities[i]>0)
2311 for (
unsigned int i=0; i<fes.size(); i++)
2313 if (multiplicities[i]>0)
2317 (clone_base_elements[ind].return_value()),
2323 Assert(ind>0, ExcInternalError());
2346 template <
int dim,
int spacedim>
2354 for (
unsigned int base_el=0; base_el<this->
n_base_elements(); ++base_el)
2371 ExcInternalError());
2378 template <
int dim,
int spacedim>
2396 for (
unsigned int base_el=0; base_el<this->
n_base_elements(); ++base_el)
2416 ExcInternalError());
2419 =
base_element(base_i).unit_face_support_points[index_in_base];
2425 template <
int dim,
int spacedim>
2440 unsigned int index = 0;
2444 = this->
base_element(b).adjust_quad_dof_index_for_face_orientation_table;
2447 for (
unsigned int i=0; i<temp.
size(0); ++i)
2448 for (
unsigned int j=0; j<8; ++j)
2451 index += temp.
size(0);
2455 ExcInternalError());
2463 const std::vector<int> &temp2
2464 = this->
base_element(b).adjust_line_dof_index_for_line_orientation_table;
2467 std::copy(temp2.begin(), temp2.end(),
2470 index += temp2.size();
2474 ExcInternalError());
2479 template <
int dim,
int spacedim>
2485 if (
base_element(b).hp_constraints_are_implemented() ==
false)
2493 template <
int dim,
int spacedim>
2504 ExcInterpolationNotImplemented());
2507 ExcDimensionMismatch (interpolation_matrix.
n(),
2510 ExcDimensionMismatch (interpolation_matrix.
m(),
2527 interpolation_matrix = 0;
2531 unsigned int base_index = 0,
2532 base_index_other = 0;
2533 unsigned int multiplicity = 0,
2534 multiplicity_other = 0;
2542 &base_other = fe_other_system->
base_element(base_index_other);
2545 ExcNotImplemented());
2548 base_to_base_interpolation.
reinit (base_other.dofs_per_face,
2551 base_to_base_interpolation);
2559 std::make_pair (base_index, multiplicity))
2560 for (
unsigned int j=0; j<fe_other_system->
dofs_per_face; ++j)
2563 std::make_pair (base_index_other, multiplicity_other))
2564 interpolation_matrix(j, i)
2577 ++multiplicity_other;
2578 if (multiplicity_other ==
2581 multiplicity_other = 0;
2590 ExcInternalError());
2597 ExcInternalError());
2603 template <
int dim,
int spacedim>
2607 const unsigned int subface,
2615 ExcInterpolationNotImplemented());
2618 ExcDimensionMismatch (interpolation_matrix.
n(),
2621 ExcDimensionMismatch (interpolation_matrix.
m(),
2638 interpolation_matrix = 0;
2642 unsigned int base_index = 0,
2643 base_index_other = 0;
2644 unsigned int multiplicity = 0,
2645 multiplicity_other = 0;
2653 &base_other = fe_other_system->
base_element(base_index_other);
2656 ExcNotImplemented());
2659 base_to_base_interpolation.
reinit (base_other.dofs_per_face,
2663 base_to_base_interpolation);
2671 std::make_pair (base_index, multiplicity))
2672 for (
unsigned int j=0; j<fe_other_system->
dofs_per_face; ++j)
2675 std::make_pair (base_index_other, multiplicity_other))
2676 interpolation_matrix(j, i)
2689 ++multiplicity_other;
2690 if (multiplicity_other ==
2693 multiplicity_other = 0;
2702 ExcInternalError());
2709 ExcInternalError());
2715 template <
int dim,
int spacedim>
2716 template <
int structdim>
2717 std::vector<std::pair<unsigned int, unsigned int> >
2737 unsigned int base_index = 0,
2738 base_index_other = 0;
2739 unsigned int multiplicity = 0,
2740 multiplicity_other = 0;
2744 unsigned int dof_offset = 0,
2745 dof_offset_other = 0;
2747 std::vector<std::pair<unsigned int, unsigned int> > identities;
2753 &base_other = fe_other_system->base_element(base_index_other);
2756 ExcNotImplemented());
2760 std::vector<std::pair<unsigned int, unsigned int> > base_identities;
2773 Assert (
false, ExcNotImplemented());
2776 for (
unsigned int i=0; i<base_identities.size(); ++i)
2777 identities.push_back
2778 (std::make_pair (base_identities[i].first + dof_offset,
2779 base_identities[i].second + dof_offset_other));
2782 dof_offset += base.template n_dofs_per_object<structdim>();
2783 dof_offset_other += base_other.template n_dofs_per_object<structdim>();
2794 ++multiplicity_other;
2795 if (multiplicity_other ==
2796 fe_other_system->element_multiplicity(base_index_other))
2798 multiplicity_other = 0;
2806 Assert (base_index_other == fe_other_system->n_base_elements(),
2807 ExcInternalError());
2813 Assert (base_index_other != fe_other_system->n_base_elements(),
2814 ExcInternalError());
2821 Assert (
false, ExcNotImplemented());
2822 return std::vector<std::pair<unsigned int, unsigned int> >();
2828 template <
int dim,
int spacedim>
2829 std::vector<std::pair<unsigned int, unsigned int> >
2832 return hp_object_dof_identities<0> (fe_other);
2835 template <
int dim,
int spacedim>
2836 std::vector<std::pair<unsigned int, unsigned int> >
2839 return hp_object_dof_identities<1> (fe_other);
2844 template <
int dim,
int spacedim>
2845 std::vector<std::pair<unsigned int, unsigned int> >
2848 return hp_object_dof_identities<2> (fe_other);
2853 template <
int dim,
int spacedim>
2864 ExcNotImplemented());
2866 ExcNotImplemented());
2869 domination = FiniteElementDomination::no_requirements;
2875 fe_sys_other->base_element(b).n_components(),
2876 ExcNotImplemented());
2878 fe_sys_other->element_multiplicity(b),
2879 ExcNotImplemented());
2885 .compare_for_face_domination (fe_sys_other->base_element(b)));
2886 domination = domination & base_domination;
2892 Assert (
false, ExcNotImplemented());
2893 return FiniteElementDomination::neither_element_dominates;
2898 template <
int dim,
int spacedim>
2909 template <
int dim,
int spacedim>
2912 const unsigned int face_index)
const 2915 .has_support_on_face(this->system_to_base_index(shape_index).second,
2921 template <
int dim,
int spacedim>
2926 ExcIndexRange (index, 0, this->dofs_per_cell));
2930 typename FEL::ExcFEHasNoSupportPoints ());
2944 template <
int dim,
int spacedim>
2949 ExcIndexRange (index, 0, this->dofs_per_face));
2953 typename FEL::ExcFEHasNoSupportPoints ());
2967 template <
int dim,
int spacedim>
2968 std::pair<Table<2,bool>, std::vector<unsigned int> >
2978 std::pair<Table<2,bool>, std::vector<unsigned int> >
2986 const unsigned int comp = components.size();
2989 Table<2,bool> new_constant_modes(comp+base_table.first.n_rows()*
2991 constant_modes.n_cols());
2992 for (
unsigned int r=0; r<comp; ++r)
2994 new_constant_modes(r,c) = constant_modes(r,c);
2995 constant_modes.
swap(new_constant_modes);
3002 std::pair<std::pair<unsigned int,unsigned int>,
unsigned int> ind
3004 if (ind.first.first == i)
3005 for (
unsigned int c=0; c<base_table.first.n_rows(); ++c)
3006 constant_modes(comp+ind.first.second*base_table.first.n_rows()+c,k)
3007 = base_table.first(c,ind.second);
3010 for (
unsigned int c=0; c<base_table.second.size(); ++c)
3011 components.push_back(comp+r*this->base_elements[i].first->n_components()
3012 +base_table.second[c]);
3015 return std::pair<Table<2,bool>, std::vector<unsigned int> >(constant_modes,
3021 template <
int dim,
int spacedim>
3039 #include "fe_system.inst" 3041 DEAL_II_NAMESPACE_CLOSE
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
void initialize_unit_support_points()
bool is_primitive() const
static const unsigned int invalid_unsigned_int
std::vector< std::vector< FullMatrix< double > > > restriction
void build_interface_constraints()
const unsigned int components
#define AssertDimension(dim1, dim2)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
void initialize(const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities)
void swap(TableBase< N, T > &v)
std::vector< std::pair< unsigned int, unsigned int > > hp_object_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
FullMatrix< double > interface_constraints
std::vector< std::pair< std_cxx11::shared_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
unsigned int n_nonzero_components(const unsigned int i) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual bool hp_constraints_are_implemented() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
::ExceptionBase & ExcMessage(std::string arg1)
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
const unsigned int dofs_per_quad
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
const unsigned int degree
virtual FiniteElement< dim, spacedim >::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
void initialize_quad_dof_index_permutation()
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
std::vector< internal::FEValues::FiniteElementRelatedData< dim, spacedim > > base_fe_output_objects
void set_fe_data(const unsigned int base_no, typename FiniteElement< dim, spacedim >::InternalDataBase *)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual Point< dim > unit_support_point(const unsigned int index) const
#define AssertThrow(cond, exc)
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
const unsigned int dofs_per_line
size_type n_elements() const
InternalData(const unsigned int n_base_elements)
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
std::vector< Point< dim > > unit_support_points
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual std::string get_name() const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual FiniteElement< dim, spacedim > * clone() const
virtual std::size_t memory_consumption() const
std::vector< std::vector< FullMatrix< double > > > prolongation
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
Third derivatives of shape functions.
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
unsigned int element_multiplicity(const unsigned int index) const
static const unsigned int invalid_face_number
Abstract base class for mapping classes.
std::vector< Point< dim-1 > > unit_face_support_points
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
const ComponentMask & get_nonzero_components(const unsigned int i) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
void initialize_unit_face_support_points()
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual std::string get_name() const =0
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
void compute_fill(const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim_1 > &quadrature, const CellSimilarity::Similarity cell_similarity, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_data, const internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
unsigned int size() const
const unsigned int dofs_per_cell
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual std::size_t memory_consumption() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
std::vector< typename FiniteElement< dim, spacedim >::InternalDataBase * > base_fe_datas
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index) const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
unsigned int n_components() const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
FiniteElement< dim, spacedim >::InternalDataBase & get_fe_data(const unsigned int base_no) const
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
Shape function gradients.
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
const Conformity conforming_space
FESystem(const FiniteElement< dim, spacedim > &fe, const unsigned int n_elements)
void push_back(const size_type size)
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
const unsigned int dofs_per_face
const std::vector< ComponentMask > nonzero_components
virtual FiniteElement< dim, spacedim >::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
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
const unsigned int dofs_per_vertex
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
unsigned int size(const unsigned int i) const
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
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
BlockIndices base_to_block_indices
internal::FEValues::FiniteElementRelatedData< dim, spacedim > & get_fe_output_object(const unsigned int base_no) const
unsigned int n_base_elements() const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const BlockIndices & block_indices() const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
TableIndices< 2 > interface_constraints_size() const
Task< RT > new_task(const std_cxx11::function< RT()> &function)
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
unsigned int tensor_degree() const