39 count_nonzeros(
const std::vector<unsigned int> &vec)
41 return std::count_if(vec.begin(), vec.end(), [](
const unsigned int i) {
49 template <
int dim,
int spacedim>
51 const unsigned int n_base_elements)
52 : base_fe_datas(n_base_elements)
53 , base_fe_output_objects(n_base_elements)
58 template <
int dim,
int spacedim>
67 template <
int dim,
int spacedim>
70 const unsigned int base_no)
const 78 template <
int dim,
int spacedim>
81 const unsigned int base_no,
90 template <
int dim,
int spacedim>
93 const unsigned int base_no)
const 104 template <
int dim,
int spacedim>
108 template <
int dim,
int spacedim>
110 const unsigned int n_elements)
118 std::vector<const FiniteElement<dim, spacedim> *> fes;
120 std::vector<unsigned int> multiplicities;
121 multiplicities.push_back(n_elements);
127 template <
int dim,
int spacedim>
129 const unsigned int n1,
131 const unsigned int n2)
141 std::vector<const FiniteElement<dim, spacedim> *> fes;
144 std::vector<unsigned int> multiplicities;
145 multiplicities.push_back(n1);
146 multiplicities.push_back(n2);
152 template <
int dim,
int spacedim>
154 const unsigned int n1,
156 const unsigned int n2,
158 const unsigned int n3)
175 std::vector<const FiniteElement<dim, spacedim> *> fes;
179 std::vector<unsigned int> multiplicities;
180 multiplicities.push_back(n1);
181 multiplicities.push_back(n2);
182 multiplicities.push_back(n3);
188 template <
int dim,
int spacedim>
190 const unsigned int n1,
192 const unsigned int n2,
194 const unsigned int n3,
196 const unsigned int n4)
224 std::vector<const FiniteElement<dim, spacedim> *> fes;
229 std::vector<unsigned int> multiplicities;
230 multiplicities.push_back(n1);
231 multiplicities.push_back(n2);
232 multiplicities.push_back(n3);
233 multiplicities.push_back(n4);
239 template <
int dim,
int spacedim>
241 const unsigned int n1,
243 const unsigned int n2,
245 const unsigned int n3,
247 const unsigned int n4,
249 const unsigned int n5)
252 multiply_dof_numbers(&fe1, n1, &fe2, n2, &fe3, n3, &fe4, n4, &fe5, n5),
273 ,
base_elements((n1 > 0) + (n2 > 0) + (n3 > 0) + (n4 > 0) + (n5 > 0))
275 std::vector<const FiniteElement<dim, spacedim> *> fes;
281 std::vector<unsigned int> multiplicities;
282 multiplicities.push_back(n1);
283 multiplicities.push_back(n2);
284 multiplicities.push_back(n3);
285 multiplicities.push_back(n4);
286 multiplicities.push_back(n5);
292 template <
int dim,
int spacedim>
295 const std::vector<unsigned int> & multiplicities)
309 template <
int dim,
int spacedim>
320 std::ostringstream namebuf;
333 return namebuf.str();
338 template <
int dim,
int spacedim>
339 std::unique_ptr<FiniteElement<dim, spacedim>>
342 std::vector<const FiniteElement<dim, spacedim> *> fes;
343 std::vector<unsigned int> multiplicities;
350 return std_cxx14::make_unique<FESystem<dim, spacedim>>(fes, multiplicities);
355 template <
int dim,
int spacedim>
358 const unsigned int first_component,
359 const unsigned int n_selected_components)
const 362 ExcMessage(
"Invalid arguments (not a part of this FiniteElement)."));
364 const unsigned int base_index =
366 const unsigned int component_in_base =
368 const unsigned int base_components =
373 if (n_selected_components <= base_components)
375 .get_sub_fe(component_in_base, n_selected_components);
378 ExcMessage(
"You can not select a part of a FiniteElement."));
384 template <
int dim,
int spacedim>
395 .shape_value(this->system_to_base_table[i].second, p));
400 template <
int dim,
int spacedim>
403 const unsigned int i,
405 const unsigned int component)
const 421 const unsigned int component_in_base =
436 template <
int dim,
int spacedim>
447 .shape_grad(this->system_to_base_table[i].second, p));
452 template <
int dim,
int spacedim>
455 const unsigned int i,
457 const unsigned int component)
const 469 const unsigned int component_in_base =
481 template <
int dim,
int spacedim>
492 .shape_grad_grad(this->system_to_base_table[i].second, p));
497 template <
int dim,
int spacedim>
500 const unsigned int i,
502 const unsigned int component)
const 514 const unsigned int component_in_base =
526 template <
int dim,
int spacedim>
537 .shape_3rd_derivative(this->system_to_base_table[i].second, p));
542 template <
int dim,
int spacedim>
545 const unsigned int i,
547 const unsigned int component)
const 559 const unsigned int component_in_base =
571 template <
int dim,
int spacedim>
582 .shape_4th_derivative(this->system_to_base_table[i].second, p));
587 template <
int dim,
int spacedim>
590 const unsigned int i,
592 const unsigned int component)
const 604 const unsigned int component_in_base =
616 template <
int dim,
int spacedim>
640 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
667 std::vector<FullMatrix<double>> base_matrices(this->
n_base_elements());
679 interpolation_matrix = 0;
684 interpolation_matrix(i, j) =
692 template <
int dim,
int spacedim>
695 const unsigned int child,
702 "Restriction matrices are only available for refined cells!"));
706 if (this->
restriction[refinement_case - 1][child].n() == 0)
708 std::lock_guard<std::mutex> lock(this->
mutex);
711 if (this->
restriction[refinement_case - 1][child].n() ==
713 return this->
restriction[refinement_case - 1][child];
716 bool do_restriction =
true;
719 std::vector<const FullMatrix<double> *> base_matrices(
725 &
base_element(i).get_restriction_matrix(child, refinement_case);
727 do_restriction =
false;
759 const unsigned int base =
762 const unsigned int base_index_i =
770 (*base_matrices[base])(base_index_i, base_index_j);
774 this->restriction[refinement_case - 1][child]));
778 return this->
restriction[refinement_case - 1][child];
783 template <
int dim,
int spacedim>
786 const unsigned int child,
793 "Restriction matrices are only available for refined cells!"));
798 if (this->
prolongation[refinement_case - 1][child].n() == 0)
800 std::lock_guard<std::mutex> lock(this->
mutex);
802 if (this->
prolongation[refinement_case - 1][child].n() ==
806 bool do_prolongation =
true;
807 std::vector<const FullMatrix<double> *> base_matrices(
812 &
base_element(i).get_prolongation_matrix(child, refinement_case);
814 do_prolongation =
false;
830 const unsigned int base =
833 const unsigned int base_index_i =
838 (*base_matrices[base])(base_index_i, base_index_j);
849 template <
int dim,
int spacedim>
852 const unsigned int face,
853 const bool face_orientation,
854 const bool face_flip,
855 const bool face_rotation)
const 860 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
863 const unsigned int base_face_to_cell_index =
865 .face_to_cell_index(face_base_index.second,
876 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> target =
877 std::make_pair(face_base_index.first, base_face_to_cell_index);
894 template <
int dim,
int spacedim>
901 for (
unsigned int base_no = 0; base_no < this->
n_base_elements(); ++base_no)
902 out |=
base_element(base_no).requires_update_flags(flags);
908 template <
int dim,
int spacedim>
909 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
925 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
926 data_ptr = std_cxx14::make_unique<InternalData>(this->
n_base_elements());
942 for (
unsigned int base_no = 0; base_no < this->
n_base_elements(); ++base_no)
945 &base_fe_output_object = data.get_fe_output_object(base_no);
949 flags |
base_element(base_no).requires_update_flags(flags));
957 auto base_fe_data =
base_element(base_no).get_data(flags,
960 base_fe_output_object);
962 data.set_fe_data(base_no, std::move(base_fe_data));
971 template <
int dim,
int spacedim>
972 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
988 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
989 data_ptr = std_cxx14::make_unique<InternalData>(this->
n_base_elements());
1005 for (
unsigned int base_no = 0; base_no < this->
n_base_elements(); ++base_no)
1008 &base_fe_output_object = data.get_fe_output_object(base_no);
1012 flags |
base_element(base_no).requires_update_flags(flags));
1020 auto base_fe_data =
base_element(base_no).get_face_data(
1021 flags, mapping, quadrature, base_fe_output_object);
1023 data.set_fe_data(base_no, std::move(base_fe_data));
1034 template <
int dim,
int spacedim>
1035 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1051 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
1052 data_ptr = std_cxx14::make_unique<InternalData>(this->
n_base_elements());
1069 for (
unsigned int base_no = 0; base_no < this->
n_base_elements(); ++base_no)
1072 &base_fe_output_object = data.get_fe_output_object(base_no);
1076 flags |
base_element(base_no).requires_update_flags(flags));
1084 auto base_fe_data =
base_element(base_no).get_subface_data(
1085 flags, mapping, quadrature, base_fe_output_object);
1087 data.set_fe_data(base_no, std::move(base_fe_data));
1095 template <
int dim,
int spacedim>
1103 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1125 template <
int dim,
int spacedim>
1129 const unsigned int face_no,
1133 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1155 template <
int dim,
int spacedim>
1159 const unsigned int face_no,
1160 const unsigned int sub_no,
1164 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
1186 template <
int dim,
int spacedim>
1187 template <
int dim_1>
1192 const unsigned int face_no,
1193 const unsigned int sub_no,
1207 Assert(dynamic_cast<const InternalData *>(&fe_internal) !=
nullptr,
1230 for (
unsigned int base_no = 0; base_no < this->
n_base_elements(); ++base_no)
1242 const Quadrature<dim - 1> *face_quadrature =
nullptr;
1243 const unsigned int n_q_points = quadrature.
size();
1247 const Subscriptor *quadrature_base_pointer = &quadrature;
1253 quadrature_base_pointer) !=
nullptr,
1263 quadrature_base_pointer) !=
nullptr,
1267 static_cast<const Quadrature<dim - 1
> *>(quadrature_base_pointer);
1327 for (
unsigned int system_index = 0; system_index < this->
dofs_per_cell;
1331 const unsigned int base_index =
1341 unsigned int out_index = 0;
1342 for (
unsigned int i = 0; i < system_index; ++i)
1344 unsigned int in_index = 0;
1345 for (
unsigned int i = 0; i < base_index; ++i)
1354 for (
unsigned int s = 0;
1357 for (
unsigned int q = 0; q < n_q_points; ++q)
1359 base_data.shape_values(in_index + s, q);
1362 for (
unsigned int s = 0;
1365 for (
unsigned int q = 0; q < n_q_points; ++q)
1367 base_data.shape_gradients[in_index + s][q];
1370 for (
unsigned int s = 0;
1373 for (
unsigned int q = 0; q < n_q_points; ++q)
1375 base_data.shape_hessians[in_index + s][q];
1378 for (
unsigned int s = 0;
1381 for (
unsigned int q = 0; q < n_q_points; ++q)
1383 base_data.shape_3rd_derivatives[in_index + s][q];
1389 template <
int dim,
int spacedim>
1398 if (
base_element(base).constraints_are_implemented() ==
false)
1422 const std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1427 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> m_index;
1449 const unsigned int index_in_line =
1451 const unsigned int sub_line =
1459 const unsigned int tmp1 =
1460 2 * this->dofs_per_vertex + index_in_line;
1476 const unsigned int tmp2 =
1502 if (m < 5 * this->dofs_per_vertex + 12 * this->
dofs_per_line)
1505 const unsigned int index_in_line =
1507 const unsigned int sub_line =
1511 const unsigned int tmp1 =
1512 4 * this->dofs_per_vertex + index_in_line;
1519 const unsigned int tmp2 =
1526 5 *
base_element(m_index.first.first).dofs_per_vertex +
1535 const unsigned int index_in_quad =
1536 (m - 5 * this->dofs_per_vertex -
1541 const unsigned int sub_quad =
1542 ((m - 5 * this->dofs_per_vertex -
1544 this->dofs_per_quad);
1547 const unsigned int tmp1 = 4 * this->dofs_per_vertex +
1556 4 *
base_element(m_index.first.first).dofs_per_vertex +
1559 const unsigned int tmp2 =
1561 4 *
base_element(m_index.first.first).dofs_per_vertex -
1567 5 *
base_element(m_index.first.first).dofs_per_vertex +
1584 if (n_index.first == m_index.first)
1587 .constraints()(m_index.second, n_index.second));
1593 template <
int dim,
int spacedim>
1597 const std::vector<unsigned int> & multiplicities)
1599 Assert(fes.size() == multiplicities.size(),
1602 ExcMessage(
"Need to pass at least one finite element."));
1603 Assert(count_nonzeros(multiplicities) > 0,
1604 ExcMessage(
"You only passed FiniteElements with multiplicity 0."));
1611 for (
unsigned int i = 0; i < fes.size(); i++)
1612 if (multiplicities[i] > 0)
1618 unsigned int ind = 0;
1619 for (
unsigned int i = 0; i < fes.size(); i++)
1620 if (multiplicities[i] > 0)
1665 for (
unsigned int base_el = 0; base_el < this->
n_base_elements(); ++base_el)
1717 const unsigned int base_i =
1719 const unsigned int index_in_base =
1727 base_element(base_i).unit_face_support_points[index_in_base];
1752 if (!
base_element(base).has_generalized_support_points())
1755 for (
const auto &
point :
1787 const auto &points =
1789 for (
unsigned int j = 0; j < points.size(); ++j)
1810 unsigned int index = 0;
1815 .adjust_quad_dof_index_for_face_orientation_table;
1818 for (
unsigned int i = 0; i < temp.
size(0); ++i)
1819 for (
unsigned int j = 0; j < 8; ++j)
1821 index + i, j) = temp(i, j);
1822 index += temp.
size(0);
1834 const std::vector<int> &temp2 =
1836 .adjust_line_dof_index_for_line_orientation_table;
1844 index += temp2.size();
1851 init_tasks.join_all();
1856 template <
int dim,
int spacedim>
1861 if (
base_element(
b).hp_constraints_are_implemented() ==
false)
1869 template <
int dim,
int spacedim>
1891 if (
const auto *fe_other_system =
1895 interpolation_matrix = 0;
1899 unsigned int base_index = 0, base_index_other = 0;
1900 unsigned int multiplicity = 0, multiplicity_other = 0;
1908 fe_other_system->base_element(
1915 base_to_base_interpolation.
reinit(base_other.dofs_per_face,
1918 base_to_base_interpolation);
1925 std::make_pair(base_index, multiplicity))
1926 for (
unsigned int j = 0; j < fe_other_system->dofs_per_face; ++j)
1927 if (fe_other_system->face_system_to_base_index(j).first ==
1928 std::make_pair(base_index_other, multiplicity_other))
1929 interpolation_matrix(j, i) = base_to_base_interpolation(
1930 fe_other_system->face_system_to_base_index(j).second,
1942 ++multiplicity_other;
1943 if (multiplicity_other ==
1944 fe_other_system->element_multiplicity(base_index_other))
1946 multiplicity_other = 0;
1954 Assert(base_index_other == fe_other_system->n_base_elements(),
1961 Assert(base_index_other != fe_other_system->n_base_elements(),
1978 template <
int dim,
int spacedim>
1982 const unsigned int subface,
1986 (x_source_fe.
get_name().find(
"FESystem<") == 0) ||
2008 if (fe_other_system !=
nullptr)
2011 interpolation_matrix = 0;
2015 unsigned int base_index = 0, base_index_other = 0;
2016 unsigned int multiplicity = 0, multiplicity_other = 0;
2031 base_to_base_interpolation.
reinit(base_other.dofs_per_face,
2035 base_to_base_interpolation);
2042 std::make_pair(base_index, multiplicity))
2043 for (
unsigned int j = 0; j < fe_other_system->
dofs_per_face; ++j)
2045 std::make_pair(base_index_other, multiplicity_other))
2046 interpolation_matrix(j, i) = base_to_base_interpolation(
2059 ++multiplicity_other;
2060 if (multiplicity_other ==
2063 multiplicity_other = 0;
2086 fe_other_system !=
nullptr,
2094 template <
int dim,
int spacedim>
2095 template <
int structdim>
2096 std::vector<std::pair<unsigned int, unsigned int>>
2117 unsigned int base_index = 0, base_index_other = 0;
2118 unsigned int multiplicity = 0, multiplicity_other = 0;
2122 unsigned int dof_offset = 0, dof_offset_other = 0;
2124 std::vector<std::pair<unsigned int, unsigned int>> identities;
2130 fe_other_system->base_element(
2138 std::vector<std::pair<unsigned int, unsigned int>> base_identities;
2154 for (
const auto &base_identity : base_identities)
2155 identities.emplace_back(base_identity.first + dof_offset,
2156 base_identity.second + dof_offset_other);
2159 dof_offset += base.template n_dofs_per_object<structdim>();
2161 base_other.template n_dofs_per_object<structdim>();
2172 ++multiplicity_other;
2173 if (multiplicity_other ==
2174 fe_other_system->element_multiplicity(base_index_other))
2176 multiplicity_other = 0;
2184 Assert(base_index_other == fe_other_system->n_base_elements(),
2191 Assert(base_index_other != fe_other_system->n_base_elements(),
2200 return std::vector<std::pair<unsigned int, unsigned int>>();
2206 template <
int dim,
int spacedim>
2207 std::vector<std::pair<unsigned int, unsigned int>>
2211 return hp_object_dof_identities<0>(fe_other);
2214 template <
int dim,
int spacedim>
2215 std::vector<std::pair<unsigned int, unsigned int>>
2219 return hp_object_dof_identities<1>(fe_other);
2224 template <
int dim,
int spacedim>
2225 std::vector<std::pair<unsigned int, unsigned int>>
2229 return hp_object_dof_identities<2>(fe_other);
2234 template <
int dim,
int spacedim>
2238 const unsigned int codim)
const 2261 fe_sys_other->base_element(
b).n_components(),
2264 fe_sys_other->element_multiplicity(
b),
2271 fe_sys_other->base_element(
b), codim));
2272 domination = domination & base_domination;
2284 template <
int dim,
int spacedim>
2294 template <
int dim,
int spacedim>
2297 const unsigned int shape_index,
2298 const unsigned int face_index)
const 2301 .has_support_on_face(this->system_to_base_index(shape_index).second,
2307 template <
int dim,
int spacedim>
2328 template <
int dim,
int spacedim>
2344 .unit_face_support_point(
2350 template <
int dim,
int spacedim>
2351 std::pair<Table<2, bool>, std::vector<unsigned int>>
2361 const std::pair<Table<2, bool>, std::vector<unsigned int>> base_table =
2369 const unsigned int comp = components.size();
2370 if (constant_modes.n_rows() <
2373 Table<2, bool> new_constant_modes(comp + base_table.first.n_rows() *
2375 constant_modes.n_cols());
2376 for (
unsigned int r = 0; r < comp; ++r)
2378 new_constant_modes(r, c) = constant_modes(r, c);
2379 constant_modes.
swap(new_constant_modes);
2386 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int> ind =
2388 if (ind.first.first == i)
2389 for (
unsigned int c = 0; c < base_table.first.n_rows(); ++c)
2390 constant_modes(comp +
2391 ind.first.second * base_table.first.n_rows() + c,
2392 k) = base_table.first(c, ind.second);
2395 for (
const unsigned int c : base_table.second)
2396 components.push_back(
2400 return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
2406 template <
int dim,
int spacedim>
2410 std::vector<double> & dof_values)
const 2413 ExcMessage(
"The FESystem does not have generalized support points"));
2419 std::vector<double> base_dof_values;
2420 std::vector<Vector<double>> base_point_values;
2425 unsigned int current_vector_component = 0;
2426 for (
unsigned int base = 0; base <
base_elements.size(); ++base)
2433 const unsigned int n_base_dofs =
base_element.dofs_per_cell;
2434 const unsigned int n_base_components =
base_element.n_components();
2439 if (n_base_dofs == 0)
2441 current_vector_component += multiplicity * n_base_components;
2447 const unsigned int n_base_points =
2450 base_dof_values.resize(n_base_dofs);
2451 base_point_values.resize(n_base_points);
2453 for (
unsigned int m = 0; m < multiplicity;
2454 ++m, current_vector_component += n_base_components)
2458 for (
unsigned int j = 0; j < base_point_values.size(); ++j)
2460 base_point_values[j].reinit(n_base_components,
false);
2468 std::begin(point_values[n]) + current_vector_component;
2469 const auto end =
begin + n_base_components;
2474 .convert_generalized_support_point_values_to_dof_values(
2475 base_point_values, base_dof_values);
2486 std::make_pair(base, m))
2500 for (
unsigned int m = 0; m < multiplicity; ++m)
2503 std::make_pair(base, m))
2504 dof_values[i] = std::numeric_limits<double>::signaling_NaN();
2506 current_vector_component += multiplicity * n_base_components;
2513 template <
int dim,
int spacedim>
2530 #include "fe_system.inst"
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual std::size_t memory_consumption() const override
virtual Point< dim > unit_support_point(const unsigned int index) const override
std::vector< internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > > base_fe_output_objects
static const unsigned int invalid_unsigned_int
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
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)
virtual std::string get_name() const override
std::vector< Point< dim > > generalized_support_points
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
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
unsigned int n_nonzero_components(const unsigned int i) const
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const override
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
Task< RT > new_task(const std::function< RT()> &function)
const std::vector< Point< dim > > & get_generalized_support_points() const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double >> &support_point_values, std::vector< double > &dof_values) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
#define AssertIndexRange(index, range)
const unsigned int dofs_per_quad
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const override
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
bool has_generalized_support_points() const
std::vector< Point< dim - 1 > > unit_face_support_points
#define AssertThrow(cond, exc)
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
bool is_primitive() 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 override
static ::ExceptionBase & ExcInterpolationNotImplemented()
virtual const FiniteElement< dim, spacedim > & get_sub_fe(const unsigned int first_component, const unsigned int n_selected_components) const override
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
const unsigned int dofs_per_line
InternalData(const unsigned int n_base_elements)
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< Point< dim > > unit_support_points
std::vector< std::vector< std::size_t > > generalized_support_points_index_table
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
virtual std::size_t memory_consumption() const
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
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.
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
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
Abstract base class for mapping classes.
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
#define DEAL_II_NAMESPACE_CLOSE
void set_fe_data(const unsigned int base_no, std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase >)
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::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 override
VectorType::value_type * end(VectorType &V)
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::string get_name() const =0
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.
virtual bool hp_constraints_are_implemented() const override
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
std::string dim_string(const int dim, const int spacedim)
unsigned int size() const
const unsigned int dofs_per_cell
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::vector< std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > > base_fe_datas
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
size_type size(const unsigned int i) const
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index) 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 const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
#define DEAL_II_NAMESPACE_OPEN
VectorType::value_type * begin(VectorType &V)
Shape function gradients.
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
FESystem(const FiniteElement< dim, spacedim > &fe, const unsigned int n_elements)
void push_back(const size_type size)
const unsigned int dofs_per_face
const std::vector< ComponentMask > nonzero_components
internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > & get_fe_output_object(const unsigned int base_no) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
static ::ExceptionBase & ExcNotImplemented()
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
const unsigned int dofs_per_vertex
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
BlockIndices base_to_block_indices
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
unsigned int n_base_elements() const
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
std::vector< int > adjust_line_dof_index_for_line_orientation_table
void copy(const T *begin, const T *end, U *dest)
TableIndices< 2 > interface_constraints_size() const
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::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 override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
std::vector< std::pair< std::unique_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index) const override
static ::ExceptionBase & ExcInternalError()