Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
fe_abf.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2003 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
19#include <deal.II/base/table.h>
21
23
24#include <deal.II/fe/fe.h>
25#include <deal.II/fe/fe_abf.h>
26#include <deal.II/fe/fe_tools.h>
28#include <deal.II/fe/mapping.h>
29
30#include <deal.II/grid/tria.h>
32
33#include <iostream>
34#include <memory>
35#include <sstream>
36
38
39// TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
40// adjust_line_dof_index_for_line_orientation_table fields, and write tests
41// similar to bits/face_orientation_and_fe_q_*
42
43template <int dim>
44FE_ABF<dim>::FE_ABF(const unsigned int deg)
45 : FE_PolyTensor<dim>(
46 PolynomialsABF<dim>(deg),
47 FiniteElementData<dim>(get_dpo_vector(deg),
48 dim,
49 deg + 2,
50 FiniteElementData<dim>::Hdiv),
51 std::vector<bool>(PolynomialsABF<dim>::n_polynomials(deg), true),
52 std::vector<ComponentMask>(PolynomialsABF<dim>::n_polynomials(deg),
53 ComponentMask(std::vector<bool>(dim, true))))
54 , rt_order(deg)
55{
56 Assert(dim >= 2, ExcImpossibleInDim(dim));
57 const unsigned int n_dofs = this->n_dofs_per_cell();
58
60 // First, initialize the
61 // generalized support points and
62 // quadrature weights, since they
63 // are required for interpolation.
65
66 // Now compute the inverse node matrix, generating the correct
67 // basis functions from the raw ones. For a discussion of what
68 // exactly happens here, see FETools::compute_node_matrix.
70 this->inverse_node_matrix.reinit(n_dofs, n_dofs);
72 // From now on, the shape functions provided by FiniteElement::shape_value
73 // and similar functions will be the correct ones, not
74 // the raw shape functions from the polynomial space anymore.
75
76 // Reinit the vectors of
77 // restriction and prolongation
78 // matrices to the right sizes.
79 // Restriction only for isotropic
80 // refinement
82 // Fill prolongation matrices with embedding operators
83 FETools::compute_embedding_matrices(*this, this->prolongation, false, 1.e-10);
84
86
87 // TODO: the implementation makes the assumption that all faces have the
88 // same number of dofs
90 const unsigned int face_no = 0;
91
92 // TODO[TL]: for anisotropic refinement we will probably need a table of
93 // submatrices with an array for each refine case
94 std::vector<FullMatrix<double>> face_embeddings(
95 1 << (dim - 1),
97 this->n_dofs_per_face(face_no)));
98 // TODO: Something goes wrong there. The error of the least squares fit
99 // is to large ...
100 // FETools::compute_face_embedding_matrices(*this, face_embeddings.data(), 0,
101 // 0);
102 this->interface_constraints.reinit((1 << (dim - 1)) *
103 this->n_dofs_per_face(face_no),
104 this->n_dofs_per_face(face_no));
105 unsigned int target_row = 0;
106 for (const auto &face_embedding : face_embeddings)
107 for (unsigned int i = 0; i < face_embedding.m(); ++i)
108 {
109 for (unsigned int j = 0; j < face_embedding.n(); ++j)
110 this->interface_constraints(target_row, j) = face_embedding(i, j);
111 ++target_row;
112 }
113
114 // We need to initialize the dof permutation table and the one for the sign
115 // change.
117}
118
119
120template <int dim>
121void
123{
124 // for 1d and 2d, do nothing
125 if (dim < 3)
126 return;
127
128 // TODO: Implement this for this class
129 return;
130}
131
132
133template <int dim>
134std::string
136{
137 // note that the
138 // FETools::get_fe_by_name
139 // function depends on the
140 // particular format of the string
141 // this function returns, so they
142 // have to be kept in synch
143
144 std::ostringstream namebuf;
145
146 namebuf << "FE_ABF<" << dim << ">(" << rt_order << ")";
147
148 return namebuf.str();
149}
150
151
152
153template <int dim>
154std::unique_ptr<FiniteElement<dim, dim>>
156{
157 return std::make_unique<FE_ABF<dim>>(rt_order);
158}
159
160
161//---------------------------------------------------------------------------
162// Auxiliary and internal functions
163//---------------------------------------------------------------------------
164
165
166
167// Version for 2d and higher. See above for 1d version
168template <int dim>
169void
171{
172 const QGauss<dim> cell_quadrature(deg + 2);
173 const unsigned int n_interior_points = cell_quadrature.size();
174
175 // TODO: the implementation makes the assumption that all faces have the
176 // same number of dofs
177 AssertDimension(this->n_unique_faces(), 1);
178 const unsigned int face_no = 0;
179
180 unsigned int n_face_points = (dim > 1) ? 1 : 0;
181 // compute (deg+1)^(dim-1)
182 for (unsigned int d = 1; d < dim; ++d)
183 n_face_points *= deg + 1;
184
185 this->generalized_support_points.resize(
186 GeometryInfo<dim>::faces_per_cell * n_face_points + n_interior_points);
187 this->generalized_face_support_points[face_no].resize(n_face_points);
188
189
190 // These might be required when the faces contribution is computed
191 // Therefore they will be initialized at this point.
192 std::array<std::unique_ptr<AnisotropicPolynomials<dim>>, dim> polynomials_abf;
193
194 // Generate x_1^{i} x_2^{r+1} ...
195 for (unsigned int dd = 0; dd < dim; ++dd)
196 {
197 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
198 for (unsigned int d = 0; d < dim; ++d)
199 poly[d].push_back(Polynomials::Monomial<double>(deg + 1));
201
202 polynomials_abf[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
203 }
204
205 // Number of the point being entered
206 unsigned int current = 0;
207
208 if (dim > 1)
209 {
210 const QGauss<dim - 1> face_points(deg + 1);
211 TensorProductPolynomials<dim - 1> legendre =
213
214 boundary_weights.reinit(n_face_points, legendre.n());
215
216 // Assert (face_points.size() == this->n_dofs_per_face(),
217 // ExcInternalError());
218
219 for (unsigned int k = 0; k < n_face_points; ++k)
220 {
221 this->generalized_face_support_points[face_no][k] =
222 face_points.point(k);
223 // Compute its quadrature
224 // contribution for each
225 // moment.
226 for (unsigned int i = 0; i < legendre.n(); ++i)
227 {
228 boundary_weights(k, i) =
229 face_points.weight(k) *
230 legendre.compute_value(i, face_points.point(k));
231 }
232 }
233
234 Quadrature<dim> faces =
235 QProjector<dim>::project_to_all_faces(this->reference_cell(),
236 face_points);
237 for (; current < GeometryInfo<dim>::faces_per_cell * n_face_points;
238 ++current)
239 {
240 // Enter the support point
241 // into the vector
242 this->generalized_support_points[current] = faces.point(current);
243 }
244
245
246 // Now initialize edge interior weights for the ABF elements.
247 // These are completely independent from the usual edge moments. They
248 // stem from applying the Gauss theorem to the nodal values, which
249 // was necessary to cast the ABF elements into the deal.II framework
250 // for vector valued elements.
251 boundary_weights_abf.reinit(faces.size(), polynomials_abf[0]->n() * dim);
252 for (unsigned int k = 0; k < faces.size(); ++k)
253 {
254 for (unsigned int i = 0; i < polynomials_abf[0]->n() * dim; ++i)
255 {
256 boundary_weights_abf(k, i) =
257 polynomials_abf[i % dim]->compute_value(i / dim,
258 faces.point(k)) *
259 faces.weight(k);
260 }
261 }
262 }
263
264 // Create Legendre basis for the
265 // space D_xi Q_k
266 if (deg > 0)
267 {
268 std::array<std::unique_ptr<AnisotropicPolynomials<dim>>, dim> polynomials;
269
270 for (unsigned int dd = 0; dd < dim; ++dd)
271 {
272 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
273 for (unsigned int d = 0; d < dim; ++d)
276
277 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
278 }
279
280 interior_weights.reinit(
281 TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
282
283 for (unsigned int k = 0; k < cell_quadrature.size(); ++k)
284 {
285 for (unsigned int i = 0; i < polynomials[0]->n(); ++i)
286 for (unsigned int d = 0; d < dim; ++d)
287 interior_weights(k, i, d) =
288 cell_quadrature.weight(k) *
289 polynomials[d]->compute_value(i, cell_quadrature.point(k));
290 }
291 }
292
293
294 // Decouple the creation of the generalized support points
295 // from computation of interior weights.
296 for (unsigned int k = 0; k < cell_quadrature.size(); ++k)
297 this->generalized_support_points[current++] = cell_quadrature.point(k);
298
299 // Additional functionality for the ABF elements
300 // TODO: Here the canonical extension of the principle
301 // behind the ABF elements is implemented. It is unclear,
302 // if this really leads to the ABF spaces in 3d!
303 interior_weights_abf.reinit(TableIndices<3>(cell_quadrature.size(),
304 polynomials_abf[0]->n() * dim,
305 dim));
306 Tensor<1, dim> poly_grad;
307
308 for (unsigned int k = 0; k < cell_quadrature.size(); ++k)
309 {
310 for (unsigned int i = 0; i < polynomials_abf[0]->n() * dim; ++i)
311 {
312 poly_grad =
313 polynomials_abf[i % dim]->compute_grad(i / dim,
314 cell_quadrature.point(k)) *
315 cell_quadrature.weight(k);
316 // The minus sign comes from the use of the Gauss theorem to replace
317 // the divergence.
318 for (unsigned int d = 0; d < dim; ++d)
319 interior_weights_abf(k, i, d) = -poly_grad[d];
320 }
321 }
322
323 Assert(current == this->generalized_support_points.size(),
325}
326
327
328
329// This function is the same Raviart-Thomas interpolation performed by
330// interpolate. Still, we cannot use interpolate, since it was written
331// for smooth functions. The functions interpolated here are not
332// smooth, maybe even not continuous. Therefore, we must double the
333// number of quadrature points in each direction in order to integrate
334// only smooth functions.
335
336// Then again, the interpolated function is chosen such that the
337// moments coincide with the function to be interpolated.
338
339template <int dim>
340void
342{
343 if (dim == 1)
344 {
345 unsigned int iso = RefinementCase<dim>::isotropic_refinement - 1;
346 for (unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_cell;
347 ++i)
348 this->restriction[iso][i].reinit(0, 0);
349 return;
350 }
351 unsigned int iso = RefinementCase<dim>::isotropic_refinement - 1;
352 const QGauss<dim - 1> q_base(rt_order + 1);
353 const unsigned int n_face_points = q_base.size();
354 // First, compute interpolation on
355 // subfaces
356 for (const unsigned int face : GeometryInfo<dim>::face_indices())
357 {
358 // The shape functions of the
359 // child cell are evaluated
360 // in the quadrature points
361 // of a full face.
362 Quadrature<dim> q_face =
363 QProjector<dim>::project_to_face(this->reference_cell(), q_base, face);
364 // Store shape values, since the
365 // evaluation suffers if not
366 // ordered by point
367 Table<2, double> cached_values_face(this->n_dofs_per_cell(),
368 q_face.size());
369 for (unsigned int k = 0; k < q_face.size(); ++k)
370 for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
371 cached_values_face(i, k) = this->shape_value_component(
373
374 for (unsigned int sub = 0; sub < GeometryInfo<dim>::max_children_per_face;
375 ++sub)
376 {
377 // The weight functions for
378 // the coarse face are
379 // evaluated on the subface
380 // only.
382 this->reference_cell(), q_base, face, sub);
383 const unsigned int child = GeometryInfo<dim>::child_cell_on_face(
385
386 // On a certain face, we must
387 // compute the moments of ALL
388 // fine level functions with
389 // the coarse level weight
390 // functions belonging to
391 // that face. Due to the
392 // orthogonalization process
393 // when building the shape
394 // functions, these weights
395 // are equal to the
396 // corresponding shape
397 // functions.
398 for (unsigned int k = 0; k < n_face_points; ++k)
399 for (unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
400 ++i_child)
401 for (unsigned int i_face = 0;
402 i_face < this->n_dofs_per_face(face);
403 ++i_face)
404 {
405 // The quadrature
406 // weights on the
407 // subcell are NOT
408 // transformed, so we
409 // have to do it here.
410 this->restriction[iso][child](
411 face * this->n_dofs_per_face(face) + i_face, i_child) +=
413 cached_values_face(i_child, k) *
414 this->shape_value_component(
415 face * this->n_dofs_per_face(face) + i_face,
416 q_sub.point(k),
418 }
419 }
420 }
421
422 if (rt_order == 0)
423 return;
424
425 // Create Legendre basis for the
426 // space D_xi Q_k. Here, we cannot
427 // use the shape functions
428 std::array<std::unique_ptr<AnisotropicPolynomials<dim>>, dim> polynomials;
429 for (unsigned int dd = 0; dd < dim; ++dd)
430 {
431 std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
432 for (unsigned int d = 0; d < dim; ++d)
434 poly[dd] = Polynomials::Legendre::generate_complete_basis(rt_order - 1);
435
436 polynomials[dd] = std::make_unique<AnisotropicPolynomials<dim>>(poly);
437 }
438
439 // TODO: the implementation makes the assumption that all faces have the
440 // same number of dofs
441 AssertDimension(this->n_unique_faces(), 1);
442 const unsigned int face_no = 0;
443
444 const QGauss<dim> q_cell(rt_order + 1);
445 const unsigned int start_cell_dofs =
446 GeometryInfo<dim>::faces_per_cell * this->n_dofs_per_face(face_no);
447
448 // Store shape values, since the
449 // evaluation suffers if not
450 // ordered by point
451 Table<3, double> cached_values_cell(this->n_dofs_per_cell(),
452 q_cell.size(),
453 dim);
454 for (unsigned int k = 0; k < q_cell.size(); ++k)
455 for (unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
456 for (unsigned int d = 0; d < dim; ++d)
457 cached_values_cell(i, k, d) =
458 this->shape_value_component(i, q_cell.point(k), d);
459
460 for (unsigned int child = 0; child < GeometryInfo<dim>::max_children_per_cell;
461 ++child)
462 {
463 Quadrature<dim> q_sub =
464 QProjector<dim>::project_to_child(this->reference_cell(),
465 q_cell,
466 child);
467
468 for (unsigned int k = 0; k < q_sub.size(); ++k)
469 for (unsigned int i_child = 0; i_child < this->n_dofs_per_cell();
470 ++i_child)
471 for (unsigned int d = 0; d < dim; ++d)
472 for (unsigned int i_weight = 0; i_weight < polynomials[d]->n();
473 ++i_weight)
474 {
475 this->restriction[iso][child](start_cell_dofs + i_weight * dim +
476 d,
477 i_child) +=
478 q_sub.weight(k) * cached_values_cell(i_child, k, d) *
479 polynomials[d]->compute_value(i_weight, q_sub.point(k));
480 }
481 }
482}
483
484
485
486template <int dim>
487std::vector<unsigned int>
488FE_ABF<dim>::get_dpo_vector(const unsigned int rt_order)
489{
490 if (dim == 1)
491 {
492 Assert(false, ExcImpossibleInDim(1));
493 return std::vector<unsigned int>();
494 }
495
496 // the element is face-based (not
497 // to be confused with George
498 // W. Bush's Faith Based
499 // Initiative...), and we have
500 // (rt_order+1)^(dim-1) DoFs per face
501 unsigned int dofs_per_face = 1;
502 for (unsigned int d = 0; d < dim - 1; ++d)
503 dofs_per_face *= rt_order + 1;
504
505 // and then there are interior dofs
506 const unsigned int interior_dofs = dim * (rt_order + 1) * dofs_per_face;
507
508 std::vector<unsigned int> dpo(dim + 1);
509 dpo[dim - 1] = dofs_per_face;
510 dpo[dim] = interior_dofs;
511
512 return dpo;
513}
514
515//---------------------------------------------------------------------------
516// Data field initialization
517//---------------------------------------------------------------------------
518
519template <int dim>
520bool
521FE_ABF<dim>::has_support_on_face(const unsigned int shape_index,
522 const unsigned int face_index) const
523{
524 AssertIndexRange(shape_index, this->n_dofs_per_cell());
526
527 // Return computed values if we
528 // know them easily. Otherwise, it
529 // is always safe to return true.
530 switch (rt_order)
531 {
532 case 0:
533 {
534 switch (dim)
535 {
536 case 2:
537 {
538 // only on the one
539 // non-adjacent face
540 // are the values
541 // actually zero. list
542 // these in a table
543 return (face_index !=
545 }
546
547 default:
548 return true;
549 }
550 }
551
552 default: // other rt_order
553 return true;
554 }
555
556 return true;
557}
558
559
560
561template <int dim>
562void
564 const std::vector<Vector<double>> &support_point_values,
565 std::vector<double> &nodal_values) const
566{
567 Assert(support_point_values.size() == this->generalized_support_points.size(),
568 ExcDimensionMismatch(support_point_values.size(),
569 this->generalized_support_points.size()));
570 Assert(support_point_values[0].size() == this->n_components(),
571 ExcDimensionMismatch(support_point_values[0].size(),
572 this->n_components()));
573 Assert(nodal_values.size() == this->n_dofs_per_cell(),
574 ExcDimensionMismatch(nodal_values.size(), this->n_dofs_per_cell()));
575
576 std::fill(nodal_values.begin(), nodal_values.end(), 0.);
577
578 const unsigned int n_face_points = boundary_weights.size(0);
579 for (const unsigned int face : GeometryInfo<dim>::face_indices())
580 for (unsigned int k = 0; k < n_face_points; ++k)
581 for (unsigned int i = 0; i < boundary_weights.size(1); ++i)
582 {
583 nodal_values[i + face * this->n_dofs_per_face(face)] +=
584 boundary_weights(k, i) *
585 support_point_values[face * n_face_points + k][GeometryInfo<
586 dim>::unit_normal_direction[face]];
587 }
588
589 // TODO: the implementation makes the assumption that all faces have the
590 // same number of dofs
591 AssertDimension(this->n_unique_faces(), 1);
592 const unsigned int face_no = 0;
593
594 const unsigned int start_cell_dofs =
595 GeometryInfo<dim>::faces_per_cell * this->n_dofs_per_face(face_no);
596 const unsigned int start_cell_points =
597 GeometryInfo<dim>::faces_per_cell * n_face_points;
598
599 for (unsigned int k = 0; k < interior_weights.size(0); ++k)
600 for (unsigned int i = 0; i < interior_weights.size(1); ++i)
601 for (unsigned int d = 0; d < dim; ++d)
602 nodal_values[start_cell_dofs + i * dim + d] +=
603 interior_weights(k, i, d) *
604 support_point_values[k + start_cell_points][d];
605
606 const unsigned int start_abf_dofs =
607 start_cell_dofs + interior_weights.size(1) * dim;
608
609 // Cell integral of ABF terms
610 for (unsigned int k = 0; k < interior_weights_abf.size(0); ++k)
611 for (unsigned int i = 0; i < interior_weights_abf.size(1); ++i)
612 for (unsigned int d = 0; d < dim; ++d)
613 nodal_values[start_abf_dofs + i] +=
614 interior_weights_abf(k, i, d) *
615 support_point_values[k + start_cell_points][d];
616
617 // Face integral of ABF terms
618 for (const unsigned int face : GeometryInfo<dim>::face_indices())
619 {
620 const double n_orient = GeometryInfo<dim>::unit_normal_orientation[face];
621 for (unsigned int fp = 0; fp < n_face_points; ++fp)
622 {
623 // TODO: Check what the face_orientation, face_flip and face_rotation
624 // have to be in 3d
626 this->reference_cell(),
627 face,
629 n_face_points);
630 for (unsigned int i = 0; i < boundary_weights_abf.size(1); ++i)
631 nodal_values[start_abf_dofs + i] +=
632 n_orient * boundary_weights_abf(k + fp, i) *
633 support_point_values[face * n_face_points + fp][GeometryInfo<
634 dim>::unit_normal_direction[face]];
635 }
636 }
637
638 // TODO: Check if this "correction" can be removed.
639 for (unsigned int i = 0; i < boundary_weights_abf.size(1); ++i)
640 if (std::fabs(nodal_values[start_abf_dofs + i]) < 1.0e-16)
641 nodal_values[start_abf_dofs + i] = 0.0;
642}
643
644
645
646template <int dim>
647std::size_t
649{
651 return 0;
652}
653
654
655
656/*-------------- Explicit Instantiations -------------------------------*/
657#include "fe_abf.inst"
658
void initialize_quad_dof_index_permutation_and_sign_change()
Definition fe_abf.cc:122
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
Definition fe_abf.cc:155
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const override
Definition fe_abf.cc:563
virtual std::size_t memory_consumption() const override
Definition fe_abf.cc:648
void initialize_restriction()
Definition fe_abf.cc:341
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
Definition fe_abf.cc:521
virtual std::string get_name() const override
Definition fe_abf.cc:135
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition fe_abf.cc:488
void initialize_support_points(const unsigned int rt_degree)
Definition fe_abf.cc:170
friend class FE_ABF
Definition fe_abf.h:255
FullMatrix< double > inverse_node_matrix
std::vector< MappingKind > mapping_kind
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_unique_faces() const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
FullMatrix< double > interface_constraints
Definition fe.h:2439
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition fe.h:2427
void invert(const FullMatrix< number2 > &M)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
static std::vector< Polynomial< number > > generate_complete_basis(const unsigned int degree)
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
static Quadrature< dim > project_to_child(const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature, const unsigned int child_no)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
static void project_to_subface(const ReferenceCell &reference_cell, 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 ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
const Point< dim > & point(const unsigned int i) const
double weight(const unsigned int i) const
unsigned int size() const
static constexpr unsigned char default_combined_face_orientation()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
@ mapping_raviart_thomas
Definition mapping.h:134
#define DEAL_II_NOT_IMPLEMENTED()
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
FullMatrix< double > compute_node_matrix(const FiniteElement< dim, spacedim > &fe)
constexpr T fixed_power(const T t)
Definition utilities.h:942
STL namespace.
static unsigned int child_cell_on_face(const RefinementCase< dim > &ref_case, const unsigned int face, const unsigned int subface, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false, const RefinementCase< dim - 1 > &face_refinement_case=RefinementCase< dim - 1 >::isotropic_refinement)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()