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