Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_raviart_thomas.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/qprojector.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/quadrature_lib.h>
20 #include <deal.II/base/std_cxx14/memory.h>
21 #include <deal.II/base/table.h>
22 #include <deal.II/base/utilities.h>
23 
24 #include <deal.II/dofs/dof_accessor.h>
25 
26 #include <deal.II/fe/fe.h>
27 #include <deal.II/fe/fe_raviart_thomas.h>
28 #include <deal.II/fe/fe_tools.h>
29 #include <deal.II/fe/fe_values.h>
30 #include <deal.II/fe/mapping.h>
31 
32 #include <deal.II/grid/tria.h>
33 #include <deal.II/grid/tria_iterator.h>
34 
35 #include <iostream>
36 #include <sstream>
37 
38 // TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
39 // adjust_line_dof_index_for_line_orientation_table fields, and write tests
40 // similar to bits/face_orientation_and_fe_q_*
41 
42 
43 DEAL_II_NAMESPACE_OPEN
44 
45 
46 template <int dim>
49  deg,
50  FiniteElementData<dim>(get_dpo_vector(deg),
51  dim,
52  deg + 1,
53  FiniteElementData<dim>::Hdiv),
54  std::vector<bool>(PolynomialsRaviartThomas<dim>::compute_n_pols(deg),
55  true),
56  std::vector<ComponentMask>(PolynomialsRaviartThomas<dim>::compute_n_pols(
57  deg),
58  std::vector<bool>(dim, true)))
59 {
60  Assert(dim >= 2, ExcImpossibleInDim(dim));
61  const unsigned int n_dofs = this->dofs_per_cell;
62 
64  // First, initialize the
65  // generalized support points and
66  // quadrature weights, since they
67  // are required for interpolation.
69 
70  // Now compute the inverse node matrix, generating the correct
71  // basis functions from the raw ones. For a discussion of what
72  // exactly happens here, see FETools::compute_node_matrix.
74  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
75  this->inverse_node_matrix.invert(M);
76  // From now on, the shape functions provided by FiniteElement::shape_value
77  // and similar functions will be the correct ones, not
78  // the raw shape functions from the polynomial space anymore.
79 
80  // Reinit the vectors of
81  // restriction and prolongation
82  // matrices to the right sizes.
83  // Restriction only for isotropic
84  // refinement
86  // Fill prolongation matrices with embedding operators
89 
90  // TODO[TL]: for anisotropic refinement we will probably need a table of
91  // submatrices with an array for each refine case
93  for (unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
94  face_embeddings[i].reinit(this->dofs_per_face, this->dofs_per_face);
95  FETools::compute_face_embedding_matrices<dim, double>(*this,
96  face_embeddings,
97  0,
98  0);
99  this->interface_constraints.reinit((1 << (dim - 1)) * this->dofs_per_face,
100  this->dofs_per_face);
101  unsigned int target_row = 0;
102  for (unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
103  for (unsigned int i = 0; i < face_embeddings[d].m(); ++i)
104  {
105  for (unsigned int j = 0; j < face_embeddings[d].n(); ++j)
106  this->interface_constraints(target_row, j) = face_embeddings[d](i, j);
107  ++target_row;
108  }
109 }
110 
111 
112 
113 template <int dim>
114 std::string
116 {
117  // note that the
118  // FETools::get_fe_by_name
119  // function depends on the
120  // particular format of the string
121  // this function returns, so they
122  // have to be kept in synch
123 
124  // note that this->degree is the maximal
125  // polynomial degree and is thus one higher
126  // than the argument given to the
127  // constructor
128  std::ostringstream namebuf;
129  namebuf << "FE_RaviartThomas<" << dim << ">(" << this->degree - 1 << ")";
130 
131  return namebuf.str();
132 }
133 
134 
135 template <int dim>
136 std::unique_ptr<FiniteElement<dim, dim>>
138 {
139  return std_cxx14::make_unique<FE_RaviartThomas<dim>>(*this);
140 }
141 
142 
143 //---------------------------------------------------------------------------
144 // Auxiliary and internal functions
145 //---------------------------------------------------------------------------
146 
147 
148 template <int dim>
149 void
151 {
152  QGauss<dim> cell_quadrature(deg + 1);
153  const unsigned int n_interior_points = (deg > 0) ? cell_quadrature.size() : 0;
154 
155  unsigned int n_face_points = (dim > 1) ? 1 : 0;
156  // compute (deg+1)^(dim-1)
157  for (unsigned int d = 1; d < dim; ++d)
158  n_face_points *= deg + 1;
159 
160 
161  this->generalized_support_points.resize(
162  GeometryInfo<dim>::faces_per_cell * n_face_points + n_interior_points);
163  this->generalized_face_support_points.resize(n_face_points);
164 
165  // Number of the point being entered
166  unsigned int current = 0;
167 
168  if (dim > 1)
169  {
170  QGauss<dim - 1> face_points(deg + 1);
171  TensorProductPolynomials<dim - 1> legendre =
173 
174  boundary_weights.reinit(n_face_points, legendre.n());
175 
176  // Assert (face_points.size() == this->dofs_per_face,
177  // ExcInternalError());
178 
179  for (unsigned int k = 0; k < n_face_points; ++k)
180  {
181  this->generalized_face_support_points[k] = face_points.point(k);
182  // Compute its quadrature
183  // contribution for each
184  // moment.
185  for (unsigned int i = 0; i < legendre.n(); ++i)
186  {
187  boundary_weights(k, i) =
188  face_points.weight(k) *
189  legendre.compute_value(i, face_points.point(k));
190  }
191  }
192 
193  Quadrature<dim> faces =
195  for (; current < GeometryInfo<dim>::faces_per_cell * n_face_points;
196  ++current)
197  {
198  // Enter the support point
199  // into the vector
200  this->generalized_support_points[current] =
202  0, true, false, false, n_face_points));
203  }
204  }
205 
206  if (deg == 0)
207  return;
208 
209  // Create Legendre basis for the space D_xi Q_k
210  std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
211  for (unsigned int dd = 0; dd < dim; ++dd)
212  {
213  std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
214  for (unsigned int d = 0; d < dim; ++d)
217 
218  polynomials[dd] =
219  std_cxx14::make_unique<AnisotropicPolynomials<dim>>(poly);
220  }
221 
222  interior_weights.reinit(
223  TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
224 
225  for (unsigned int k = 0; k < cell_quadrature.size(); ++k)
226  {
227  this->generalized_support_points[current++] = cell_quadrature.point(k);
228  for (unsigned int i = 0; i < polynomials[0]->n(); ++i)
229  for (unsigned int d = 0; d < dim; ++d)
230  interior_weights(k, i, d) =
231  cell_quadrature.weight(k) *
232  polynomials[d]->compute_value(i, cell_quadrature.point(k));
233  }
234 
235  Assert(current == this->generalized_support_points.size(),
236  ExcInternalError());
237 }
238 
239 
240 
241 template <>
242 void
244 {
245  // there is only one refinement case in 1d,
246  // which is the isotropic one (first index of
247  // the matrix array has to be 0)
248  for (unsigned int i = 0; i < GeometryInfo<1>::max_children_per_cell; ++i)
249  this->restriction[0][i].reinit(0, 0);
250 }
251 
252 
253 
254 // This function is the same Raviart-Thomas interpolation performed by
255 // interpolate. Still, we cannot use interpolate, since it was written
256 // for smooth functions. The functions interpolated here are not
257 // smooth, maybe even not continuous. Therefore, we must double the
258 // number of quadrature points in each direction in order to integrate
259 // only smooth functions.
260 
261 // Then again, the interpolated function is chosen such that the
262 // moments coincide with the function to be interpolated.
263 
264 template <int dim>
265 void
267 {
268  const unsigned int iso = RefinementCase<dim>::isotropic_refinement - 1;
269 
270  QGauss<dim - 1> q_base(this->degree);
271  const unsigned int n_face_points = q_base.size();
272  // First, compute interpolation on
273  // subfaces
274  for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
275  {
276  // The shape functions of the
277  // child cell are evaluated
278  // in the quadrature points
279  // of a full face.
280  Quadrature<dim> q_face = QProjector<dim>::project_to_face(q_base, face);
281  // Store shape values, since the
282  // evaluation suffers if not
283  // ordered by point
284  Table<2, double> cached_values_on_face(this->dofs_per_cell,
285  q_face.size());
286  for (unsigned int k = 0; k < q_face.size(); ++k)
287  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
288  cached_values_on_face(i, k) = this->shape_value_component(
289  i, q_face.point(k), GeometryInfo<dim>::unit_normal_direction[face]);
290 
291  for (unsigned int sub = 0; sub < GeometryInfo<dim>::max_children_per_face;
292  ++sub)
293  {
294  // The weight functions for
295  // the coarse face are
296  // evaluated on the subface
297  // only.
298  Quadrature<dim> q_sub =
299  QProjector<dim>::project_to_subface(q_base, face, sub);
300  const unsigned int child = GeometryInfo<dim>::child_cell_on_face(
302 
303  // On a certain face, we must
304  // compute the moments of ALL
305  // fine level functions with
306  // the coarse level weight
307  // functions belonging to
308  // that face. Due to the
309  // orthogonalization process
310  // when building the shape
311  // functions, these weights
312  // are equal to the
313  // corresponding shape
314  // functions.
315  for (unsigned int k = 0; k < n_face_points; ++k)
316  for (unsigned int i_child = 0; i_child < this->dofs_per_cell;
317  ++i_child)
318  for (unsigned int i_face = 0; i_face < this->dofs_per_face;
319  ++i_face)
320  {
321  // The quadrature
322  // weights on the
323  // subcell are NOT
324  // transformed, so we
325  // have to do it here.
326  this->restriction[iso][child](face * this->dofs_per_face +
327  i_face,
328  i_child) +=
329  Utilities::fixed_power<dim - 1>(.5) * q_sub.weight(k) *
330  cached_values_on_face(i_child, k) *
331  this->shape_value_component(
332  face * this->dofs_per_face + i_face,
333  q_sub.point(k),
335  }
336  }
337  }
338 
339  if (this->degree == 1)
340  return;
341 
342  // Create Legendre basis for the space D_xi Q_k. Here, we cannot
343  // use the shape functions
344  std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
345  for (unsigned int dd = 0; dd < dim; ++dd)
346  {
347  std::vector<std::vector<Polynomials::Polynomial<double>>> poly(dim);
348  for (unsigned int d = 0; d < dim; ++d)
349  poly[d] =
351  poly[dd] =
353 
354  polynomials[dd] =
355  std_cxx14::make_unique<AnisotropicPolynomials<dim>>(poly);
356  }
357 
358  QGauss<dim> q_cell(this->degree);
359  const unsigned int start_cell_dofs =
360  GeometryInfo<dim>::faces_per_cell * this->dofs_per_face;
361 
362  // Store shape values, since the
363  // evaluation suffers if not
364  // ordered by point
365  Table<3, double> cached_values_on_cell(this->dofs_per_cell,
366  q_cell.size(),
367  dim);
368  for (unsigned int k = 0; k < q_cell.size(); ++k)
369  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
370  for (unsigned int d = 0; d < dim; ++d)
371  cached_values_on_cell(i, k, d) =
372  this->shape_value_component(i, q_cell.point(k), d);
373 
374  for (unsigned int child = 0; child < GeometryInfo<dim>::max_children_per_cell;
375  ++child)
376  {
377  Quadrature<dim> q_sub = QProjector<dim>::project_to_child(q_cell, child);
378 
379  for (unsigned int k = 0; k < q_sub.size(); ++k)
380  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
381  for (unsigned int d = 0; d < dim; ++d)
382  for (unsigned int i_weight = 0; i_weight < polynomials[d]->n();
383  ++i_weight)
384  {
385  this->restriction[iso][child](start_cell_dofs + i_weight * dim +
386  d,
387  i_child) +=
388  q_sub.weight(k) * cached_values_on_cell(i_child, k, d) *
389  polynomials[d]->compute_value(i_weight, q_sub.point(k));
390  }
391  }
392 }
393 
394 
395 
396 template <int dim>
397 std::vector<unsigned int>
399 {
400  // the element is face-based and we have
401  // (deg+1)^(dim-1) DoFs per face
402  unsigned int dofs_per_face = 1;
403  for (unsigned int d = 1; d < dim; ++d)
404  dofs_per_face *= deg + 1;
405 
406  // and then there are interior dofs
407  const unsigned int interior_dofs = dim * deg * dofs_per_face;
408 
409  std::vector<unsigned int> dpo(dim + 1);
410  dpo[dim - 1] = dofs_per_face;
411  dpo[dim] = interior_dofs;
412 
413  return dpo;
414 }
415 
416 
417 
418 template <int dim>
419 std::pair<Table<2, bool>, std::vector<unsigned int>>
421 {
422  Table<2, bool> constant_modes(dim, this->dofs_per_cell);
423  for (unsigned int d = 0; d < dim; ++d)
424  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
425  constant_modes(d, i) = true;
426  std::vector<unsigned int> components;
427  for (unsigned int d = 0; d < dim; ++d)
428  components.push_back(d);
429  return std::pair<Table<2, bool>, std::vector<unsigned int>>(constant_modes,
430  components);
431 }
432 
433 
434 
435 //---------------------------------------------------------------------------
436 // Data field initialization
437 //---------------------------------------------------------------------------
438 
439 
440 template <int dim>
441 bool
442 FE_RaviartThomas<dim>::has_support_on_face(const unsigned int shape_index,
443  const unsigned int face_index) const
444 {
445  Assert(shape_index < this->dofs_per_cell,
446  ExcIndexRange(shape_index, 0, this->dofs_per_cell));
449 
450  // Return computed values if we
451  // know them easily. Otherwise, it
452  // is always safe to return true.
453  switch (this->degree)
454  {
455  case 1:
456  {
457  switch (dim)
458  {
459  case 2:
460  {
461  // only on the one
462  // non-adjacent face
463  // are the values
464  // actually zero. list
465  // these in a table
466  return (face_index !=
467  GeometryInfo<dim>::opposite_face[shape_index]);
468  }
469 
470  default:
471  return true;
472  }
473  }
474 
475  default: // other rt_order
476  return true;
477  }
478 
479  return true;
480 }
481 
482 
483 
484 template <int dim>
485 void
487  const std::vector<Vector<double>> &support_point_values,
488  std::vector<double> & nodal_values) const
489 {
490  Assert(support_point_values.size() == this->generalized_support_points.size(),
491  ExcDimensionMismatch(support_point_values.size(),
492  this->generalized_support_points.size()));
493  Assert(nodal_values.size() == this->dofs_per_cell,
494  ExcDimensionMismatch(nodal_values.size(), this->dofs_per_cell));
495  Assert(support_point_values[0].size() == this->n_components(),
496  ExcDimensionMismatch(support_point_values[0].size(),
497  this->n_components()));
498 
499  std::fill(nodal_values.begin(), nodal_values.end(), 0.);
500 
501  const unsigned int n_face_points = boundary_weights.size(0);
502  for (unsigned int face = 0; face < GeometryInfo<dim>::faces_per_cell; ++face)
503  for (unsigned int k = 0; k < n_face_points; ++k)
504  for (unsigned int i = 0; i < boundary_weights.size(1); ++i)
505  {
506  nodal_values[i + face * this->dofs_per_face] +=
507  boundary_weights(k, i) *
508  support_point_values[face * n_face_points + k](
510  }
511 
512  const unsigned int start_cell_dofs =
513  GeometryInfo<dim>::faces_per_cell * this->dofs_per_face;
514  const unsigned int start_cell_points =
515  GeometryInfo<dim>::faces_per_cell * n_face_points;
516 
517  for (unsigned int k = 0; k < interior_weights.size(0); ++k)
518  for (unsigned int i = 0; i < interior_weights.size(1); ++i)
519  for (unsigned int d = 0; d < dim; ++d)
520  nodal_values[start_cell_dofs + i * dim + d] +=
521  interior_weights(k, i, d) *
522  support_point_values[k + start_cell_points](d);
523 }
524 
525 
526 
527 template <int dim>
528 std::size_t
530 {
531  Assert(false, ExcNotImplemented());
532  return 0;
533 }
534 
535 
536 
537 // explicit instantiations
538 #include "fe_raviart_thomas.inst"
539 
540 
541 DEAL_II_NAMESPACE_CLOSE
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1063
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 > interface_constraints
Definition: fe.h:2468
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 std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:866
const Point< dim > & point(const unsigned int i) const
STL namespace.
FullMatrix< double > compute_node_matrix(const FiniteElement< dim, spacedim > &fe)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual std::size_t memory_consumption() const override
virtual std::string get_name() const override
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
friend class FE_RaviartThomas
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2456
void invert(const FullMatrix< number2 > &M)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:1407
void initialize_support_points(const unsigned int rt_degree)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:282
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim >> &q_points)
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
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
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)
const unsigned int dofs_per_face
Definition: fe_base.h:275
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:276
static ::ExceptionBase & ExcNotImplemented()
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const override
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()