Reference documentation for deal.II version 9.0.0
fe_raviart_thomas.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/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 #include <deal.II/base/std_cxx14/memory.h>
34 
35 //TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
36 //adjust_line_dof_index_for_line_orientation_table fields, and write tests
37 //similar to bits/face_orientation_and_fe_q_*
38 
39 
40 DEAL_II_NAMESPACE_OPEN
41 
42 
43 template <int dim>
45  :
47  deg,
48  FiniteElementData<dim>(get_dpo_vector(deg),
49  dim, deg+1, FiniteElementData<dim>::Hdiv),
50  std::vector<bool>(PolynomialsRaviartThomas<dim>::compute_n_pols(deg), true),
51  std::vector<ComponentMask>(PolynomialsRaviartThomas<dim>::compute_n_pols(deg),
52  std::vector<bool>(dim,true)))
53 {
54  Assert (dim >= 2, ExcImpossibleInDim(dim));
55  const unsigned int n_dofs = this->dofs_per_cell;
56 
58  // First, initialize the
59  // generalized support points and
60  // quadrature weights, since they
61  // are required for interpolation.
63 
64  // Now compute the inverse node matrix, generating the correct
65  // basis functions from the raw ones. For a discussion of what
66  // exactly happens here, see FETools::compute_node_matrix.
68  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
69  this->inverse_node_matrix.invert(M);
70  // From now on, the shape functions provided by FiniteElement::shape_value
71  // and similar functions will be the correct ones, not
72  // the raw shape functions from the polynomial space anymore.
73 
74  // Reinit the vectors of
75  // restriction and prolongation
76  // matrices to the right sizes.
77  // Restriction only for isotropic
78  // refinement
80  // Fill prolongation matrices with embedding operators
83 
84  // TODO[TL]: for anisotropic refinement we will probably need a table of submatrices with an array for each refine case
86  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
87  face_embeddings[i].reinit (this->dofs_per_face, this->dofs_per_face);
88  FETools::compute_face_embedding_matrices<dim,double>(*this, face_embeddings, 0, 0);
89  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
90  this->dofs_per_face);
91  unsigned int target_row=0;
92  for (unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
93  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
94  {
95  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
96  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
97  ++target_row;
98  }
99 }
100 
101 
102 
103 template <int dim>
104 std::string
106 {
107  // note that the
108  // FETools::get_fe_by_name
109  // function depends on the
110  // particular format of the string
111  // this function returns, so they
112  // have to be kept in synch
113 
114  // note that this->degree is the maximal
115  // polynomial degree and is thus one higher
116  // than the argument given to the
117  // constructor
118  std::ostringstream namebuf;
119  namebuf << "FE_RaviartThomas<" << dim << ">(" << this->degree-1 << ")";
120 
121  return namebuf.str();
122 }
123 
124 
125 template <int dim>
126 std::unique_ptr<FiniteElement<dim,dim> >
128 {
129  return std_cxx14::make_unique<FE_RaviartThomas<dim>>(*this);
130 }
131 
132 
133 //---------------------------------------------------------------------------
134 // Auxiliary and internal functions
135 //---------------------------------------------------------------------------
136 
137 
138 template <int dim>
139 void
141 {
142  QGauss<dim> cell_quadrature(deg+1);
143  const unsigned int n_interior_points
144  = (deg>0) ? cell_quadrature.size() : 0;
145 
146  unsigned int n_face_points = (dim>1) ? 1 : 0;
147  // compute (deg+1)^(dim-1)
148  for (unsigned int d=1; d<dim; ++d)
149  n_face_points *= deg+1;
150 
151 
152  this->generalized_support_points.resize (GeometryInfo<dim>::faces_per_cell*n_face_points
153  + n_interior_points);
154  this->generalized_face_support_points.resize (n_face_points);
155 
156  // Number of the point being entered
157  unsigned int current = 0;
158 
159  if (dim>1)
160  {
161  QGauss<dim-1> face_points (deg+1);
162  TensorProductPolynomials<dim-1> legendre
164 
165  boundary_weights.reinit(n_face_points, legendre.n());
166 
167 // Assert (face_points.size() == this->dofs_per_face,
168 // ExcInternalError());
169 
170  for (unsigned int k=0; k<n_face_points; ++k)
171  {
172  this->generalized_face_support_points[k] = face_points.point(k);
173  // Compute its quadrature
174  // contribution for each
175  // moment.
176  for (unsigned int i=0; i<legendre.n(); ++i)
177  {
178  boundary_weights(k, i)
179  = face_points.weight(k)
180  * legendre.compute_value(i, face_points.point(k));
181  }
182  }
183 
185  for (; current<GeometryInfo<dim>::faces_per_cell*n_face_points;
186  ++current)
187  {
188  // Enter the support point
189  // into the vector
190  this->generalized_support_points[current] = faces.point(current+QProjector<dim>::DataSetDescriptor::face(0,true,false,false,n_face_points));
191  }
192  }
193 
194  if (deg==0) return;
195 
196  // Create Legendre basis for the space D_xi Q_k
197  std::unique_ptr<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] = std_cxx14::make_unique<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  Assert (current == this->generalized_support_points.size(),
220  ExcInternalError());
221 }
222 
223 
224 
225 template <>
226 void
228 {
229  // there is only one refinement case in 1d,
230  // which is the isotropic one (first index of
231  // the matrix array has to be 0)
232  for (unsigned int i=0; i<GeometryInfo<1>::max_children_per_cell; ++i)
233  this->restriction[0][i].reinit(0,0);
234 }
235 
236 
237 
238 // This function is the same Raviart-Thomas interpolation performed by
239 // interpolate. Still, we cannot use interpolate, since it was written
240 // for smooth functions. The functions interpolated here are not
241 // smooth, maybe even not continuous. Therefore, we must double the
242 // number of quadrature points in each direction in order to integrate
243 // only smooth functions.
244 
245 // Then again, the interpolated function is chosen such that the
246 // moments coincide with the function to be interpolated.
247 
248 template <int dim>
249 void
251 {
252  const unsigned int iso=RefinementCase<dim>::isotropic_refinement-1;
253 
254  QGauss<dim-1> q_base (this->degree);
255  const unsigned int n_face_points = q_base.size();
256  // First, compute interpolation on
257  // subfaces
258  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
259  {
260  // The shape functions of the
261  // child cell are evaluated
262  // in the quadrature points
263  // of a full face.
264  Quadrature<dim> q_face
265  = QProjector<dim>::project_to_face(q_base, face);
266  // Store shape values, since the
267  // evaluation suffers if not
268  // ordered by point
269  Table<2,double> cached_values_on_face(this->dofs_per_cell, q_face.size());
270  for (unsigned int k=0; k<q_face.size(); ++k)
271  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
272  cached_values_on_face(i,k)
273  = this->shape_value_component(i, q_face.point(k),
275 
276  for (unsigned int sub=0; sub<GeometryInfo<dim>::max_children_per_face; ++sub)
277  {
278  // The weight functions for
279  // the coarse face are
280  // evaluated on the subface
281  // only.
282  Quadrature<dim> q_sub
283  = QProjector<dim>::project_to_subface(q_base, face, sub);
284  const unsigned int child
287 
288  // On a certain face, we must
289  // compute the moments of ALL
290  // fine level functions with
291  // the coarse level weight
292  // functions belonging to
293  // that face. Due to the
294  // orthogonalization process
295  // when building the shape
296  // functions, these weights
297  // are equal to the
298  // corresponding shape
299  // functions.
300  for (unsigned int k=0; k<n_face_points; ++k)
301  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
302  for (unsigned int i_face = 0; i_face < this->dofs_per_face; ++i_face)
303  {
304  // The quadrature
305  // weights on the
306  // subcell are NOT
307  // transformed, so we
308  // have to do it here.
309  this->restriction[iso][child](face*this->dofs_per_face+i_face,
310  i_child)
311  += Utilities::fixed_power<dim-1>(.5) * q_sub.weight(k)
312  * cached_values_on_face(i_child, k)
313  * this->shape_value_component(face*this->dofs_per_face+i_face,
314  q_sub.point(k),
316  }
317  }
318  }
319 
320  if (this->degree == 1) return;
321 
322  // Create Legendre basis for the space D_xi Q_k. Here, we cannot
323  // use the shape functions
324  std::unique_ptr<AnisotropicPolynomials<dim>> polynomials[dim];
325  for (unsigned int dd=0; dd<dim; ++dd)
326  {
327  std::vector<std::vector<Polynomials::Polynomial<double> > > poly(dim);
328  for (unsigned int d=0; d<dim; ++d)
329  poly[d] = Polynomials::Legendre::generate_complete_basis(this->degree-1);
330  poly[dd] = Polynomials::Legendre::generate_complete_basis(this->degree-2);
331 
332  polynomials[dd] = std_cxx14::make_unique<AnisotropicPolynomials<dim>>(poly);
333  }
334 
335  QGauss<dim> q_cell(this->degree);
336  const unsigned int start_cell_dofs
337  = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
338 
339  // Store shape values, since the
340  // evaluation suffers if not
341  // ordered by point
342  Table<3,double> cached_values_on_cell(this->dofs_per_cell, q_cell.size(), dim);
343  for (unsigned int k=0; k<q_cell.size(); ++k)
344  for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
345  for (unsigned int d=0; d<dim; ++d)
346  cached_values_on_cell(i,k,d) = this->shape_value_component(i, q_cell.point(k), d);
347 
348  for (unsigned int child=0; child<GeometryInfo<dim>::max_children_per_cell; ++child)
349  {
350  Quadrature<dim> q_sub = QProjector<dim>::project_to_child(q_cell, child);
351 
352  for (unsigned int k=0; k<q_sub.size(); ++k)
353  for (unsigned int i_child = 0; i_child < this->dofs_per_cell; ++i_child)
354  for (unsigned int d=0; d<dim; ++d)
355  for (unsigned int i_weight=0; i_weight<polynomials[d]->n(); ++i_weight)
356  {
357  this->restriction[iso][child](start_cell_dofs+i_weight*dim+d,
358  i_child)
359  += q_sub.weight(k)
360  * cached_values_on_cell(i_child, k, d)
361  * polynomials[d]->compute_value(i_weight, q_sub.point(k));
362  }
363  }
364 }
365 
366 
367 
368 template <int dim>
369 std::vector<unsigned int>
371 {
372  // the element is face-based and we have
373  // (deg+1)^(dim-1) DoFs per face
374  unsigned int dofs_per_face = 1;
375  for (unsigned int d=1; d<dim; ++d)
376  dofs_per_face *= deg+1;
377 
378  // and then there are interior dofs
379  const unsigned int
380  interior_dofs = dim*deg*dofs_per_face;
381 
382  std::vector<unsigned int> dpo(dim+1);
383  dpo[dim-1] = dofs_per_face;
384  dpo[dim] = interior_dofs;
385 
386  return dpo;
387 }
388 
389 
390 
391 template <int dim>
392 std::pair<Table<2,bool>, std::vector<unsigned int> >
394 {
395  Table<2,bool> constant_modes(dim, this->dofs_per_cell);
396  for (unsigned int d=0; d<dim; ++d)
397  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
398  constant_modes(d,i) = true;
399  std::vector<unsigned int> components;
400  for (unsigned int d=0; d<dim; ++d)
401  components.push_back(d);
402  return std::pair<Table<2,bool>, std::vector<unsigned int> >
403  (constant_modes, components);
404 }
405 
406 
407 
408 //---------------------------------------------------------------------------
409 // Data field initialization
410 //---------------------------------------------------------------------------
411 
412 
413 template <int dim>
414 bool
416  const unsigned int shape_index,
417  const unsigned int face_index) const
418 {
419  Assert (shape_index < this->dofs_per_cell,
420  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
423 
424  // Return computed values if we
425  // know them easily. Otherwise, it
426  // is always safe to return true.
427  switch (this->degree)
428  {
429  case 1:
430  {
431  switch (dim)
432  {
433  case 2:
434  {
435  // only on the one
436  // non-adjacent face
437  // are the values
438  // actually zero. list
439  // these in a table
440  return (face_index != GeometryInfo<dim>::opposite_face[shape_index]);
441  }
442 
443  default:
444  return true;
445  };
446  };
447 
448  default: // other rt_order
449  return true;
450  };
451 
452  return true;
453 }
454 
455 
456 
457 template <int dim>
458 void
461  std::vector<double> &nodal_values) const
462 {
463  Assert (support_point_values.size() == this->generalized_support_points.size(),
464  ExcDimensionMismatch(support_point_values.size(), this->generalized_support_points.size()));
465  Assert (nodal_values.size() == this->dofs_per_cell,
466  ExcDimensionMismatch(nodal_values.size(),this->dofs_per_cell));
467  Assert (support_point_values[0].size() == this->n_components(),
468  ExcDimensionMismatch(support_point_values[0].size(),this->n_components()));
469 
470  std::fill(nodal_values.begin(), nodal_values.end(), 0.);
471 
472  const unsigned int n_face_points = boundary_weights.size(0);
473  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
474  for (unsigned int k=0; k<n_face_points; ++k)
475  for (unsigned int i=0; i<boundary_weights.size(1); ++i)
476  {
477  nodal_values[i+face*this->dofs_per_face] += boundary_weights(k,i)
478  * support_point_values[face*n_face_points+k](GeometryInfo<dim>::unit_normal_direction[face]);
479  }
480 
481  const unsigned int start_cell_dofs = GeometryInfo<dim>::faces_per_cell*this->dofs_per_face;
482  const unsigned int start_cell_points = GeometryInfo<dim>::faces_per_cell*n_face_points;
483 
484  for (unsigned int k=0; k<interior_weights.size(0); ++k)
485  for (unsigned int i=0; i<interior_weights.size(1); ++i)
486  for (unsigned int d=0; d<dim; ++d)
487  nodal_values[start_cell_dofs+i*dim+d] += interior_weights(k,i,d) * support_point_values[k+start_cell_points](d);
488 }
489 
490 
491 
492 template <int dim>
493 std::size_t
495 {
496  Assert (false, ExcNotImplemented ());
497  return 0;
498 }
499 
500 
501 
502 // explicit instantiations
503 #include "fe_raviart_thomas.inst"
504 
505 
506 DEAL_II_NAMESPACE_CLOSE
virtual std::string get_name() const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
static Quadrature< dim > project_to_child(const Quadrature< dim > &quadrature, const unsigned int child_no)
Definition: quadrature.cc:1096
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
FullMatrix< double > interface_constraints
Definition: fe.h:2401
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:874
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)
FullMatrix< double > compute_node_matrix(const FiniteElement< dim, spacedim > &fe)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
friend class FE_RaviartThomas
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
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: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)
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)
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:274
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()