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