Reference documentation for deal.II version 8.4.2
fe_bdm.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/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  //Now compute the inverse node
61  //matrix, generating the correct
62  //basis functions from the raw
63  //ones.
64 
65  // We use an auxiliary matrix in
66  // this function. Therefore,
67  // inverse_node_matrix is still
68  // empty and shape_value_component
69  // returns the 'raw' shape values.
70  FullMatrix<double> M(n_dofs, n_dofs);
72 
73 // std::cout << std::endl;
74 // M.print_formatted(std::cout, 2, true);
75 
76  this->inverse_node_matrix.reinit(n_dofs, n_dofs);
77  this->inverse_node_matrix.invert(M);
78  // From now on, the shape functions
79  // will be the correct ones, not
80  // the raw shape functions anymore.
81 
82  // Embedding errors become pretty large, so we just replace the
83  // regular threshold in both "computing_..." functions by 1.
85  FETools::compute_embedding_matrices (*this, this->prolongation, true, 1.);
86 
88  for (unsigned int i=0; i<GeometryInfo<dim>::max_children_per_face; ++i)
89  face_embeddings[i].reinit (this->dofs_per_face, this->dofs_per_face);
90  FETools::compute_face_embedding_matrices(*this, face_embeddings, 0, 0, 1.);
91  this->interface_constraints.reinit((1<<(dim-1)) * this->dofs_per_face,
92  this->dofs_per_face);
93  unsigned int target_row=0;
94  for (unsigned int d=0; d<GeometryInfo<dim>::max_children_per_face; ++d)
95  for (unsigned int i=0; i<face_embeddings[d].m(); ++i)
96  {
97  for (unsigned int j=0; j<face_embeddings[d].n(); ++j)
98  this->interface_constraints(target_row,j) = face_embeddings[d](i,j);
99  ++target_row;
100  }
101 }
102 
103 
104 
105 template <int dim>
106 std::string
108 {
109  // note that the
110  // FETools::get_fe_from_name
111  // function depends on the
112  // particular format of the string
113  // this function returns, so they
114  // have to be kept in synch
115 
116  // note that this->degree is the maximal
117  // polynomial degree and is thus one higher
118  // than the argument given to the
119  // constructor
120  std::ostringstream namebuf;
121  namebuf << "FE_BDM<" << dim << ">(" << this->degree-1 << ")";
122 
123  return namebuf.str();
124 }
125 
126 
127 template <int dim>
130 {
131  return new FE_BDM<dim>(*this);
132 }
133 
134 
135 
136 template <int dim>
137 void
139  std::vector<double> &,
140  const std::vector<double> &) const
141 {
142  Assert(false, ExcNotImplemented());
143 }
144 
145 
146 template <int dim>
147 void
149  std::vector<double> &,
150  const std::vector<Vector<double> > &,
151  unsigned int) const
152 {
153  Assert(false, ExcNotImplemented());
154 }
155 
156 
157 
158 template <int dim>
159 void
161  std::vector<double> &local_dofs,
162  const VectorSlice<const std::vector<std::vector<double> > > &values) const
163 {
164  AssertDimension (values.size(), dim);
165  Assert (values[0].size() == this->generalized_support_points.size(),
166  ExcDimensionMismatch(values.size(), this->generalized_support_points.size()));
167  Assert (local_dofs.size() == this->dofs_per_cell,
168  ExcDimensionMismatch(local_dofs.size(),this->dofs_per_cell));
169 
170  // First do interpolation on faces. There, the component evaluated
171  // depends on the face direction and orientation.
172 
173  // The index of the first dof on this face or the cell
174  unsigned int dbase = 0;
175  // The index of the first generalized support point on this face or the cell
176  unsigned int pbase = 0;
177  for (unsigned int f = 0; f<GeometryInfo<dim>::faces_per_cell; ++f)
178  {
179  // Old version with no moments in 2D. See comment below in
180  // initialize_support_points()
181  if (test_values_face.size() == 0)
182  {
183  for (unsigned int i=0; i<this->dofs_per_face; ++i)
184  local_dofs[dbase+i] = values[GeometryInfo<dim>::unit_normal_direction[f]][pbase+i];
185  pbase += this->dofs_per_face;
186  }
187  else
188  {
189  for (unsigned int i=0; i<this->dofs_per_face; ++i)
190  {
191  double s = 0.;
192  for (unsigned int k=0; k<test_values_face.size(); ++k)
193  s += values[GeometryInfo<dim>::unit_normal_direction[f]][pbase+k] * test_values_face[k][i];
194  local_dofs[dbase+i] = s;
195  }
196  pbase += test_values_face.size();
197  }
198  dbase += this->dofs_per_face;
199  }
200 
202  AssertDimension (pbase, this->generalized_support_points.size() - test_values_cell.size());
203 
204  // Done for BDM1
205  if (dbase == this->dofs_per_cell) return;
206 
207  // What's missing are the interior
208  // degrees of freedom. In each
209  // point, we take all components of
210  // the solution.
211  Assert ((this->dofs_per_cell - dbase) % dim == 0, ExcInternalError());
212 
213  for (unsigned int d=0; d<dim; ++d, dbase += test_values_cell[0].size())
214  {
215  for (unsigned int i=0; i<test_values_cell[0].size(); ++i)
216  {
217  double s = 0.;
218  for (unsigned int k=0; k<test_values_cell.size(); ++k)
219  s += values[d][pbase+k] * test_values_cell[k][i];
220  local_dofs[dbase+i] = s;
221  }
222  }
223 
224  Assert (dbase == this->dofs_per_cell, ExcInternalError());
225 }
226 
227 
228 
229 
230 template <int dim>
231 std::vector<unsigned int>
232 FE_BDM<dim>::get_dpo_vector (const unsigned int deg)
233 {
234  // the element is face-based and we have as many degrees of freedom
235  // on the faces as there are polynomials of degree up to
236  // deg. Observe the odd convention of
237  // PolynomialSpace::compute_n_pols()!
239 
240  // and then there are interior dofs, namely the number of
241  // polynomials up to degree deg-2 in dim dimensions.
242  unsigned int interior_dofs = 0;
243  if (deg>1)
244  interior_dofs = dim * PolynomialSpace<dim>::compute_n_pols(deg-1);
245 
246  std::vector<unsigned int> dpo(dim+1);
247  dpo[dim-1] = dofs_per_face;
248  dpo[dim] = interior_dofs;
249 
250  return dpo;
251 }
252 
253 
254 
255 template <int dim>
256 std::vector<bool>
257 FE_BDM<dim>::get_ria_vector (const unsigned int deg)
258 {
259  if (dim==1)
260  {
261  Assert (false, ExcImpossibleInDim(1));
262  return std::vector<bool>();
263  }
264 
265  const unsigned int dofs_per_cell = PolynomialsBDM<dim>::compute_n_pols(deg);
266  const unsigned int dofs_per_face = PolynomialSpace<dim-1>::compute_n_pols(deg);
267 
268  Assert(GeometryInfo<dim>::faces_per_cell*dofs_per_face < dofs_per_cell,
269  ExcInternalError());
270 
271  // all dofs need to be
272  // non-additive, since they have
273  // continuity requirements.
274  // however, the interior dofs are
275  // made additive
276  std::vector<bool> ret_val(dofs_per_cell,false);
277  for (unsigned int i=GeometryInfo<dim>::faces_per_cell*dofs_per_face;
278  i < dofs_per_cell; ++i)
279  ret_val[i] = true;
280 
281  return ret_val;
282 }
283 
284 
285 namespace
286 {
287  // This function sets up the values of the polynomials we want to
288  // take moments with in the quadrature points. In fact, we multiply
289  // thos by the weights, such that the sum of function values and
290  // test_values over quadrature points yields the interpolated degree
291  // of freedom.
292  template <int dim>
293  void
294  initialize_test_values (std::vector<std::vector<double> > &test_values,
295  const Quadrature<dim> &quadrature,
296  const unsigned int deg)
297  {
298  PolynomialsP<dim> poly(deg);
299  std::vector<Tensor<1,dim> > dummy1;
300  std::vector<Tensor<2,dim> > dummy2;
301  std::vector<Tensor<3,dim> > dummy3;
302  std::vector<Tensor<4,dim> > dummy4;
303 
304  test_values.resize(quadrature.size());
305 
306  for (unsigned int k=0; k<quadrature.size(); ++k)
307  {
308  test_values[k].resize(poly.n());
309  poly.compute(quadrature.point(k), test_values[k], dummy1, dummy2,
310  dummy3, dummy4);
311  for (unsigned int i=0; i < poly.n(); ++i)
312  {
313  test_values[k][i] *= quadrature.weight(k);
314  }
315  }
316  }
317 
318  // This specialization only serves to avoid error messages. Nothing
319  // useful can be computed in dimension zero and thus the vector
320  // length stays zero.
321  template <>
322  void
323  initialize_test_values (std::vector<std::vector<double> > &,
324  const Quadrature<0> &,
325  const unsigned int)
326  {}
327 }
328 
329 
330 template <int dim>
331 void
333 {
334  // Our support points are quadrature points on faces and inside the
335  // cell. First on the faces, we have to test polynomials of degree
336  // up to deg, which means we need dg+1 points in each direction. The
337  // fact that we do not have tensor product polynomials will be
338  // considered later. In 2D, we can use point values.
339  QGauss<dim-1> face_points (deg+1);
340 
341  // Copy the quadrature formula to the face points.
342  this->generalized_face_support_points.resize (face_points.size());
343  for (unsigned int k=0; k<face_points.size(); ++k)
344  this->generalized_face_support_points[k] = face_points.point(k);
345 
346  // In the interior, we only test with polynomials of degree up to
347  // deg-2, thus we use deg points. Note that deg>=1 and the lowest
348  // order element has no points in the cell, such that we have to
349  // distinguish this case.
350  QGauss<dim> cell_points(deg==1 ? 0 : deg);
351 
352  // Compute the size of the whole support point set
353  const unsigned int npoints
354  = cell_points.size() + GeometryInfo<dim>::faces_per_cell * face_points.size();
355 
356  this->generalized_support_points.resize (npoints);
357 
359  for (unsigned int k=0; k < face_points.size()*GeometryInfo<dim>::faces_per_cell; ++k)
361  = faces.point(k+QProjector<dim>
362  ::DataSetDescriptor::face(0, true, false, false,
363  this->dofs_per_face));
364 
365  // Currently, for backward compatibility, we do not use moments, but
366  // point values on faces in 2D. In 3D, this is impossible, since the
367  // moments are only taken with respect to PolynomialsP.
368  if (dim>2)
369  initialize_test_values(test_values_face, face_points, deg);
370 
371  if (deg<=1) return;
372 
373  // Remember where interior points start
374  const unsigned int ibase = face_points.size()*GeometryInfo<dim>::faces_per_cell;
375  for (unsigned int k=0; k<cell_points.size(); ++k)
376  {
377  this->generalized_support_points[ibase+k] = cell_points.point(k);
378  }
379  // Finally, compute the values of
380  // the test functions in the
381  // interior quadrature points
382 
383  initialize_test_values(test_values_cell, cell_points, deg-2);
384 }
385 
386 
387 
388 /*-------------- Explicit Instantiations -------------------------------*/
389 #include "fe_bdm.inst"
390 
391 DEAL_II_NAMESPACE_CLOSE
392 
size_type m() const
std::vector< std::vector< double > > test_values_cell
Definition: fe_bdm.h:118
void initialize_support_points(const unsigned int bdm_degree)
Definition: fe_bdm.cc:332
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2157
FullMatrix< double > interface_constraints
Definition: fe.h:2132
::ExceptionBase & ExcMessage(std::string arg1)
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:914
const unsigned int degree
Definition: fe_base.h:299
const Point< dim > & point(const unsigned int i) 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
Definition: fe_bdm.h:56
static std::vector< bool > get_ria_vector(const unsigned int degree)
Definition: fe_bdm.cc:257
FE_BDM(const unsigned int p)
Definition: fe_bdm.cc:37
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2120
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:294
void compute(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const
unsigned int n() 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
virtual FiniteElement< dim > * clone() const
Definition: fe_bdm.cc:129
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:107
std::vector< std::vector< double > > test_values_face
Definition: fe_bdm.h:112
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_bdm.cc:232
static unsigned int compute_n_pols(unsigned int degree)
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 void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
Definition: fe_bdm.cc:138
void compute_node_matrix(FullMatrix< double > &M, const FiniteElement< dim, spacedim > &fe)
Definition: fe_tools.cc:633
double weight(const unsigned int i) const