Reference documentation for deal.II version Git 7026f387cc 2019-10-15 14:19:01 -0400
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_raviart_thomas_nodal.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2018 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_lib.h>
19 #include <deal.II/base/std_cxx14/memory.h>
20 
21 #include <deal.II/dofs/dof_accessor.h>
22 
23 #include <deal.II/fe/fe.h>
24 #include <deal.II/fe/fe_nothing.h>
25 #include <deal.II/fe/fe_raviart_thomas.h>
26 #include <deal.II/fe/fe_tools.h>
27 #include <deal.II/fe/fe_values.h>
28 #include <deal.II/fe/mapping.h>
29 
30 #include <deal.II/grid/tria.h>
31 #include <deal.II/grid/tria_iterator.h>
32 
33 #include <sstream>
34 
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 
39 template <int dim>
41  : FE_PolyTensor<dim>(PolynomialsRaviartThomas<dim>(deg),
42  FiniteElementData<dim>(get_dpo_vector(deg),
43  dim,
44  deg + 1,
45  FiniteElementData<dim>::Hdiv),
46  get_ria_vector(deg),
47  std::vector<ComponentMask>(
48  PolynomialsRaviartThomas<dim>::n_polynomials(deg),
49  std::vector<bool>(dim, true)))
50 {
51  Assert(dim >= 2, ExcImpossibleInDim(dim));
52  const unsigned int n_dofs = this->dofs_per_cell;
53 
55  // First, initialize the
56  // generalized support points and
57  // quadrature weights, since they
58  // are required for interpolation.
60 
61  // Now compute the inverse node matrix, generating the correct
62  // basis functions from the raw ones. For a discussion of what
63  // exactly happens here, see FETools::compute_node_matrix.
65  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
66  this->inverse_node_matrix.invert(M);
67  // From now on, the shape functions provided by FiniteElement::shape_value
68  // and similar functions will be the correct ones, not
69  // the raw shape functions from the polynomial space anymore.
70 
71  // Reinit the vectors of
72  // prolongation matrices to the
73  // right sizes. There are no
74  // restriction matrices implemented
75  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
76  ref_case < RefinementCase<dim>::isotropic_refinement + 1;
77  ++ref_case)
78  {
79  const unsigned int nc =
81 
82  for (unsigned int i = 0; i < nc; ++i)
83  this->prolongation[ref_case - 1][i].reinit(n_dofs, n_dofs);
84  }
85  // Fill prolongation matrices with embedding operators
87  // TODO[TL]: for anisotropic refinement we will probably need a table of
88  // submatrices with an array for each refine case
90  for (unsigned int i = 0; i < GeometryInfo<dim>::max_children_per_face; ++i)
91  face_embeddings[i].reinit(this->dofs_per_face, this->dofs_per_face);
92  FETools::compute_face_embedding_matrices<dim, double>(*this,
93  face_embeddings,
94  0,
95  0);
96  this->interface_constraints.reinit((1 << (dim - 1)) * this->dofs_per_face,
97  this->dofs_per_face);
98  unsigned int target_row = 0;
99  for (unsigned int d = 0; d < GeometryInfo<dim>::max_children_per_face; ++d)
100  for (unsigned int i = 0; i < face_embeddings[d].m(); ++i)
101  {
102  for (unsigned int j = 0; j < face_embeddings[d].n(); ++j)
103  this->interface_constraints(target_row, j) = face_embeddings[d](i, j);
104  ++target_row;
105  }
106 }
107 
108 
109 
110 template <int dim>
111 std::string
113 {
114  // note that the
115  // FETools::get_fe_by_name
116  // function depends on the
117  // particular format of the string
118  // this function returns, so they
119  // have to be kept in synch
120 
121  // note that this->degree is the maximal
122  // polynomial degree and is thus one higher
123  // than the argument given to the
124  // constructor
125  std::ostringstream namebuf;
126  namebuf << "FE_RaviartThomasNodal<" << dim << ">(" << this->degree - 1 << ")";
127 
128  return namebuf.str();
129 }
130 
131 
132 template <int dim>
133 std::unique_ptr<FiniteElement<dim, dim>>
135 {
136  return std_cxx14::make_unique<FE_RaviartThomasNodal<dim>>(*this);
137 }
138 
139 
140 //---------------------------------------------------------------------------
141 // Auxiliary and internal functions
142 //---------------------------------------------------------------------------
143 
144 
145 
146 template <int dim>
147 void
149 {
150  this->generalized_support_points.resize(this->dofs_per_cell);
152 
153  // Number of the point being entered
154  unsigned int current = 0;
155 
156  // On the faces, we choose as many
157  // Gauss points as necessary to
158  // determine the normal component
159  // uniquely. This is the deg of
160  // the Raviart-Thomas element plus
161  // one.
162  if (dim > 1)
163  {
164  QGauss<dim - 1> face_points(deg + 1);
165  Assert(face_points.size() == this->dofs_per_face, ExcInternalError());
166  for (unsigned int k = 0; k < this->dofs_per_face; ++k)
167  this->generalized_face_support_points[k] = face_points.point(k);
168  Quadrature<dim> faces =
170  for (unsigned int k = 0;
171  k < this->dofs_per_face * GeometryInfo<dim>::faces_per_cell;
172  ++k)
173  this->generalized_support_points[k] =
175  0, true, false, false, this->dofs_per_face));
176 
177  current = this->dofs_per_face * GeometryInfo<dim>::faces_per_cell;
178  }
179 
180  if (deg == 0)
181  return;
182  // In the interior, we need
183  // anisotropic Gauss quadratures,
184  // different for each direction.
185  QGauss<1> high(deg + 1);
186  QGauss<1> low(deg);
187 
188  for (unsigned int d = 0; d < dim; ++d)
189  {
190  std::unique_ptr<QAnisotropic<dim>> quadrature;
191  switch (dim)
192  {
193  case 1:
194  quadrature = std_cxx14::make_unique<QAnisotropic<dim>>(high);
195  break;
196  case 2:
197  quadrature = std_cxx14::make_unique<QAnisotropic<dim>>(
198  ((d == 0) ? low : high), ((d == 1) ? low : high));
199  break;
200  case 3:
201  quadrature =
202  std_cxx14::make_unique<QAnisotropic<dim>>(((d == 0) ? low : high),
203  ((d == 1) ? low : high),
204  ((d == 2) ? low :
205  high));
206  break;
207  default:
208  Assert(false, ExcNotImplemented());
209  }
210 
211  for (unsigned int k = 0; k < quadrature->size(); ++k)
212  this->generalized_support_points[current++] = quadrature->point(k);
213  }
214  Assert(current == this->dofs_per_cell, ExcInternalError());
215 }
216 
217 
218 
219 template <int dim>
220 std::vector<unsigned int>
222 {
223  // the element is face-based and we have
224  // (deg+1)^(dim-1) DoFs per face
225  unsigned int dofs_per_face = 1;
226  for (unsigned int d = 1; d < dim; ++d)
227  dofs_per_face *= deg + 1;
228 
229  // and then there are interior dofs
230  const unsigned int interior_dofs = dim * deg * dofs_per_face;
231 
232  std::vector<unsigned int> dpo(dim + 1);
233  dpo[dim - 1] = dofs_per_face;
234  dpo[dim] = interior_dofs;
235 
236  return dpo;
237 }
238 
239 
240 
241 template <>
242 std::vector<bool>
244 {
245  Assert(false, ExcImpossibleInDim(1));
246  return std::vector<bool>();
247 }
248 
249 
250 
251 template <int dim>
252 std::vector<bool>
254 {
255  const unsigned int dofs_per_cell =
257  unsigned int dofs_per_face = deg + 1;
258  for (unsigned int d = 2; d < dim; ++d)
259  dofs_per_face *= deg + 1;
260  // all face dofs need to be
261  // non-additive, since they have
262  // continuity requirements.
263  // however, the interior dofs are
264  // made additive
265  std::vector<bool> ret_val(dofs_per_cell, false);
266  for (unsigned int i = GeometryInfo<dim>::faces_per_cell * dofs_per_face;
267  i < dofs_per_cell;
268  ++i)
269  ret_val[i] = true;
270 
271  return ret_val;
272 }
273 
274 
275 template <int dim>
276 bool
278  const unsigned int shape_index,
279  const unsigned int face_index) const
280 {
281  Assert(shape_index < this->dofs_per_cell,
282  ExcIndexRange(shape_index, 0, this->dofs_per_cell));
285 
286  // The first degrees of freedom are
287  // on the faces and each face has
288  // degree degrees.
289  const unsigned int support_face = shape_index / this->degree;
290 
291  // The only thing we know for sure
292  // is that shape functions with
293  // support on one face are zero on
294  // the opposite face.
295  if (support_face < GeometryInfo<dim>::faces_per_cell)
296  return (face_index != GeometryInfo<dim>::opposite_face[support_face]);
297 
298  // In all other cases, return true,
299  // which is safe
300  return true;
301 }
302 
303 
304 
305 template <int dim>
306 void
309  const std::vector<Vector<double>> &support_point_values,
310  std::vector<double> & nodal_values) const
311 {
312  Assert(support_point_values.size() == this->generalized_support_points.size(),
313  ExcDimensionMismatch(support_point_values.size(),
314  this->generalized_support_points.size()));
315  Assert(nodal_values.size() == this->dofs_per_cell,
316  ExcDimensionMismatch(nodal_values.size(), this->dofs_per_cell));
317  Assert(support_point_values[0].size() == this->n_components(),
318  ExcDimensionMismatch(support_point_values[0].size(),
319  this->n_components()));
320 
321  // First do interpolation on
322  // faces. There, the component
323  // evaluated depends on the face
324  // direction and orientation.
325  unsigned int fbase = 0;
326  unsigned int f = 0;
327  for (; f < GeometryInfo<dim>::faces_per_cell;
328  ++f, fbase += this->dofs_per_face)
329  {
330  for (unsigned int i = 0; i < this->dofs_per_face; ++i)
331  {
332  nodal_values[fbase + i] = support_point_values[fbase + i](
334  }
335  }
336 
337  // The remaining points form dim
338  // chunks, one for each component.
339  const unsigned int istep = (this->dofs_per_cell - fbase) / dim;
340  Assert((this->dofs_per_cell - fbase) % dim == 0, ExcInternalError());
341 
342  f = 0;
343  while (fbase < this->dofs_per_cell)
344  {
345  for (unsigned int i = 0; i < istep; ++i)
346  {
347  nodal_values[fbase + i] = support_point_values[fbase + i](f);
348  }
349  fbase += istep;
350  ++f;
351  }
352  Assert(fbase == this->dofs_per_cell, ExcInternalError());
353 }
354 
355 
356 
357 // TODO: There are tests that check that the following few functions don't
358 // produce assertion failures, but none that actually check whether they do the
359 // right thing. one example for such a test would be to project a function onto
360 // an hp space and make sure that the convergence order is correct with regard
361 // to the lowest used polynomial degree
362 
363 template <int dim>
364 bool
366 {
367  return true;
368 }
369 
370 
371 template <int dim>
372 std::vector<std::pair<unsigned int, unsigned int>>
374  const FiniteElement<dim> &fe_other) const
375 {
376  // we can presently only compute these
377  // identities if both FEs are
378  // FE_RaviartThomasNodals or the other is FE_Nothing.
379  // In either case, no dofs are assigned on the vertex,
380  // so we shouldn't be getting here at all.
381  if (dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe_other) != nullptr)
382  return std::vector<std::pair<unsigned int, unsigned int>>();
383  else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
384  return std::vector<std::pair<unsigned int, unsigned int>>();
385  else
386  {
387  Assert(false, ExcNotImplemented());
388  return std::vector<std::pair<unsigned int, unsigned int>>();
389  }
390 }
391 
392 
393 
394 template <int dim>
395 std::vector<std::pair<unsigned int, unsigned int>>
397  const FiniteElement<dim> &fe_other) const
398 {
399  // we can presently only compute
400  // these identities if both FEs are
401  // FE_RaviartThomasNodals or if the other
402  // one is FE_Nothing
403  if (const FE_RaviartThomasNodal<dim> *fe_q_other =
404  dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe_other))
405  {
406  // dofs are located on faces; these are
407  // only lines in 2d
408  if (dim != 2)
409  return std::vector<std::pair<unsigned int, unsigned int>>();
410 
411  // dofs are located along lines, so two
412  // dofs are identical only if in the
413  // following two cases (remember that
414  // the face support points are Gauss
415  // points):
416  // 1. this->degree = fe_q_other->degree,
417  // in the case, all the dofs on
418  // the line are identical
419  // 2. this->degree-1 and fe_q_other->degree-1
420  // are both even, i.e. this->dof_per_line
421  // and fe_q_other->dof_per_line are both odd,
422  // there exists only one point (the middle one)
423  // such that dofs are identical on this point
424  //
425  // to understand this, note that
426  // this->degree is the *maximal*
427  // polynomial degree, and is thus one
428  // higher than the argument given to
429  // the constructor
430  const unsigned int p = this->degree - 1;
431  const unsigned int q = fe_q_other->degree - 1;
432 
433  std::vector<std::pair<unsigned int, unsigned int>> identities;
434 
435  if (p == q)
436  for (unsigned int i = 0; i < p + 1; ++i)
437  identities.emplace_back(i, i);
438 
439  else if (p % 2 == 0 && q % 2 == 0)
440  identities.emplace_back(p / 2, q / 2);
441 
442  return identities;
443  }
444  else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
445  {
446  // the FE_Nothing has no degrees of freedom, so there are no
447  // equivalencies to be recorded
448  return std::vector<std::pair<unsigned int, unsigned int>>();
449  }
450  else
451  {
452  Assert(false, ExcNotImplemented());
453  return std::vector<std::pair<unsigned int, unsigned int>>();
454  }
455 }
456 
457 
458 template <int dim>
459 std::vector<std::pair<unsigned int, unsigned int>>
461  const FiniteElement<dim> &fe_other) const
462 {
463  // we can presently only compute
464  // these identities if both FEs are
465  // FE_RaviartThomasNodals or if the other
466  // one is FE_Nothing
467  if (const FE_RaviartThomasNodal<dim> *fe_q_other =
468  dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe_other))
469  {
470  // dofs are located on faces; these are
471  // only quads in 3d
472  if (dim != 3)
473  return std::vector<std::pair<unsigned int, unsigned int>>();
474 
475  // this works exactly like the line
476  // case above
477  const unsigned int p = this->dofs_per_quad;
478  const unsigned int q = fe_q_other->dofs_per_quad;
479 
480  std::vector<std::pair<unsigned int, unsigned int>> identities;
481 
482  if (p == q)
483  for (unsigned int i = 0; i < p; ++i)
484  identities.emplace_back(i, i);
485 
486  else if (p % 2 != 0 && q % 2 != 0)
487  identities.emplace_back(p / 2, q / 2);
488 
489  return identities;
490  }
491  else if (dynamic_cast<const FE_Nothing<dim> *>(&fe_other) != nullptr)
492  {
493  // the FE_Nothing has no degrees of freedom, so there are no
494  // equivalencies to be recorded
495  return std::vector<std::pair<unsigned int, unsigned int>>();
496  }
497  else
498  {
499  Assert(false, ExcNotImplemented());
500  return std::vector<std::pair<unsigned int, unsigned int>>();
501  }
502 }
503 
504 
505 template <int dim>
508  const FiniteElement<dim> &fe_other,
509  const unsigned int codim) const
510 {
511  Assert(codim <= dim, ExcImpossibleInDim(dim));
512  (void)codim;
513 
514  // vertex/line/face/cell domination
515  // --------------------------------
516  if (const FE_RaviartThomasNodal<dim> *fe_rt_nodal_other =
517  dynamic_cast<const FE_RaviartThomasNodal<dim> *>(&fe_other))
518  {
519  if (this->degree < fe_rt_nodal_other->degree)
521  else if (this->degree == fe_rt_nodal_other->degree)
523  else
525  }
526  else if (const FE_Nothing<dim> *fe_nothing =
527  dynamic_cast<const FE_Nothing<dim> *>(&fe_other))
528  {
529  if (fe_nothing->is_dominating())
531  else
532  // the FE_Nothing has no degrees of freedom and it is typically used
533  // in a context where we don't require any continuity along the
534  // interface
536  }
537 
538  Assert(false, ExcNotImplemented());
540 }
541 
542 
543 
544 template <>
545 void
547  const FiniteElement<1, 1> & /*x_source_fe*/,
548  FullMatrix<double> & /*interpolation_matrix*/) const
549 {
550  Assert(false, ExcImpossibleInDim(1));
551 }
552 
553 
554 template <>
555 void
557  const FiniteElement<1, 1> & /*x_source_fe*/,
558  const unsigned int /*subface*/,
559  FullMatrix<double> & /*interpolation_matrix*/) const
560 {
561  Assert(false, ExcImpossibleInDim(1));
562 }
563 
564 
565 
566 template <int dim>
567 void
569  const FiniteElement<dim> &x_source_fe,
570  FullMatrix<double> & interpolation_matrix) const
571 {
572  // this is only implemented, if the
573  // source FE is also a
574  // RaviartThomasNodal element
575  AssertThrow((x_source_fe.get_name().find("FE_RaviartThomasNodal<") == 0) ||
576  (dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
577  &x_source_fe) != nullptr),
579 
580  Assert(interpolation_matrix.n() == this->dofs_per_face,
581  ExcDimensionMismatch(interpolation_matrix.n(), this->dofs_per_face));
582  Assert(interpolation_matrix.m() == x_source_fe.dofs_per_face,
583  ExcDimensionMismatch(interpolation_matrix.m(),
584  x_source_fe.dofs_per_face));
585 
586  // ok, source is a RaviartThomasNodal element, so
587  // we will be able to do the work
588  const FE_RaviartThomasNodal<dim> &source_fe =
589  dynamic_cast<const FE_RaviartThomasNodal<dim> &>(x_source_fe);
590 
591  // Make sure, that the element,
592  // for which the DoFs should be
593  // constrained is the one with
594  // the higher polynomial degree.
595  // Actually the procedure will work
596  // also if this assertion is not
597  // satisfied. But the matrices
598  // produced in that case might
599  // lead to problems in the
600  // hp procedures, which use this
601  // method.
602  Assert(this->dofs_per_face <= source_fe.dofs_per_face,
604 
605  // generate a quadrature
606  // with the generalized support points.
607  // This is later based as a
608  // basis for the QProjector,
609  // which returns the support
610  // points on the face.
611  Quadrature<dim - 1> quad_face_support(
613 
614  // Rule of thumb for FP accuracy,
615  // that can be expected for a
616  // given polynomial degree.
617  // This value is used to cut
618  // off values close to zero.
619  double eps = 2e-13 * this->degree * (dim - 1);
620 
621  // compute the interpolation
622  // matrix by simply taking the
623  // value at the support points.
624  const Quadrature<dim> face_projection =
625  QProjector<dim>::project_to_face(quad_face_support, 0);
626 
627  for (unsigned int i = 0; i < source_fe.dofs_per_face; ++i)
628  {
629  const Point<dim> &p = face_projection.point(i);
630 
631  for (unsigned int j = 0; j < this->dofs_per_face; ++j)
632  {
633  double matrix_entry =
634  this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
635 
636  // Correct the interpolated
637  // value. I.e. if it is close
638  // to 1 or 0, make it exactly
639  // 1 or 0. Unfortunately, this
640  // is required to avoid problems
641  // with higher order elements.
642  if (std::fabs(matrix_entry - 1.0) < eps)
643  matrix_entry = 1.0;
644  if (std::fabs(matrix_entry) < eps)
645  matrix_entry = 0.0;
646 
647  interpolation_matrix(i, j) = matrix_entry;
648  }
649  }
650 
651  // make sure that the row sum of
652  // each of the matrices is 1 at
653  // this point. this must be so
654  // since the shape functions sum up
655  // to 1
656  for (unsigned int j = 0; j < source_fe.dofs_per_face; ++j)
657  {
658  double sum = 0.;
659 
660  for (unsigned int i = 0; i < this->dofs_per_face; ++i)
661  sum += interpolation_matrix(j, i);
662 
663  Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
664  ExcInternalError());
665  }
666 }
667 
668 
669 template <int dim>
670 void
672  const FiniteElement<dim> &x_source_fe,
673  const unsigned int subface,
674  FullMatrix<double> & interpolation_matrix) const
675 {
676  // this is only implemented, if the
677  // source FE is also a
678  // RaviartThomasNodal element
679  AssertThrow((x_source_fe.get_name().find("FE_RaviartThomasNodal<") == 0) ||
680  (dynamic_cast<const FE_RaviartThomasNodal<dim> *>(
681  &x_source_fe) != nullptr),
683 
684  Assert(interpolation_matrix.n() == this->dofs_per_face,
685  ExcDimensionMismatch(interpolation_matrix.n(), this->dofs_per_face));
686  Assert(interpolation_matrix.m() == x_source_fe.dofs_per_face,
687  ExcDimensionMismatch(interpolation_matrix.m(),
688  x_source_fe.dofs_per_face));
689 
690  // ok, source is a RaviartThomasNodal element, so
691  // we will be able to do the work
692  const FE_RaviartThomasNodal<dim> &source_fe =
693  dynamic_cast<const FE_RaviartThomasNodal<dim> &>(x_source_fe);
694 
695  // Make sure, that the element,
696  // for which the DoFs should be
697  // constrained is the one with
698  // the higher polynomial degree.
699  // Actually the procedure will work
700  // also if this assertion is not
701  // satisfied. But the matrices
702  // produced in that case might
703  // lead to problems in the
704  // hp procedures, which use this
705  // method.
706  Assert(this->dofs_per_face <= source_fe.dofs_per_face,
708 
709  // generate a quadrature
710  // with the generalized support points.
711  // This is later based as a
712  // basis for the QProjector,
713  // which returns the support
714  // points on the face.
715  Quadrature<dim - 1> quad_face_support(
717 
718  // Rule of thumb for FP accuracy,
719  // that can be expected for a
720  // given polynomial degree.
721  // This value is used to cut
722  // off values close to zero.
723  double eps = 2e-13 * this->degree * (dim - 1);
724 
725  // compute the interpolation
726  // matrix by simply taking the
727  // value at the support points.
728 
729  const Quadrature<dim> subface_projection =
730  QProjector<dim>::project_to_subface(quad_face_support, 0, subface);
731 
732  for (unsigned int i = 0; i < source_fe.dofs_per_face; ++i)
733  {
734  const Point<dim> &p = subface_projection.point(i);
735 
736  for (unsigned int j = 0; j < this->dofs_per_face; ++j)
737  {
738  double matrix_entry =
739  this->shape_value_component(this->face_to_cell_index(j, 0), p, 0);
740 
741  // Correct the interpolated
742  // value. I.e. if it is close
743  // to 1 or 0, make it exactly
744  // 1 or 0. Unfortunately, this
745  // is required to avoid problems
746  // with higher order elements.
747  if (std::fabs(matrix_entry - 1.0) < eps)
748  matrix_entry = 1.0;
749  if (std::fabs(matrix_entry) < eps)
750  matrix_entry = 0.0;
751 
752  interpolation_matrix(i, j) = matrix_entry;
753  }
754  }
755 
756  // make sure that the row sum of
757  // each of the matrices is 1 at
758  // this point. this must be so
759  // since the shape functions sum up
760  // to 1
761  for (unsigned int j = 0; j < source_fe.dofs_per_face; ++j)
762  {
763  double sum = 0.;
764 
765  for (unsigned int i = 0; i < this->dofs_per_face; ++i)
766  sum += interpolation_matrix(j, i);
767 
768  Assert(std::fabs(sum - 1) < 2e-13 * this->degree * (dim - 1),
769  ExcInternalError());
770  }
771 }
772 
773 
774 
775 // explicit instantiations
776 #include "fe_raviart_thomas_nodal.inst"
777 
778 
779 DEAL_II_NAMESPACE_CLOSE
size_type m() const
static unsigned int n_polynomials(const unsigned int degree)
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2493
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)
FE_RaviartThomasNodal(const unsigned int p)
const unsigned int dofs_per_quad
Definition: fe_base.h:237
const unsigned int degree
Definition: fe_base.h:298
const Point< dim > & point(const unsigned int i) const
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const override
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
void invert(const FullMatrix< number2 > &M)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1519
FullMatrix< double > inverse_node_matrix
FullMatrix< double > compute_node_matrix(const FiniteElement< dim, spacedim > &fe)
std::vector< Point< dim - 1 > > generalized_face_support_points
Definition: fe.h:2499
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual bool hp_constraints_are_implemented() const override
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2456
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:1407
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void initialize_support_points(const unsigned int rt_degree)
virtual std::string get_name() const =0
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
Definition: fe.cc:550
const unsigned int dofs_per_cell
Definition: fe_base.h:282
virtual std::string get_name() const override
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim >> &q_points)
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim > &fe_other, const unsigned int codim=0) const override final
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
unsigned int n_components() const
const unsigned int dofs_per_face
Definition: fe_base.h:275
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInternalError()
static std::vector< bool > get_ria_vector(const unsigned int degree)
const std::vector< Point< dim - 1 > > & get_generalized_face_support_points() const
Definition: fe.cc:1093
std::vector< MappingKind > mapping_kind