Reference documentation for deal.II version 9.0.0
fe_bdm.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/quadrature_lib.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/polynomials_p.h>
20 #include <deal.II/base/table.h>
21 #include <deal.II/grid/tria.h>
22 #include <deal.II/grid/tria_iterator.h>
23 #include <deal.II/dofs/dof_accessor.h>
24 #include <deal.II/fe/fe.h>
25 #include <deal.II/fe/mapping.h>
26 #include <deal.II/fe/fe_bdm.h>
27 #include <deal.II/fe/fe_values.h>
28 #include <deal.II/fe/fe_tools.h>
29 
30 #include <iostream>
31 #include <sstream>
32 #include <deal.II/base/std_cxx14/memory.h>
33 
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 template <int dim>
38 FE_BDM<dim>::FE_BDM (const unsigned int deg)
39  :
40  FE_PolyTensor<PolynomialsBDM<dim>, dim> (
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>(PolynomialsBDM<dim>::compute_n_pols(deg),
46  std::vector<bool>(dim,true)))
47 {
48  Assert (dim >= 2, ExcImpossibleInDim(dim));
49  Assert (deg > 0, ExcMessage("Lowest order BDM element are degree 1, but you asked for degree 0"));
50 
51  const unsigned int n_dofs = this->dofs_per_cell;
52 
53  this->mapping_type = mapping_bdm;
54  // These must be done first, since
55  // they change the evaluation of
56  // basis functions
57 
58  // Set up the generalized support
59  // points
61 
62  // Now compute the inverse node matrix, generating the correct
63  // basis functions from the raw ones. For a discussion of what
64  // exactly happens here, see FETools::compute_node_matrix.
66  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
67  this->inverse_node_matrix.invert(M);
68  // From now on, the shape functions provided by FiniteElement::shape_value
69  // and similar functions will be the correct ones, not
70  // the raw shape functions from the polynomial space anymore.
71 
72  // Embedding errors become pretty large, so we just replace the
73  // regular threshold in both "computing_..." functions by 1.
75  FETools::compute_embedding_matrices (*this, this->prolongation, true, 1.);
76 
78  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
79  face_embeddings[i].reinit (this->dofs_per_face, this->dofs_per_face);
80  FETools::compute_face_embedding_matrices(*this, face_embeddings, 0, 0, 1.);
81  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
82  this->dofs_per_face);
83  unsigned int target_row=0;
84  for (unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
85  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
86  {
87  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
88  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
89  ++target_row;
90  }
91 }
92 
93 
94 
95 template <int dim>
96 std::string
98 {
99  // note that the
100  // FETools::get_fe_by_name
101  // function depends on the
102  // particular format of the string
103  // this function returns, so they
104  // have to be kept in synch
105 
106  // note that this->degree is the maximal
107  // polynomial degree and is thus one higher
108  // than the argument given to the
109  // constructor
110  std::ostringstream namebuf;
111  namebuf << "FE_BDM<" << dim << ">(" << this->degree-1 << ")";
112 
113  return namebuf.str();
114 }
115 
116 
117 template <int dim>
118 std::unique_ptr<FiniteElement<dim,dim> >
120 {
121  return std_cxx14::make_unique<FE_BDM<dim>>(*this);
122 }
123 
124 
125 
126 template <int dim>
127 void
130  std::vector<double> &nodal_values) const
131 {
132  Assert (support_point_values.size() == this->generalized_support_points.size(),
133  ExcDimensionMismatch(support_point_values.size(), this->generalized_support_points.size()));
134  AssertDimension (support_point_values[0].size(), dim);
135  Assert (nodal_values.size() == this->dofs_per_cell,
136  ExcDimensionMismatch(nodal_values.size(),this->dofs_per_cell));
137 
138  // First do interpolation on faces. There, the component evaluated
139  // depends on the face direction and orientation.
140 
141  // The index of the first dof on this face or the cell
142  unsigned int dbase = 0;
143  // The index of the first generalized support point on this face or the cell
144  unsigned int pbase = 0;
145  for (unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
146  {
147  // Old version with no moments in 2D. See comment below in
148  // initialize_support_points()
149  if (test_values_face.size() == 0)
150  {
151  for (unsigned int i=0; i<this->dofs_per_face; ++i)
152  nodal_values[dbase+i] = support_point_values[pbase+i][GeometryInfo<dim>::unit_normal_direction[f]];
153  pbase += this->dofs_per_face;
154  }
155  else
156  {
157  for (unsigned int i=0; i<this->dofs_per_face; ++i)
158  {
159  double s = 0.;
160  for (unsigned int k=0; k<test_values_face.size(); ++k)
161  s += support_point_values[pbase+k][GeometryInfo<dim>::unit_normal_direction[f]] * test_values_face[k][i];
162  nodal_values[dbase+i] = s;
163  }
164  pbase += test_values_face.size();
165  }
166  dbase += this->dofs_per_face;
167  }
168 
169  AssertDimension (dbase, this->dofs_per_face * GeometryInfo<dim>::faces_per_cell);
170  AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
171 
172  // Done for BDM1
173  if (dbase == this->dofs_per_cell) return;
174 
175  // What's missing are the interior
176  // degrees of freedom. In each
177  // point, we take all components of
178  // the solution.
179  Assert ((this->dofs_per_cell - dbase) % dim == 0, ExcInternalError());
180 
181  for (unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
182  {
183  for (unsigned int i=0; i<test_values_cell[0].size(); ++i)
184  {
185  double s = 0.;
186  for (unsigned int k=0; k<test_values_cell.size(); ++k)
187  s += support_point_values[pbase+k][d] * test_values_cell[k][i];
188  nodal_values[dbase+i] = s;
189  }
190  }
191 
192  Assert (dbase == this->dofs_per_cell, ExcInternalError());
193 }
194 
195 
196 
197 template <int dim>
198 std::vector<unsigned int>
199 FE_BDM<dim>::get_dpo_vector (const unsigned int deg)
200 {
201  // compute the number of unknowns per cell interior/face/edge
202  //
203  // for the number of interior dofs, this is the number of
204  // polynomials up to degree deg-2 in dim dimensions.
205  //
206  // the element is face-based and we have as many degrees of freedom
207  // on the faces as there are polynomials of degree up to
208  // deg. Observe the odd convention of
209  // PolynomialSpace::compute_n_pols()!
210 
211  std::vector<unsigned int> dpo(dim+1, 0u);
212  dpo[dim] = (deg > 1 ?
214  0u);
215  dpo[dim-1] = PolynomialSpace<dim-1>::compute_n_pols(deg+1);
216 
217  return dpo;
218 }
219 
220 
221 
222 template <int dim>
223 std::vector<bool>
224 FE_BDM<dim>::get_ria_vector (const unsigned int deg)
225 {
226  if (dim==1)
227  {
228  Assert (false, ExcImpossibleInDim(1));
229  return std::vector<bool>();
230  }
231 
232  const unsigned int dofs_per_cell = PolynomialsBDM<dim>::compute_n_pols(deg);
233  const unsigned int dofs_per_face = PolynomialSpace<dim-1>::compute_n_pols(deg+1);
234 
235  Assert(GeometryInfo<dim>::faces_per_cell*dofs_per_face <= dofs_per_cell,
236  ExcInternalError());
237 
238  // all dofs need to be
239  // non-additive, since they have
240  // continuity requirements.
241  // however, the interior dofs are
242  // made additive
243  std::vector<bool> ret_val(dofs_per_cell,false);
244  for (unsigned int i=GeometryInfo<dim>::faces_per_cell*dofs_per_face;
245  i < dofs_per_cell; ++i)
246  ret_val[i] = true;
247 
248  return ret_val;
249 }
250 
251 
252 namespace internal
253 {
254  namespace FE_BDM
255  {
256  namespace
257  {
258  // This function sets up the values of the polynomials we want to
259  // take moments with in the quadrature points. In fact, we multiply
260  // thos by the weights, such that the sum of function values and
261  // test_values over quadrature points yields the interpolated degree
262  // of freedom.
263  template <int dim>
264  void
265  initialize_test_values (std::vector<std::vector<double> > &test_values,
266  const Quadrature<dim> &quadrature,
267  const unsigned int deg)
268  {
269  PolynomialsP<dim> poly(deg);
270  std::vector<Tensor<1,dim> > dummy1;
271  std::vector<Tensor<2,dim> > dummy2;
272  std::vector<Tensor<3,dim> > dummy3;
273  std::vector<Tensor<4,dim> > dummy4;
274 
275  test_values.resize(quadrature.size());
276 
277  for (unsigned int k=0; k<quadrature.size(); ++k)
278  {
279  test_values[k].resize(poly.n());
280  poly.compute(quadrature.point(k), test_values[k], dummy1, dummy2,
281  dummy3, dummy4);
282  for (unsigned int i=0; i < poly.n(); ++i)
283  {
284  test_values[k][i] *= quadrature.weight(k);
285  }
286  }
287  }
288 
289  // This specialization only serves to avoid error messages. Nothing
290  // useful can be computed in dimension zero and thus the vector
291  // length stays zero.
292  template <>
293  void
294  initialize_test_values (std::vector<std::vector<double> > &,
295  const Quadrature<0> &,
296  const unsigned int)
297  {}
298  }
299  }
300 }
301 
302 
303 template <int dim>
304 void
306 {
307  // Our support points are quadrature points on faces and inside the
308  // cell. First on the faces, we have to test polynomials of degree
309  // up to deg, which means we need dg+1 points in each direction. The
310  // fact that we do not have tensor product polynomials will be
311  // considered later. In 2D, we can use point values.
312  QGauss<dim-1> face_points (deg+1);
313 
314  // Copy the quadrature formula to the face points.
315  this->generalized_face_support_points.resize (face_points.size());
316  for (unsigned int k=0; k<face_points.size(); ++k)
317  this->generalized_face_support_points[k] = face_points.point(k);
318 
319  // In the interior, we only test with polynomials of degree up to
320  // deg-2, thus we use deg points. Note that deg>=1 and the lowest
321  // order element has no points in the cell, such that we have to
322  // distinguish this case.
323  QGauss<dim> cell_points(deg==1 ? 0 : deg);
324 
325  // Compute the size of the whole support point set
326  const unsigned int npoints
327  = cell_points.size() + GeometryInfo<dim>::faces_per_cell * face_points.size();
328 
329  this->generalized_support_points.resize (npoints);
330 
332  for (unsigned int k=0; k < face_points.size()*GeometryInfo<dim>::faces_per_cell; ++k)
333  this->generalized_support_points[k]
334  = faces.point(k+QProjector<dim>
335  ::DataSetDescriptor::face(0, true, false, false,
336  this->dofs_per_face));
337 
338  // Currently, for backward compatibility, we do not use moments, but
339  // point values on faces in 2D. In 3D, this is impossible, since the
340  // moments are only taken with respect to PolynomialsP.
341  if (dim>2)
342  internal::FE_BDM::initialize_test_values(test_values_face, face_points, deg);
343 
344  if (deg<=1) return;
345 
346  // Remember where interior points start
347  const unsigned int ibase = face_points.size()*GeometryInfo<dim>::faces_per_cell;
348  for (unsigned int k=0; k<cell_points.size(); ++k)
349  {
350  this->generalized_support_points[ibase+k] = cell_points.point(k);
351  }
352  // Finally, compute the values of
353  // the test functions in the
354  // interior quadrature points
355 
356  internal::FE_BDM::initialize_test_values(test_values_cell, cell_points, deg-2);
357 }
358 
359 
360 
361 /*-------------- Explicit Instantiations -------------------------------*/
362 #include "fe_bdm.inst"
363 
364 DEAL_II_NAMESPACE_CLOSE
size_type m() const
void initialize_support_points(const unsigned int bdm_degree)
Definition: fe_bdm.cc:305
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
FullMatrix< double > interface_constraints
Definition: fe.h:2401
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
Definition: fe_bdm.cc:119
void compute_face_embedding_matrices(const FiniteElement< dim, spacedim > &fe, FullMatrix< number >(&matrices)[GeometryInfo< dim >::max_children_per_face], const unsigned int face_coarse, const unsigned int face_fine, const double threshold=1.e-12)
const Point< dim > & point(const unsigned int i) 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)
Definition: fe_bdm.h:56
static std::vector< bool > get_ria_vector(const unsigned int degree)
Definition: fe_bdm.cc:224
FE_BDM(const unsigned int p)
Definition: fe_bdm.cc:38
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_bdm.cc:129
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
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 ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
static unsigned int compute_n_pols(const unsigned int n)
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual std::string get_name() const
Definition: fe_bdm.cc:97
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_bdm.cc:199
static unsigned int compute_n_pols(unsigned int degree)
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
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()