Reference documentation for deal.II version 8.5.1
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 
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 template <int dim>
37 FE_BDM<dim>::FE_BDM (const unsigned int deg)
38  :
39  FE_PolyTensor<PolynomialsBDM<dim>, dim> (
40  deg,
41  FiniteElementData<dim>(get_dpo_vector(deg),
42  dim, deg+1, FiniteElementData<dim>::Hdiv),
43  get_ria_vector (deg),
44  std::vector<ComponentMask>(PolynomialsBDM<dim>::compute_n_pols(deg),
45  std::vector<bool>(dim,true)))
46 {
47  Assert (dim >= 2, ExcImpossibleInDim(dim));
48  Assert (deg > 0, ExcMessage("Lowest order BDM element are degree 1, but you asked for degree 0"));
49 
50  const unsigned int n_dofs = this->dofs_per_cell;
51 
52  this->mapping_type = mapping_bdm;
53  // These must be done first, since
54  // they change the evaluation of
55  // basis functions
56 
57  // Set up the generalized support
58  // points
60 
61  // Now compute the inverse node matrix, generating the correct
62  // basis functions from the raw ones. For a discussion of what
63  // exactly happens here, see FETools::compute_node_matrix.
65  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
66  this->inverse_node_matrix.invert(M);
67  // From now on, the shape functions provided by FiniteElement::shape_value
68  // and similar functions will be the correct ones, not
69  // the raw shape functions from the polynomial space anymore.
70 
71  // Embedding errors become pretty large, so we just replace the
72  // regular threshold in both "computing_..." functions by 1.
74  FETools::compute_embedding_matrices (*this, this->prolongation, true, 1.);
75 
77  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
78  face_embeddings[i].reinit (this->dofs_per_face, this->dofs_per_face);
79  FETools::compute_face_embedding_matrices(*this, face_embeddings, 0, 0, 1.);
80  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
81  this->dofs_per_face);
82  unsigned int target_row=0;
83  for (unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
84  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
85  {
86  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
87  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
88  ++target_row;
89  }
90 }
91 
92 
93 
94 template <int dim>
95 std::string
97 {
98  // note that the
99  // FETools::get_fe_by_name
100  // function depends on the
101  // particular format of the string
102  // this function returns, so they
103  // have to be kept in synch
104 
105  // note that this->degree is the maximal
106  // polynomial degree and is thus one higher
107  // than the argument given to the
108  // constructor
109  std::ostringstream namebuf;
110  namebuf << "FE_BDM<" << dim << ">(" << this->degree-1 << ")";
111 
112  return namebuf.str();
113 }
114 
115 
116 template <int dim>
119 {
120  return new FE_BDM<dim>(*this);
121 }
122 
123 
124 
125 template <int dim>
126 void
129  std::vector<double> &nodal_values) const
130 {
131  Assert (support_point_values.size() == this->generalized_support_points.size(),
132  ExcDimensionMismatch(support_point_values.size(), this->generalized_support_points.size()));
133  AssertDimension (support_point_values[0].size(), dim);
134  Assert (nodal_values.size() == this->dofs_per_cell,
135  ExcDimensionMismatch(nodal_values.size(),this->dofs_per_cell));
136 
137  // First do interpolation on faces. There, the component evaluated
138  // depends on the face direction and orientation.
139 
140  // The index of the first dof on this face or the cell
141  unsigned int dbase = 0;
142  // The index of the first generalized support point on this face or the cell
143  unsigned int pbase = 0;
144  for (unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
145  {
146  // Old version with no moments in 2D. See comment below in
147  // initialize_support_points()
148  if (test_values_face.size() == 0)
149  {
150  for (unsigned int i=0; i<this->dofs_per_face; ++i)
151  nodal_values[dbase+i] = support_point_values[pbase+i][GeometryInfo<dim>::unit_normal_direction[f]];
152  pbase += this->dofs_per_face;
153  }
154  else
155  {
156  for (unsigned int i=0; i<this->dofs_per_face; ++i)
157  {
158  double s = 0.;
159  for (unsigned int k=0; k<test_values_face.size(); ++k)
160  s += support_point_values[pbase+k][GeometryInfo<dim>::unit_normal_direction[f]] * test_values_face[k][i];
161  nodal_values[dbase+i] = s;
162  }
163  pbase += test_values_face.size();
164  }
165  dbase += this->dofs_per_face;
166  }
167 
168  AssertDimension (dbase, this->dofs_per_face * GeometryInfo<dim>::faces_per_cell);
169  AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
170 
171  // Done for BDM1
172  if (dbase == this->dofs_per_cell) return;
173 
174  // What's missing are the interior
175  // degrees of freedom. In each
176  // point, we take all components of
177  // the solution.
178  Assert ((this->dofs_per_cell - dbase) % dim == 0, ExcInternalError());
179 
180  for (unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
181  {
182  for (unsigned int i=0; i<test_values_cell[0].size(); ++i)
183  {
184  double s = 0.;
185  for (unsigned int k=0; k<test_values_cell.size(); ++k)
186  s += support_point_values[pbase+k][d] * test_values_cell[k][i];
187  nodal_values[dbase+i] = s;
188  }
189  }
190 
191  Assert (dbase == this->dofs_per_cell, ExcInternalError());
192 }
193 
194 
195 
196 template <int dim>
197 void
199  std::vector<double> &,
200  const std::vector<double> &) const
201 {
202  Assert(false, ExcNotImplemented());
203 }
204 
205 
206 template <int dim>
207 void
209  std::vector<double> &,
210  const std::vector<Vector<double> > &,
211  unsigned int) const
212 {
213  Assert(false, ExcNotImplemented());
214 }
215 
216 
217 
218 template <int dim>
219 void
221  std::vector<double> &local_dofs,
222  const VectorSlice<const std::vector<std::vector<double> > > &values) const
223 {
224  AssertDimension (values.size(), dim);
225  Assert (values[0].size() == this->generalized_support_points.size(),
226  ExcDimensionMismatch(values.size(), this->generalized_support_points.size()));
227  Assert (local_dofs.size() == this->dofs_per_cell,
228  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
229 
230  // First do interpolation on faces. There, the component evaluated
231  // depends on the face direction and orientation.
232 
233  // The index of the first dof on this face or the cell
234  unsigned int dbase = 0;
235  // The index of the first generalized support point on this face or the cell
236  unsigned int pbase = 0;
237  for (unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
238  {
239  // Old version with no moments in 2D. See comment below in
240  // initialize_support_points()
241  if (test_values_face.size() == 0)
242  {
243  for (unsigned int i=0; i<this->dofs_per_face; ++i)
244  local_dofs[dbase+i] = values[GeometryInfo<dim>::unit_normal_direction[f]][pbase+i];
245  pbase += this->dofs_per_face;
246  }
247  else
248  {
249  for (unsigned int i=0; i<this->dofs_per_face; ++i)
250  {
251  double s = 0.;
252  for (unsigned int k=0; k<test_values_face.size(); ++k)
253  s += values[GeometryInfo<dim>::unit_normal_direction[f]][pbase+k] * test_values_face[k][i];
254  local_dofs[dbase+i] = s;
255  }
256  pbase += test_values_face.size();
257  }
258  dbase += this->dofs_per_face;
259  }
260 
261  AssertDimension (dbase, this->dofs_per_face * GeometryInfo<dim>::faces_per_cell);
262  AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
263 
264  // Done for BDM1
265  if (dbase == this->dofs_per_cell) return;
266 
267  // What's missing are the interior
268  // degrees of freedom. In each
269  // point, we take all components of
270  // the solution.
271  Assert ((this->dofs_per_cell - dbase) % dim == 0, ExcInternalError());
272 
273  for (unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
274  {
275  for (unsigned int i=0; i<test_values_cell[0].size(); ++i)
276  {
277  double s = 0.;
278  for (unsigned int k=0; k<test_values_cell.size(); ++k)
279  s += values[d][pbase+k] * test_values_cell[k][i];
280  local_dofs[dbase+i] = s;
281  }
282  }
283 
284  Assert (dbase == this->dofs_per_cell, ExcInternalError());
285 }
286 
287 
288 
289 
290 template <int dim>
291 std::vector<unsigned int>
292 FE_BDM<dim>::get_dpo_vector (const unsigned int deg)
293 {
294  // compute the number of unknowns per cell interior/face/edge
295  //
296  // for the number of interior dofs, this is the number of
297  // polynomials up to degree deg-2 in dim dimensions.
298  //
299  // the element is face-based and we have as many degrees of freedom
300  // on the faces as there are polynomials of degree up to
301  // deg. Observe the odd convention of
302  // PolynomialSpace::compute_n_pols()!
303 
304  std::vector<unsigned int> dpo(dim+1, 0u);
305  dpo[dim] = (deg > 1 ?
307  0u);
308  dpo[dim-1] = PolynomialSpace<dim-1>::compute_n_pols(deg+1);
309 
310  return dpo;
311 }
312 
313 
314 
315 template <int dim>
316 std::vector<bool>
317 FE_BDM<dim>::get_ria_vector (const unsigned int deg)
318 {
319  if (dim==1)
320  {
321  Assert (false, ExcImpossibleInDim(1));
322  return std::vector<bool>();
323  }
324 
325  const unsigned int dofs_per_cell = PolynomialsBDM<dim>::compute_n_pols(deg);
326  const unsigned int dofs_per_face = PolynomialSpace<dim-1>::compute_n_pols(deg+1);
327 
328  Assert(GeometryInfo<dim>::faces_per_cell*dofs_per_face <= dofs_per_cell,
329  ExcInternalError());
330 
331  // all dofs need to be
332  // non-additive, since they have
333  // continuity requirements.
334  // however, the interior dofs are
335  // made additive
336  std::vector<bool> ret_val(dofs_per_cell,false);
337  for (unsigned int i=GeometryInfo<dim>::faces_per_cell*dofs_per_face;
338  i < dofs_per_cell; ++i)
339  ret_val[i] = true;
340 
341  return ret_val;
342 }
343 
344 
345 namespace
346 {
347  // This function sets up the values of the polynomials we want to
348  // take moments with in the quadrature points. In fact, we multiply
349  // thos by the weights, such that the sum of function values and
350  // test_values over quadrature points yields the interpolated degree
351  // of freedom.
352  template <int dim>
353  void
354  initialize_test_values (std::vector<std::vector<double> > &test_values,
355  const Quadrature<dim> &quadrature,
356  const unsigned int deg)
357  {
358  PolynomialsP<dim> poly(deg);
359  std::vector<Tensor<1,dim> > dummy1;
360  std::vector<Tensor<2,dim> > dummy2;
361  std::vector<Tensor<3,dim> > dummy3;
362  std::vector<Tensor<4,dim> > dummy4;
363 
364  test_values.resize(quadrature.size());
365 
366  for (unsigned int k=0; k<quadrature.size(); ++k)
367  {
368  test_values[k].resize(poly.n());
369  poly.compute(quadrature.point(k), test_values[k], dummy1, dummy2,
370  dummy3, dummy4);
371  for (unsigned int i=0; i < poly.n(); ++i)
372  {
373  test_values[k][i] *= quadrature.weight(k);
374  }
375  }
376  }
377 
378  // This specialization only serves to avoid error messages. Nothing
379  // useful can be computed in dimension zero and thus the vector
380  // length stays zero.
381  template <>
382  void
383  initialize_test_values (std::vector<std::vector<double> > &,
384  const Quadrature<0> &,
385  const unsigned int)
386  {}
387 }
388 
389 
390 template <int dim>
391 void
393 {
394  // Our support points are quadrature points on faces and inside the
395  // cell. First on the faces, we have to test polynomials of degree
396  // up to deg, which means we need dg+1 points in each direction. The
397  // fact that we do not have tensor product polynomials will be
398  // considered later. In 2D, we can use point values.
399  QGauss<dim-1> face_points (deg+1);
400 
401  // Copy the quadrature formula to the face points.
402  this->generalized_face_support_points.resize (face_points.size());
403  for (unsigned int k=0; k<face_points.size(); ++k)
404  this->generalized_face_support_points[k] = face_points.point(k);
405 
406  // In the interior, we only test with polynomials of degree up to
407  // deg-2, thus we use deg points. Note that deg>=1 and the lowest
408  // order element has no points in the cell, such that we have to
409  // distinguish this case.
410  QGauss<dim> cell_points(deg==1 ? 0 : deg);
411 
412  // Compute the size of the whole support point set
413  const unsigned int npoints
414  = cell_points.size() + GeometryInfo<dim>::faces_per_cell * face_points.size();
415 
416  this->generalized_support_points.resize (npoints);
417 
419  for (unsigned int k=0; k < face_points.size()*GeometryInfo<dim>::faces_per_cell; ++k)
420  this->generalized_support_points[k]
421  = faces.point(k+QProjector<dim>
422  ::DataSetDescriptor::face(0, true, false, false,
423  this->dofs_per_face));
424 
425  // Currently, for backward compatibility, we do not use moments, but
426  // point values on faces in 2D. In 3D, this is impossible, since the
427  // moments are only taken with respect to PolynomialsP.
428  if (dim>2)
429  initialize_test_values(test_values_face, face_points, deg);
430 
431  if (deg<=1) return;
432 
433  // Remember where interior points start
434  const unsigned int ibase = face_points.size()*GeometryInfo<dim>::faces_per_cell;
435  for (unsigned int k=0; k<cell_points.size(); ++k)
436  {
437  this->generalized_support_points[ibase+k] = cell_points.point(k);
438  }
439  // Finally, compute the values of
440  // the test functions in the
441  // interior quadrature points
442 
443  initialize_test_values(test_values_cell, cell_points, deg-2);
444 }
445 
446 
447 
448 /*-------------- Explicit Instantiations -------------------------------*/
449 #include "fe_bdm.inst"
450 
451 DEAL_II_NAMESPACE_CLOSE
452 
size_type m() const
void initialize_support_points(const unsigned int bdm_degree)
Definition: fe_bdm.cc:392
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
FullMatrix< double > interface_constraints
Definition: fe.h:2211
virtual void convert_generalized_support_point_values_to_nodal_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe_bdm.cc:128
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)
Definition: fe_tools.cc:1920
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)
Definition: fe_tools.cc:1897
Definition: fe_bdm.h:56
static std::vector< bool > get_ria_vector(const unsigned int degree)
Definition: fe_bdm.cc:317
FE_BDM(const unsigned int p)
Definition: fe_bdm.cc:37
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:2199
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:313
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
virtual FiniteElement< dim > * clone() const
Definition: fe_bdm.cc:118
static unsigned int compute_n_pols(const unsigned int n)
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
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const 1
Definition: fe_bdm.cc:198
virtual std::string get_name() const
Definition: fe_bdm.cc:96
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_bdm.cc:292
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:277
static ::ExceptionBase & ExcNotImplemented()
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()