54 zero_indices(
unsigned int (&indices)[dim])
56 for (
unsigned int d = 0; d < dim; ++d)
67 increment_indices(
unsigned int (&indices)[dim],
const unsigned int dofs1d)
70 for (
unsigned int d = 0; d < dim - 1; ++d)
71 if (indices[d] == dofs1d)
87template <
int xdim,
int xspacedim>
94 template <
int spacedim>
103 template <
int spacedim>
108 const unsigned int dim = 2;
110 unsigned int q_deg = fe.
degree;
159 std::vector<
Point<dim - 1>> constraint_points;
161 constraint_points.emplace_back(0.5);
165 const unsigned int n = q_deg - 1;
166 const double step = 1. / q_deg;
168 for (
unsigned int i = 1; i <= n; ++i)
169 constraint_points.push_back(
173 for (
unsigned int i = 1; i <= n; ++i)
174 constraint_points.push_back(
187 const std::vector<unsigned int> &index_map_inverse =
189 const std::vector<unsigned int> face_index_map =
196 for (
unsigned int i = 0; i < constraint_points.size(); ++i)
197 for (
unsigned int j = 0; j < q_deg + 1; ++j)
200 p[0] = constraint_points[i](0);
202 fe.
poly_space->compute_value(index_map_inverse[j], p);
214 template <
int spacedim>
219 const unsigned int dim = 3;
221 unsigned int q_deg = fe.
degree;
241 std::vector<
Point<dim - 1>> constraint_points;
244 constraint_points.emplace_back(0.5, 0.5);
247 constraint_points.emplace_back(0, 0.5);
248 constraint_points.emplace_back(1, 0.5);
249 constraint_points.emplace_back(0.5, 0);
250 constraint_points.emplace_back(0.5, 1);
254 const unsigned int n = q_deg - 1;
255 const double step = 1. / q_deg;
256 std::vector<
Point<dim - 2>> line_support_points(n);
257 for (
unsigned int i = 0; i < n; ++i)
258 line_support_points[i](0) = (i + 1) * step;
259 Quadrature<dim - 2> qline(line_support_points);
262 std::vector<
Point<dim - 1>> p_line(n);
268 ReferenceCells::get_hypercube<dim - 1>(), qline, 0, 0, p_line);
269 for (
unsigned int i = 0; i < n; ++i)
273 ReferenceCells::get_hypercube<dim - 1>(), qline, 0, 1, p_line);
274 for (
unsigned int i = 0; i < n; ++i)
278 ReferenceCells::get_hypercube<dim - 1>(), qline, 2, 0, p_line);
279 for (
unsigned int i = 0; i < n; ++i)
283 ReferenceCells::get_hypercube<dim - 1>(), qline, 2, 1, p_line);
284 for (
unsigned int i = 0; i < n; ++i)
288 for (
unsigned int face = 0;
291 for (
unsigned int subface = 0;
296 ReferenceCells::get_hypercube<dim - 1>(),
301 constraint_points.insert(constraint_points.end(),
307 std::vector<
Point<dim - 1>> inner_points(n * n);
308 for (
unsigned int i = 0, iy = 1; iy <= n; ++iy)
309 for (
unsigned int ix = 1; ix <= n; ++ix)
313 for (
unsigned int child = 0;
316 for (
const auto &inner_point : inner_points)
317 constraint_points.push_back(
324 const unsigned int pnts = (q_deg + 1) * (q_deg + 1);
328 const std::vector<unsigned int> &index_map_inverse =
330 const std::vector<unsigned int> face_index_map =
340 for (
unsigned int i = 0; i < constraint_points.size(); ++i)
342 const double interval =
static_cast<double>(q_deg * 2);
343 bool mirror[dim - 1];
355 for (
unsigned int k = 0; k < dim - 1; ++k)
357 const int coord_int =
358 static_cast<int>(constraint_points[i](k) * interval + 0.25);
359 constraint_point(k) = 1. * coord_int / interval;
381 mirror[k] = (constraint_point(k) > 0.5);
383 constraint_point(k) = 1.0 - constraint_point(k);
386 for (
unsigned int j = 0; j < pnts; ++j)
388 unsigned int indices[2] = {j % (q_deg + 1), j / (q_deg + 1)};
390 for (
unsigned int k = 0; k < 2; ++k)
392 indices[k] = q_deg - indices[k];
394 const unsigned int new_index =
395 indices[1] * (q_deg + 1) + indices[0];
398 fe.
poly_space->compute_value(index_map_inverse[new_index],
415template <
int dim,
int spacedim>
419 const std::vector<bool> & restriction_is_additive_flags)
423 restriction_is_additive_flags,
426 &poly_space) != nullptr ?
433template <
int dim,
int spacedim>
438 ExcMessage(
"The first support point has to be zero."));
439 Assert(points.back()[0] == 1,
440 ExcMessage(
"The last support point has to be one."));
445 const unsigned int q_dofs_per_cell =
446 Utilities::fixed_power<dim>(
q_degree + 1);
452 [
this, q_dofs_per_cell]() {
453 std::vector<unsigned int> renumber =
454 FETools::hierarchic_to_lexicographic_numbering<dim>(
q_degree);
455 for (
unsigned int i = q_dofs_per_cell; i < this->
n_dofs_per_cell(); ++i)
456 renumber.push_back(i);
457 auto *tensor_poly_space_ptr =
459 if (tensor_poly_space_ptr !=
nullptr)
464 auto *tensor_piecewise_poly_space_ptr =
dynamic_cast<
467 if (tensor_piecewise_poly_space_ptr !=
nullptr)
472 auto *tensor_bubbles_poly_space_ptr =
475 if (tensor_bubbles_poly_space_ptr !=
nullptr)
480 auto *tensor_const_poly_space_ptr =
483 if (tensor_const_poly_space_ptr !=
nullptr)
509template <
int dim,
int spacedim>
520 Assert(interpolation_matrix.
m() == this->n_dofs_per_cell(),
522 this->n_dofs_per_cell()));
528 const unsigned int q_dofs_per_cell =
529 Utilities::fixed_power<dim>(
q_degree + 1);
530 const unsigned int source_q_dofs_per_cell =
531 Utilities::fixed_power<dim>(source_fe->degree + 1);
537 for (
unsigned int j = 0; j < q_dofs_per_cell; ++j)
547 for (
unsigned int i = 0; i < source_q_dofs_per_cell; ++i)
548 interpolation_matrix(j, i) =
549 source_fe->poly_space->compute_value(i, p);
556 source_fe->n_dofs_per_cell());
557 for (
unsigned int i = 0; i < source_q_dofs_per_cell; ++i)
558 interpolation_matrix(q_dofs_per_cell, i) = 0.;
559 for (
unsigned int j = 0; j < q_dofs_per_cell; ++j)
560 interpolation_matrix(j, source_q_dofs_per_cell) = 0.;
561 interpolation_matrix(q_dofs_per_cell, source_q_dofs_per_cell) = 1.;
567 for (
unsigned int j = 0; j < source_fe->n_dofs_per_cell(); ++j)
568 if (std::fabs(interpolation_matrix(i, j)) <
eps)
569 interpolation_matrix(i, j) = 0.;
577 for (
unsigned int j = 0; j < source_fe->n_dofs_per_cell(); ++j)
578 sum += interpolation_matrix(i, j);
609template <
int dim,
int spacedim>
614 const unsigned int face_no)
const
618 interpolation_matrix,
624template <
int dim,
int spacedim>
628 const unsigned int subface,
630 const unsigned int face_no)
const
643 Assert(interpolation_matrix.
n() == this->n_dofs_per_face(face_no),
645 this->n_dofs_per_face(face_no)));
686 double matrix_entry =
692 if (std::fabs(matrix_entry - 1.0) < eps)
694 if (std::fabs(matrix_entry) < eps)
697 interpolation_matrix(i, j) = matrix_entry;
709 sum += interpolation_matrix(j, i);
728template <
int dim,
int spacedim>
737template <
int dim,
int spacedim>
738std::vector<std::pair<unsigned int, unsigned int>>
781template <
int dim,
int spacedim>
782std::vector<std::pair<unsigned int, unsigned int>>
799 const unsigned int p = this->
degree;
800 const unsigned int q = fe_q_other->degree;
802 std::vector<std::pair<unsigned int, unsigned int>> identities;
804 const std::vector<unsigned int> &index_map_inverse =
806 const std::vector<unsigned int> &index_map_inverse_other =
807 fe_q_other->get_poly_space_numbering_inverse();
809 for (
unsigned int i = 0; i < p - 1; ++i)
810 for (
unsigned int j = 0; j < q - 1; ++j)
813 fe_q_other->unit_support_points[index_map_inverse_other[j + 1]]
815 identities.emplace_back(i, j);
836 const std::vector<unsigned int> &index_map_inverse_q =
839 std::vector<std::pair<unsigned int, unsigned int>> identities;
841 for (
unsigned int i = 0; i < this->
degree - 1; ++i)
842 for (
unsigned int j = 0; j < fe_p_other->degree - 1; ++j)
845 fe_p_other->get_unit_support_points()[j + 3][0]) < 1e-14)
846 identities.emplace_back(i, j);
876template <
int dim,
int spacedim>
877std::vector<std::pair<unsigned int, unsigned int>>
880 const unsigned int)
const
892 const unsigned int p = this->
degree;
893 const unsigned int q = fe_q_other->degree;
895 std::vector<std::pair<unsigned int, unsigned int>> identities;
897 const std::vector<unsigned int> &index_map_inverse =
899 const std::vector<unsigned int> &index_map_inverse_other =
900 fe_q_other->get_poly_space_numbering_inverse();
902 for (
unsigned int i1 = 0; i1 < p - 1; ++i1)
903 for (
unsigned int i2 = 0; i2 < p - 1; ++i2)
904 for (
unsigned int j1 = 0; j1 < q - 1; ++j1)
905 for (
unsigned int j2 = 0; j2 < q - 1; ++j2)
916 identities.emplace_back(i1 * (p - 1) + i2, j1 * (q - 1) + j2);
924 return std::vector<std::pair<unsigned int, unsigned int>>();
935 return std::vector<std::pair<unsigned int, unsigned int>>();
940 return std::vector<std::pair<unsigned int, unsigned int>>();
952template <
int dim,
int spacedim>
955 const std::vector<
Point<1>> &points)
957 const std::vector<unsigned int> &index_map_inverse =
970 for (
unsigned int k = 0; k < support_quadrature.size(); ++k)
972 support_quadrature.point(k);
977template <
int dim,
int spacedim>
980 const std::vector<
Point<1>> &points)
985 const unsigned int face_no = 0;
988 Utilities::fixed_power<dim - 1>(
q_degree + 1));
996 const std::vector<unsigned int> face_index_map =
1004 const Quadrature<dim - 1> support_quadrature(support_1d);
1009 for (
unsigned int k = 0; k < support_quadrature.size(); ++k)
1011 support_quadrature.point(k);
1016template <
int dim,
int spacedim>
1027 const unsigned int face_no = 0;
1033 const unsigned int n =
q_degree - 1;
1054 for (
unsigned int local = 0; local < this->
n_dofs_per_quad(face_no); ++local)
1058 unsigned int i = local % n, j = local / n;
1067 i + (n - 1 - j) * n - local;
1071 (n - 1 - j) + (n - 1 - i) * n - local;
1075 (n - 1 - i) + j * n - local;
1082 j + (n - 1 - i) * n - local;
1086 (n - 1 - i) + (n - 1 - j) * n - local;
1090 (n - 1 - j) + i * n - local;
1101template <
int dim,
int spacedim>
1104 const unsigned int face,
1105 const bool face_orientation,
1106 const bool face_flip,
1107 const bool face_rotation)
const
1129 const unsigned int dof_index_on_vertex =
1135 face, face_vertex, face_orientation, face_flip, face_rotation) *
1137 dof_index_on_vertex);
1144 const unsigned int index =
1152 unsigned int adjusted_dof_index_on_line = 0;
1162 if (face_flip ==
false)
1163 adjusted_dof_index_on_line = dof_index_on_line;
1165 adjusted_dof_index_on_line =
1178 ((face_orientation ==
true) && (face_flip ==
false) &&
1179 (face_rotation ==
false)),
1181 adjusted_dof_index_on_line = dof_index_on_line;
1190 face, face_line, face_orientation, face_flip, face_rotation) *
1192 adjusted_dof_index_on_line);
1200 const unsigned int index =
1207 ((face_orientation ==
true) && (face_flip ==
false) &&
1208 (face_rotation ==
false)),
1216template <
int dim,
int spacedim>
1217std::vector<unsigned int>
1222 std::vector<unsigned int> dpo(dim + 1, 1U);
1223 for (
unsigned int i = 1; i < dpo.size(); ++i)
1224 dpo[i] = dpo[i - 1] * (
degree - 1);
1230template <
int dim,
int spacedim>
1233 const std::vector<
Point<1>> &points)
1240template <
int dim,
int spacedim>
1243 const unsigned int child,
1250 "Prolongation matrices are only available for refined cells!"));
1254 if (this->
prolongation[refinement_case - 1][child].n() == 0)
1256 std::lock_guard<std::mutex> lock(this->
mutex);
1259 if (this->
prolongation[refinement_case - 1][child].n() ==
1264 const unsigned int q_dofs_per_cell =
1265 Utilities::fixed_power<dim>(
q_degree + 1);
1278 for (
unsigned int i = 0; i < q_dofs_per_cell; ++i)
1281 i, this->unit_support_points[i])) < eps,
1283 "to one or zero in a nodal point. "
1284 "This typically indicates that the "
1285 "polynomial interpolation is "
1286 "ill-conditioned such that round-off "
1287 "prevents the sum to be one."));
1288 for (
unsigned int j = 0; j < q_dofs_per_cell; ++j)
1291 i, this->unit_support_points[j])) < eps,
1293 "The Lagrange polynomial does not evaluate "
1294 "to one or zero in a nodal point. "
1295 "This typically indicates that the "
1296 "polynomial interpolation is "
1297 "ill-conditioned such that round-off "
1298 "prevents the sum to be one."));
1306 const unsigned int dofs1d =
q_degree + 1;
1307 std::vector<Table<2, double>> subcell_evaluations(
1310 const std::vector<unsigned int> &index_map_inverse =
1315 unsigned int step_size_diag = 0;
1317 unsigned int factor = 1;
1318 for (
unsigned int d = 0;
d < dim; ++
d)
1320 step_size_diag += factor;
1330 for (
unsigned int j = 0; j < dofs1d; ++j)
1332 const unsigned int diag_comp = index_map_inverse[j * step_size_diag];
1338 for (
unsigned int i = 0; i < dofs1d; ++i)
1339 for (
unsigned int d = 0;
d < dim; ++
d)
1344 const double cell_value =
1345 this->
poly_space->compute_value(index_map_inverse[i], point);
1364 if (std::fabs(cell_value) < eps)
1365 subcell_evaluations[
d](j, i) = 0;
1367 subcell_evaluations[
d](j, i) = cell_value;
1373 unsigned int j_indices[dim];
1374 internal::FE_Q_Base::zero_indices<dim>(j_indices);
1375 for (
unsigned int j = 0; j < q_dofs_per_cell; j += dofs1d)
1377 unsigned int i_indices[dim];
1378 internal::FE_Q_Base::zero_indices<dim>(i_indices);
1379 for (
unsigned int i = 0; i < q_dofs_per_cell; i += dofs1d)
1381 double val_extra_dim = 1.;
1382 for (
unsigned int d = 1;
d < dim; ++
d)
1384 subcell_evaluations[d](j_indices[d - 1], i_indices[d - 1]);
1388 for (
unsigned int jj = 0; jj < dofs1d; ++jj)
1390 const unsigned int j_ind = index_map_inverse[j + jj];
1391 for (
unsigned int ii = 0; ii < dofs1d; ++ii)
1392 prolongate(j_ind, index_map_inverse[i + ii]) =
1393 val_extra_dim * subcell_evaluations[0](jj, ii);
1398 internal::FE_Q_Base::increment_indices<dim>(i_indices, dofs1d);
1401 internal::FE_Q_Base::increment_indices<dim>(j_indices, dofs1d);
1407 prolongate(q_dofs_per_cell, q_dofs_per_cell) = 1.;
1416 sum += prolongate(row, col);
1417 Assert(std::fabs(sum - 1.) <
1420 "prolongation matrix do not add to one. "
1421 "This typically indicates that the "
1422 "polynomial interpolation is "
1423 "ill-conditioned such that round-off "
1424 "prevents the sum to be one."));
1439template <
int dim,
int spacedim>
1442 const unsigned int child,
1449 "Restriction matrices are only available for refined cells!"));
1453 if (this->
restriction[refinement_case - 1][child].n() == 0)
1455 std::lock_guard<std::mutex> lock(this->
mutex);
1458 if (this->
restriction[refinement_case - 1][child].n() ==
1460 return this->
restriction[refinement_case - 1][child];
1465 const unsigned int q_dofs_per_cell =
1466 Utilities::fixed_power<dim>(
q_degree + 1);
1487 const std::vector<unsigned int> &index_map_inverse =
1490 const unsigned int dofs1d =
q_degree + 1;
1491 std::vector<Tensor<1, dim>> evaluations1d(dofs1d);
1495 for (
unsigned int i = 0; i < q_dofs_per_cell; ++i)
1497 unsigned int mother_dof = index_map_inverse[i];
1512 for (
unsigned int j = 0; j < dofs1d; ++j)
1513 for (
unsigned int d = 0;
d < dim; ++
d)
1517 evaluations1d[j][
d] =
1518 this->
poly_space->compute_value(index_map_inverse[j],
1521 unsigned int j_indices[dim];
1522 internal::FE_Q_Base::zero_indices<dim>(j_indices);
1524 double sum_check = 0;
1526 for (
unsigned int j = 0; j < q_dofs_per_cell; j += dofs1d)
1528 double val_extra_dim = 1.;
1529 for (
unsigned int d = 1;
d < dim; ++
d)
1530 val_extra_dim *= evaluations1d[j_indices[d - 1]][d];
1531 for (
unsigned int jj = 0; jj < dofs1d; ++jj)
1539 const double val = val_extra_dim * evaluations1d[jj][0];
1540 const unsigned int child_dof = index_map_inverse[j + jj];
1541 if (std::fabs(val - 1.) < eps)
1542 my_restriction(mother_dof, child_dof) = 1.;
1543 else if (std::fabs(val) > eps)
1544 my_restriction(mother_dof, child_dof) = val;
1549 internal::FE_Q_Base::increment_indices<dim>(j_indices,
1552 Assert(std::fabs(sum_check - 1) <
1556 "restriction matrix do not add to one. "
1557 "This typically indicates that the "
1558 "polynomial interpolation is "
1559 "ill-conditioned such that round-off "
1560 "prevents the sum to be one."));
1577 return this->
restriction[refinement_case - 1][child];
1587template <
int dim,
int spacedim>
1590 const unsigned int shape_index,
1591 const unsigned int face_index)
const
1600 return (((shape_index == 0) && (face_index == 0)) ||
1601 ((shape_index == 1) && (face_index == 1)));
1616 const unsigned int vertex_no = shape_index;
1620 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_face; ++v)
1630 const unsigned int line_index =
1637 return (line_index == face_index);
1641 const unsigned int lines_per_face =
1644 for (
unsigned int l = 0;
l < lines_per_face; ++
l)
1657 const unsigned int quad_index =
1660 Assert(
static_cast<signed int>(quad_index) <
1670 return (quad_index == face_index);
1689template <
int dim,
int spacedim>
1690std::pair<Table<2, bool>, std::vector<unsigned int>>
1698 for (
unsigned int i = 0; i < Utilities::fixed_power<dim>(
q_degree + 1); ++i)
1699 constant_modes(0, i) =
true;
1700 return std::pair<Table<2, bool>, std::vector<unsigned int>>(
1701 constant_modes, std::vector<unsigned int>(1, 0));
1707#include "fe_q_base.inst"
std::vector< unsigned int > get_poly_space_numbering_inverse() const
const ScalarPolynomialsBase< dim > & get_poly_space() const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
const std::unique_ptr< ScalarPolynomialsBase< dim > > poly_space
void initialize_unit_face_support_points(const std::vector< Point< 1 > > &points)
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) 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_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
const unsigned int q_degree
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
void initialize_quad_dof_index_permutation()
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
void initialize(const std::vector< Point< 1 > > &support_points_1d)
void initialize_unit_support_points(const std::vector< Point< 1 > > &points)
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
void initialize_constraints(const std::vector< Point< 1 > > &points)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
FE_Q_Base(const ScalarPolynomialsBase< dim > &poly_space, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags)
unsigned int get_first_line_index() const
unsigned int n_dofs_per_vertex() const
const unsigned int degree
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_components() const
ReferenceCell reference_cell() const
unsigned int get_first_face_line_index(const unsigned int face_no=0) const
unsigned int get_first_face_quad_index(const unsigned int face_no=0) const
unsigned int n_unique_faces() const
unsigned int n_dofs_per_quad(unsigned int face_no=0) const
unsigned int get_first_hex_index() const
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
bool has_face_support_points(const unsigned int face_no=0) const
std::vector< std::vector< FullMatrix< double > > > restriction
std::vector< Table< 2, int > > adjust_quad_dof_index_for_face_orientation_table
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
TableIndices< 2 > interface_constraints_size() const
std::vector< Point< dim > > unit_support_points
FullMatrix< double > interface_constraints
std::vector< std::vector< FullMatrix< double > > > prolongation
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const Point< dim > & point(const unsigned int i) const
void set_numbering(const std::vector< unsigned int > &renumber)
void set_numbering(const std::vector< unsigned int > &renumber)
void set_numbering(const std::vector< unsigned int > &renumber)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcInterpolationNotImplemented()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Task< RT > new_task(const std::function< RT()> &function)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T sum(const T &t, const MPI_Comm &mpi_communicator)
static const unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
static void initialize_constraints(const std::vector< Point< 1 > > &, FE_Q_Base< 3, spacedim > &fe)
static void initialize_constraints(const std::vector< Point< 1 > > &, FE_Q_Base< 1, spacedim > &)
static void initialize_constraints(const std::vector< Point< 1 > > &, FE_Q_Base< 2, spacedim > &fe)
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static Point< dim > cell_to_child_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)