Reference documentation for deal.II version 8.5.1
fe_bernstein.cc
1 // ---------------------------------------------------------------------
2 // @f$Id: fe_q.cc 30037 2013-07-18 16:55:40Z maier @f$
3 //
4 // Copyright (C) 2000 - 2016 by the deal.II authors
5 //
6 // This file is part of the deal.II library.
7 //
8 // The deal.II library is free software; you can use it, redistribute
9 // it, and/or modify it under the terms of the GNU Lesser General
10 // Public License as published by the Free Software Foundation; either
11 // version 2.1 of the License, or (at your option) any later version.
12 // The full text of the license can be found in the file LICENSE at
13 // the top level of the deal.II distribution.
14 //
15 // ---------------------------------------------------------------------
16 
17 
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/fe/fe_nothing.h>
22 #include <deal.II/fe/fe_q.h>
23 #include <deal.II/fe/fe_tools.h>
24 #include <deal.II/fe/fe_bernstein.h>
25 #include <deal.II/base/polynomials_bernstein.h>
26 
27 #include <vector>
28 #include <sstream>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 
34 template <int dim, int spacedim>
35 FE_Bernstein<dim,spacedim>::FE_Bernstein (const unsigned int degree)
36  :
37  FE_Q_Base<TensorProductPolynomials<dim>, dim, spacedim> (
38  this->renumber_bases(degree),
39  FiniteElementData<dim>(this->get_dpo_vector(degree),
40  1, degree,
41  FiniteElementData<dim>::H1),
42  std::vector<bool> (1, false))
43 {}
44 
45 
46 
47 template <int dim, int spacedim>
48 void
51  FullMatrix<double> &) const
52 {
53  // no interpolation possible. throw exception, as documentation says
54  typedef FiniteElement<dim,spacedim> FEE;
55  AssertThrow (false,
56  typename FEE::ExcInterpolationNotImplemented());
57 }
58 
59 
60 
61 template <int dim, int spacedim>
62 const FullMatrix<double> &
64  const RefinementCase<dim> &) const
65 {
66  typedef FiniteElement<dim,spacedim> FEE;
67  AssertThrow (false,
68  typename FEE::ExcProjectionVoid());
69  // return dummy, nothing will happen because the base class FE_Q_Base
70  // implements lazy evaluation of those matrices
71  return this->restriction[0][0];
72 }
73 
74 
75 
76 template <int dim, int spacedim>
77 const FullMatrix<double> &
79  const RefinementCase<dim> &) const
80 {
81  typedef FiniteElement<dim,spacedim> FEE;
82  AssertThrow (false,
83  typename FEE::ExcEmbeddingVoid());
84  // return dummy, nothing will happen because the base class FE_Q_Base
85  // implements lazy evaluation of those matrices
86  return this->prolongation[0][0];
87 }
88 
89 
90 
91 template <int dim, int spacedim>
92 void
95  FullMatrix<double> &interpolation_matrix) const
96 {
97  Assert (dim > 1, ExcImpossibleInDim(1));
98  get_subface_interpolation_matrix (source_fe, numbers::invalid_unsigned_int,
99  interpolation_matrix);
100 }
101 
102 
103 template <int dim, int spacedim>
104 void
107  const unsigned int subface,
108  FullMatrix<double> &interpolation_matrix) const
109 {
110  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
111  ExcDimensionMismatch (interpolation_matrix.m(),
112  x_source_fe.dofs_per_face));
113 
114  // see if source is a Bernstein element
115  if (const FE_Bernstein<dim,spacedim> *source_fe
116  = dynamic_cast<const FE_Bernstein<dim,spacedim> *>(&x_source_fe))
117  {
118  // have this test in here since a table of size 2x0 reports its size as
119  // 0x0
120  Assert (interpolation_matrix.n() == this->dofs_per_face,
121  ExcDimensionMismatch (interpolation_matrix.n(),
122  this->dofs_per_face));
123 
124  // Make sure that the element for which the DoFs should be constrained
125  // is the one with the higher polynomial degree. Actually the procedure
126  // will work also if this assertion is not satisfied. But the matrices
127  // produced in that case might lead to problems in the hp procedures,
128  // which use this method.
129  Assert (this->dofs_per_face <= source_fe->dofs_per_face,
130  (typename FiniteElement<dim,spacedim>::
131  ExcInterpolationNotImplemented ()));
132 
133  const Quadrature<dim-1>
134  quad_face_support(FE_Q<dim,spacedim>(QIterated<1>(QTrapez<1>(),source_fe->degree)).get_unit_face_support_points ());
135 
136  // Rule of thumb for FP accuracy, that can be expected for a given
137  // polynomial degree. This value is used to cut off values close to
138  // zero.
139  double eps = 2e-13 * std::max(this->degree, source_fe->degree) * (dim-1);
140 
141  // compute the interpolation matrix by simply taking the value at the
142  // support points.
143 //TODO: Verify that all faces are the same with respect to
144 // these support points. Furthermore, check if something has to
145 // be done for the face orientation flag in 3D.
146  const Quadrature<dim> subface_quadrature
147  = subface == numbers::invalid_unsigned_int
148  ?
149  QProjector<dim>::project_to_face (quad_face_support, 0)
150  :
151  QProjector<dim>::project_to_subface (quad_face_support, 0, subface);
152 
153  for (unsigned int i=0; i<source_fe->dofs_per_face; ++i)
154  {
155  const Point<dim> &p = subface_quadrature.point (i);
156  for (unsigned int j=0; j<this->dofs_per_face; ++j)
157  {
158  double matrix_entry = this->shape_value (this->face_to_cell_index(j, 0), p);
159 
160  // Correct the interpolated value. I.e. if it is close to 1 or
161  // 0, make it exactly 1 or 0. Unfortunately, this is required to
162  // avoid problems with higher order elements.
163  if (std::fabs (matrix_entry - 1.0) < eps)
164  matrix_entry = 1.0;
165  if (std::fabs (matrix_entry) < eps)
166  matrix_entry = 0.0;
167 
168  interpolation_matrix(i,j) = matrix_entry;
169  }
170  }
171 
172  // make sure that the row sum of each of the matrices is 1 at this
173  // point. this must be so since the shape functions sum up to 1
174  for (unsigned int j=0; j<source_fe->dofs_per_face; ++j)
175  {
176  double sum = 0.;
177 
178  for (unsigned int i=0; i<this->dofs_per_face; ++i)
179  sum += interpolation_matrix(j,i);
180 
181  Assert (std::fabs(sum-1) < eps, ExcInternalError());
182  }
183  }
184  else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe) != 0)
185  {
186  // nothing to do here, the FE_Nothing has no degrees of freedom anyway
187  }
188  else
189  AssertThrow (false,(typename FiniteElement<dim,spacedim>::
190  ExcInterpolationNotImplemented()));
191 }
192 
193 
194 
195 template <int dim, int spacedim>
196 bool
198 {
199  return true;
200 }
201 
202 
203 template <int dim, int spacedim>
204 std::vector<std::pair<unsigned int, unsigned int> >
206 {
207  // we can presently only compute these identities if both FEs are FE_Bernsteins
208  // or if the other one is an FE_Nothing. in the first case, there should be
209  // exactly one single DoF of each FE at a vertex, and they should have
210  // identical value
211  if (dynamic_cast<const FE_Bernstein<dim,spacedim>*>(&fe_other) != 0)
212  {
213  return
214  std::vector<std::pair<unsigned int, unsigned int> >
215  (1, std::make_pair (0U, 0U));
216  }
217  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
218  {
219  // the FE_Nothing has no degrees of freedom, so there are no
220  // equivalencies to be recorded
221  return std::vector<std::pair<unsigned int, unsigned int> > ();
222  }
223  else if (fe_other.dofs_per_face == 0)
224  {
225  // if the other element has no elements on faces at all,
226  // then it would be impossible to enforce any kind of
227  // continuity even if we knew exactly what kind of element
228  // we have -- simply because the other element declares
229  // that it is discontinuous because it has no DoFs on
230  // its faces. in that case, just state that we have no
231  // constraints to declare
232  return std::vector<std::pair<unsigned int, unsigned int> > ();
233  }
234  else
235  {
236  Assert (false, ExcNotImplemented());
237  return std::vector<std::pair<unsigned int, unsigned int> > ();
238  }
239 }
240 
241 
242 template <int dim, int spacedim>
243 std::vector<std::pair<unsigned int, unsigned int> >
245 {
246  // Since this fe is not interpolatory but on the vertices, we can
247  // not identify dofs on lines and on quads even if there are dofs
248  // on lines and on quads.
249  //
250  // we also have nothing to say about interpolation to other finite
251  // elements. consequently, we never have anything to say at all
252  return std::vector<std::pair<unsigned int, unsigned int> > ();
253 }
254 
255 
256 template <int dim, int spacedim>
257 std::vector<std::pair<unsigned int, unsigned int> >
259 {
260  // Since this fe is not interpolatory but on the vertices, we can
261  // not identify dofs on lines and on quads even if there are dofs
262  // on lines and on quads.
263  //
264  // we also have nothing to say about interpolation to other finite
265  // elements. consequently, we never have anything to say at all
266  return std::vector<std::pair<unsigned int, unsigned int> > ();
267 }
268 
269 
270 template <int dim, int spacedim>
273 {
274  if (const FE_Bernstein<dim,spacedim> *fe_b_other
275  = dynamic_cast<const FE_Bernstein<dim,spacedim>*>(&fe_other))
276  {
277  if (this->degree < fe_b_other->degree)
279  else if (this->degree == fe_b_other->degree)
281  else
283  }
284  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
285  {
286  if (fe_nothing->is_dominating())
287  {
289  }
290  else
291  {
292  // the FE_Nothing has no degrees of freedom and it is typically used in
293  // a context where we don't require any continuity along the interface
295  }
296  }
297 
298  Assert (false, ExcNotImplemented());
300 }
301 
302 
303 template <int dim, int spacedim>
304 std::string
306 {
307  // note that the FETools::get_fe_by_name function depends on the
308  // particular format of the string this function returns, so they have to be
309  // kept in synch
310 
311  std::ostringstream namebuf;
312  namebuf << "FE_Bernstein<" << dim << ">(" << this->degree << ")";
313  return namebuf.str();
314 }
315 
316 
317 template <int dim, int spacedim>
320 {
321  return new FE_Bernstein<dim,spacedim>(*this);
322 }
323 
324 
328 template <int dim, int spacedim>
329 std::vector<unsigned int>
331 {
332  AssertThrow(deg>0,ExcMessage("FE_Bernstein needs to be of degree > 0."));
333  std::vector<unsigned int> dpo(dim+1, 1U);
334  for (unsigned int i=1; i<dpo.size(); ++i)
335  dpo[i]=dpo[i-1]*(deg-1);
336  return dpo;
337 }
338 
339 
340 template <int dim, int spacedim>
343 {
344  TensorProductPolynomials<dim> tpp(::generate_complete_bernstein_basis<double>(deg));
345  std::vector<unsigned int> renumber(Utilities::fixed_power<dim>(deg+1));
346  const FiniteElementData<dim> fe(this->get_dpo_vector(deg),1,
347  deg);
349  tpp.set_numbering(renumber);
350  return tpp;
351 }
352 
353 
354 // explicit instantiations
355 #include "fe_bernstein.inst"
356 
357 DEAL_II_NAMESPACE_CLOSE
size_type m() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_bernstein.cc:78
static const unsigned int invalid_unsigned_int
Definition: types.h:170
FE_Bernstein(const unsigned int p)
Definition: fe_bernstein.cc:35
const std::vector< Point< dim-1 > > & get_unit_face_support_points() const
Definition: fe.cc:1017
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_bernstein.cc:50
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_bernstein.cc:63
virtual std::string get_name() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
void hierarchic_to_lexicographic_numbering(unsigned int degree, std::vector< unsigned int > &h2l)
Definition: fe_tools.cc:2928
void set_numbering(const std::vector< unsigned int > &renumber)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual bool hp_constraints_are_implemented() const
static ::ExceptionBase & ExcMessage(std::string arg1)
Definition: fe_q.h:552
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
size_type n() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
#define Assert(cond, exc)
Definition: exceptions.h:313
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)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_bernstein.cc:94
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
TensorProductPolynomials< dim > renumber_bases(const unsigned int degree)
virtual FiniteElement< dim, spacedim > * clone() const
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
const unsigned int dofs_per_face
Definition: fe_base.h:288
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcInternalError()