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