Reference documentation for deal.II version 8.4.2
fe_raviart_thomas.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/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  // Now compute the inverse node
63  //matrix, generating the correct
64  //basis functions from the raw
65  //ones.
66 
67  // We use an auxiliary matrix in
68  // this function. Therefore,
69  // inverse_node_matrix is still
70  // empty and shape_value_component
71  // returns the 'raw' shape values.
72  FullMatrix<double> M(n_dofs, n_dofs);
74  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
75  this->inverse_node_matrix.invert(M);
76  // From now on, the shape functions
77  // will be the correct ones, not
78  // the raw shape functions anymore.
79 
80  // Reinit the vectors of
81  // restriction and prolongation
82  // matrices to the right sizes.
83  // Restriction only for isotropic
84  // refinement
86  // Fill prolongation matrices with embedding operators
89 
90  // TODO[TL]: for anisotropic refinement we will probably need a table of submatrices with an array for each refine case
92  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
93  face_embeddings[i].reinit (this->dofs_per_face, this->dofs_per_face);
94  FETools::compute_face_embedding_matrices<dim,double>(*this, face_embeddings, 0, 0);
95  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
96  this->dofs_per_face);
97  unsigned int target_row=0;
98  for (unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
99  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
100  {
101  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
102  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
103  ++target_row;
104  }
105 }
106 
107 
108 
109 template <int dim>
110 std::string
112 {
113  // note that the
114  // FETools::get_fe_from_name
115  // function depends on the
116  // particular format of the string
117  // this function returns, so they
118  // have to be kept in synch
119 
120  // note that this->degree is the maximal
121  // polynomial degree and is thus one higher
122  // than the argument given to the
123  // constructor
124  std::ostringstream namebuf;
125  namebuf << "FE_RaviartThomas<" << dim << ">(" << this->degree-1 << ")";
126 
127  return namebuf.str();
128 }
129 
130 
131 template <int dim>
134 {
135  return new FE_RaviartThomas<dim>(*this);
136 }
137 
138 
139 //---------------------------------------------------------------------------
140 // Auxiliary and internal functions
141 //---------------------------------------------------------------------------
142 
143 
144 template <int dim>
145 void
147 {
148  QGauss<dim> cell_quadrature(deg+1);
149  const unsigned int n_interior_points
150  = (deg>0) ? cell_quadrature.size() : 0;
151 
152  unsigned int n_face_points = (dim>1) ? 1 : 0;
153  // compute (deg+1)^(dim-1)
154  for (unsigned int d=1; d<dim; ++d)
155  n_face_points *= deg+1;
156 
157 
159  + n_interior_points);
160  this->generalized_face_support_points.resize (n_face_points);
161 
162  // Number of the point being entered
163  unsigned int current = 0;
164 
165  if (dim>1)
166  {
167  QGauss<dim-1> face_points (deg+1);
168  TensorProductPolynomials<dim-1> legendre
170 
171  boundary_weights.reinit(n_face_points, legendre.n());
172 
173 // Assert (face_points.size() == this->dofs_per_face,
174 // ExcInternalError());
175 
176  for (unsigned int k=0; k<n_face_points; ++k)
177  {
178  this->generalized_face_support_points[k] = face_points.point(k);
179  // Compute its quadrature
180  // contribution for each
181  // moment.
182  for (unsigned int i=0; i<legendre.n(); ++i)
183  {
184  boundary_weights(k, i)
185  = face_points.weight(k)
186  * legendre.compute_value(i, face_points.point(k));
187  }
188  }
189 
191  for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
192  ++current)
193  {
194  // Enter the support point
195  // into the vector
196  this->generalized_support_points[current] = faces.point(current+QProjector<dim>::DataSetDescriptor::face(0,true,false,false,n_face_points));
197  }
198  }
199 
200  if (deg==0) return;
201 
202  // Create Legendre basis for the
203  // space D_xi Q_k
204  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
205  for (unsigned int dd=0; dd<dim; ++dd)
206  {
207  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
208  for (unsigned int d=0; d<dim; ++d)
211 
212  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
213  }
214 
215  interior_weights.reinit(TableIndices<3>(n_interior_points, polynomials[0]->n(), dim));
216 
217  for (unsigned int k=0; k<cell_quadrature.size(); ++k)
218  {
219  this->generalized_support_points[current++] = cell_quadrature.point(k);
220  for (unsigned int i=0; i<polynomials[0]->n(); ++i)
221  for (unsigned int d=0; d<dim; ++d)
222  interior_weights(k,i,d) = cell_quadrature.weight(k)
223  * polynomials[d]->compute_value(i,cell_quadrature.point(k));
224  }
225 
226  for (unsigned int d=0; d<dim; ++d)
227  delete polynomials[d];
228 
229  Assert (current == this->generalized_support_points.size(),
230  ExcInternalError());
231 }
232 
233 
234 
235 template <>
236 void
238 {
239  // there is only one refinement case in 1d,
240  // which is the isotropic one (first index of
241  // the matrix array has to be 0)
242  for (unsigned int i=0; i<GeometryInfo<1>::max_children_per_cell; ++i)
243  this->restriction[0][i].reinit(0,0);
244 }
245 
246 
247 
248 // This function is the same Raviart-Thomas interpolation performed by
249 // interpolate. Still, we cannot use interpolate, since it was written
250 // for smooth functions. The functions interpolated here are not
251 // smooth, maybe even not continuous. Therefore, we must double the
252 // number of quadrature points in each direction in order to integrate
253 // only smooth functions.
254 
255 // Then again, the interpolated function is chosen such that the
256 // moments coincide with the function to be interpolated.
257 
258 template <int dim>
259 void
261 {
262  const unsigned int iso=RefinementCase<dim>::isotropic_refinement-1;
263 
264  QGauss<dim-1> q_base (this->degree);
265  const unsigned int n_face_points = q_base.size();
266  // First, compute interpolation on
267  // subfaces
268  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
269  {
270  // The shape functions of the
271  // child cell are evaluated
272  // in the quadrature points
273  // of a full face.
274  Quadrature<dim> q_face
275  = QProjector<dim>::project_to_face(q_base, face);
276  // Store shape values, since the
277  // evaluation suffers if not
278  // ordered by point
280  for (unsigned int k=0; k<q_face.size(); ++k)
281  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
282  cached_values(i,k)
283  = this->shape_value_component(i, q_face.point(k),
285 
286  for (unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
287  {
288  // The weight functions for
289  // the coarse face are
290  // evaluated on the subface
291  // only.
292  Quadrature<dim> q_sub
293  = QProjector<dim>::project_to_subface(q_base, face, sub);
294  const unsigned int child
297 
298  // On a certain face, we must
299  // compute the moments of ALL
300  // fine level functions with
301  // the coarse level weight
302  // functions belonging to
303  // that face. Due to the
304  // orthogonalization process
305  // when building the shape
306  // functions, these weights
307  // are equal to the
308  // corresponding shape
309  // functions.
310  for (unsigned int k=0; k<n_face_points; ++k)
311  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
312  for (unsigned int i_face = 0; i_face < this->dofs_per_face; ++i_face)
313  {
314  // The quadrature
315  // weights on the
316  // subcell are NOT
317  // transformed, so we
318  // have to do it here.
319  this->restriction[iso][child](face*this->dofs_per_face+i_face,
320  i_child)
321  += Utilities::fixed_power<dim-1>(.5) * q_sub.weight(k)
322  * cached_values(i_child, k)
323  * this->shape_value_component(face*this->dofs_per_face+i_face,
324  q_sub.point(k),
326  }
327  }
328  }
329 
330  if (this->degree == 1) return;
331 
332  // Create Legendre basis for the
333  // space D_xi Q_k. Here, we cannot
334  // use the shape functions
335  std::vector<AnisotropicPolynomials<dim>* > polynomials(dim);
336  for (unsigned int dd=0; dd<dim; ++dd)
337  {
338  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
339  for (unsigned int d=0; d<dim; ++d)
342 
343  polynomials[dd] = new AnisotropicPolynomials<dim>(poly);
344  }
345 
346  QGauss<dim> q_cell(this->degree);
347  const unsigned int start_cell_dofs
349 
350  // Store shape values, since the
351  // evaluation suffers if not
352  // ordered by point
353  Table<3,double> cached_values(this->dofs_per_cell, q_cell.size(), dim);
354  for (unsigned int k=0; k<q_cell.size(); ++k)
355  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
356  for (unsigned int d=0; d<dim; ++d)
357  cached_values(i,k,d) = this->shape_value_component(i, q_cell.point(k), d);
358 
359  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
360  {
361  Quadrature<dim> q_sub = QProjector<dim>::project_to_child(q_cell, child);
362 
363  for (unsigned int k=0; k<q_sub.size(); ++k)
364  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
365  for (unsigned int d=0; d<dim; ++d)
366  for (unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
367  {
368  this->restriction[iso][child](start_cell_dofs+i_weight*dim+d,
369  i_child)
370  += q_sub.weight(k)
371  * cached_values(i_child, k, d)
372  * polynomials[d]->compute_value(i_weight, q_sub.point(k));
373  }
374  }
375 
376  for (unsigned int d=0; d<dim; ++d)
377  delete polynomials[d];
378 }
379 
380 
381 
382 template <int dim>
383 std::vector<unsigned int>
385 {
386  // the element is face-based and we have
387  // (deg+1)^(dim-1) DoFs per face
388  unsigned int dofs_per_face = 1;
389  for (unsigned int d=1; d<dim; ++d)
390  dofs_per_face *= deg+1;
391 
392  // and then there are interior dofs
393  const unsigned int
394  interior_dofs = dim*deg*dofs_per_face;
395 
396  std::vector<unsigned int> dpo(dim+1);
397  dpo[dim-1] = dofs_per_face;
398  dpo[dim] = interior_dofs;
399 
400  return dpo;
401 }
402 
403 
404 
405 template <int dim>
406 std::pair<Table<2,bool>, std::vector<unsigned int> >
408 {
409  Table<2,bool> constant_modes(dim, this->dofs_per_cell);
410  for (unsigned int d=0; d<dim; ++d)
411  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
412  constant_modes(d,i) = true;
413  std::vector<unsigned int> components;
414  for (unsigned int d=0; d<dim; ++d)
415  components.push_back(d);
416  return std::pair<Table<2,bool>, std::vector<unsigned int> >
417  (constant_modes, components);
418 }
419 
420 
421 
422 //---------------------------------------------------------------------------
423 // Data field initialization
424 //---------------------------------------------------------------------------
425 
426 
427 template <int dim>
428 bool
430  const unsigned int shape_index,
431  const unsigned int face_index) const
432 {
433  Assert (shape_index < this->dofs_per_cell,
434  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
436  ExcIndexRange (face_index, 0, GeometryInfo<dim>::faces_per_cell));
437 
438  // Return computed values if we
439  // know them easily. Otherwise, it
440  // is always safe to return true.
441  switch (this->degree)
442  {
443  case 1:
444  {
445  switch (dim)
446  {
447  case 2:
448  {
449  // only on the one
450  // non-adjacent face
451  // are the values
452  // actually zero. list
453  // these in a table
454  return (face_index != GeometryInfo<dim>::opposite_face[shape_index]);
455  }
456 
457  default:
458  return true;
459  };
460  };
461 
462  default: // other rt_order
463  return true;
464  };
465 
466  return true;
467 }
468 
469 
470 // Since this is a vector valued element, we cannot interpolate a
471 // scalar function
472 template <int dim>
473 void
475  std::vector<double> &,
476  const std::vector<double> &) const
477 {
478  Assert(false, ExcNotImplemented());
479 }
480 
481 
482 template <int dim>
483 void
485  std::vector<double> &local_dofs,
486  const std::vector<Vector<double> > &values,
487  unsigned int offset) const
488 {
489  Assert (values.size() == this->generalized_support_points.size(),
490  ExcDimensionMismatch(values.size(), this->generalized_support_points.size()));
491  Assert (local_dofs.size() == this->dofs_per_cell,
492  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
493  Assert (values[0].size() >= offset+this->n_components(),
494  ExcDimensionMismatch(values[0].size(),offset+this->n_components()));
495 
496  std::fill(local_dofs.begin(), local_dofs.end(), 0.);
497 
498  const unsigned int n_face_points = boundary_weights.size(0);
499  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
500  for (unsigned int k=0; k<n_face_points; ++k)
501  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
502  {
503  local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
504  * values[face*n_face_points+k](GeometryInfo<dim>::unit_normal_direction[face]+offset);
505  }
506 
507  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
508  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
509 
510  for (unsigned int k=0; k<interior_weights.size(0); ++k)
511  for (unsigned int i=0; i<interior_weights.size(1); ++i)
512  for (unsigned int d=0; d<dim; ++d)
513  local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[k+start_cell_points](d+offset);
514 }
515 
516 
517 template <int dim>
518 void
520  std::vector<double> &local_dofs,
521  const VectorSlice<const std::vector<std::vector<double> > > &values) const
522 {
523  Assert (values.size() == this->n_components(),
524  ExcDimensionMismatch(values.size(), this->n_components()));
525  Assert (values[0].size() == this->generalized_support_points.size(),
526  ExcDimensionMismatch(values[0].size(), this->generalized_support_points.size()));
527  Assert (local_dofs.size() == this->dofs_per_cell,
528  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
529 
530  std::fill(local_dofs.begin(), local_dofs.end(), 0.);
531 
532  const unsigned int n_face_points = boundary_weights.size(0);
533  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
534  for (unsigned int k=0; k<n_face_points; ++k)
535  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
536  {
537  local_dofs[i+face*this->dofs_per_face] += boundary_weights(k,i)
538  * values[GeometryInfo<dim>::unit_normal_direction[face]][face*n_face_points+k];
539  }
540 
541  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
542  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
543 
544  for (unsigned int k=0; k<interior_weights.size(0); ++k)
545  for (unsigned int i=0; i<interior_weights.size(1); ++i)
546  for (unsigned int d=0; d<dim; ++d)
547  local_dofs[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * values[d][k+start_cell_points];
548 }
549 
550 
551 
552 template <int dim>
553 std::size_t
555 {
556  Assert (false, ExcNotImplemented ());
557  return 0;
558 }
559 
560 
561 
562 // explicit instantiations
563 #include "fe_raviart_thomas.inst"
564 
565 
566 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
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2106
const unsigned int components
Definition: fe_base.h:293
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2157
FullMatrix< double > interface_constraints
Definition: fe.h:2132
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:1009
const unsigned int degree
Definition: fe_base.h:299
const Point< dim > & point(const unsigned int i) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() 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
friend class FE_RaviartThomas
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2120
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: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::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: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
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:278
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
unsigned int size(const unsigned int i) const
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
void compute_node_matrix(FullMatrix< double > &M, const FiniteElement< dim, spacedim > &fe)
Definition: fe_tools.cc:633
Table< 3, double > interior_weights
double weight(const unsigned int i) const
Table< 2, double > boundary_weights