Reference documentation for deal.II version 9.0.0
fe_q_bubbles.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2012 - 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/quadrature.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/template_constraints.h>
20 #include <deal.II/fe/fe_q_bubbles.h>
21 #include <deal.II/fe/fe_nothing.h>
22 #include <deal.II/fe/fe_tools.h>
23 #include <deal.II/fe/fe_values.h>
24 #include <deal.II/fe/mapping_q1.h>
25 #include <deal.II/base/quadrature_lib.h>
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/dofs/dof_handler.h>
28 #include <deal.II/grid/tria.h>
29 
30 #include <deal.II/grid/grid_generator.h>
31 
32 
33 #include <vector>
34 #include <sstream>
35 #include <memory>
36 #include <deal.II/base/std_cxx14/memory.h>
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 namespace internal
42 {
43  namespace FE_Q_Bubbles
44  {
45  namespace
46  {
47  template <int dim, int spacedim>
48  inline
49  void
50  compute_embedding_matrices(const ::FE_Q_Bubbles<dim, spacedim> &fe,
51  std::vector<std::vector<FullMatrix<double> > > &matrices,
52  const bool isotropic_only)
53  {
54  const unsigned int dpc = fe.dofs_per_cell;
55  const unsigned int degree = fe.degree;
56 
57  // Initialize quadrature formula on fine cells
58  std::unique_ptr<Quadrature<dim> > q_fine;
59  Quadrature<1> q_dummy(std::vector<Point<1> >(1), std::vector<double> (1,1.));
60  switch (dim)
61  {
62  case 1:
63  if (spacedim==1)
64  q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
65  else if (spacedim==2)
66  q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy);
67  else
68  q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy, q_dummy);
69  break;
70  case 2:
71  if (spacedim==2)
72  q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
73  else
74  q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), QGauss<1>(degree+1), q_dummy);
75  break;
76  case 3:
77  q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
78  break;
79  default:
80  Assert(false, ExcInternalError());
81  }
82 
83  Assert(q_fine.get() != nullptr, ExcInternalError());
84  const unsigned int nq = q_fine->size();
85 
86  // loop over all possible refinement cases
87  unsigned int ref_case = (isotropic_only)
90  for (; ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
91  {
92  const unsigned int nc
94 
95  for (unsigned int i=0; i<nc; ++i)
96  {
97  Assert(matrices[ref_case-1][i].n() == dpc,
98  ExcDimensionMismatch(matrices[ref_case-1][i].n(),dpc));
99  Assert(matrices[ref_case-1][i].m() == dpc,
100  ExcDimensionMismatch(matrices[ref_case-1][i].m(),dpc));
101  }
102 
103  // create a respective refinement on the triangulation
105  GridGenerator::hyper_cube (tr, 0, 1);
106  tr.begin_active()->set_refine_flag(RefinementCase<dim>(ref_case));
108 
110  dh.distribute_dofs(fe);
111 
116 
117  const unsigned int n_dofs = dh.n_dofs();
118 
119  FullMatrix<double> fine_mass(n_dofs);
120  FullMatrix<double> coarse_rhs_matrix(n_dofs, dpc);
121 
122  std::vector<std::vector<types::global_dof_index> > child_ldi
123  (nc, std::vector<types::global_dof_index>(fe.dofs_per_cell));
124 
125  //now create the mass matrix and all the right_hand sides
126  unsigned int child_no = 0;
127  typename ::DoFHandler<dim>::active_cell_iterator cell
128  = dh.begin_active();
129  for (; cell!=dh.end(); ++cell, ++child_no)
130  {
131  fine.reinit(cell);
132  cell->get_dof_indices(child_ldi[child_no]);
133 
134  for (unsigned int q=0; q<nq; ++q)
135  for (unsigned int i=0; i<dpc; ++i)
136  for (unsigned int j=0; j<dpc; ++j)
137  {
138  const unsigned int gdi=child_ldi[child_no][i];
139  const unsigned int gdj=child_ldi[child_no][j];
140  fine_mass(gdi, gdj)+=fine.shape_value(i,q)
141  *fine.shape_value(j,q)
142  *fine.JxW(q);
143  Point<dim> quad_tmp;
144  for (unsigned int k=0; k<dim; ++k)
145  quad_tmp(k) = fine.quadrature_point(q)(k);
146  coarse_rhs_matrix(gdi, j)
147  +=fine.shape_value(i,q)
148  *fe.shape_value(j, quad_tmp)
149  *fine.JxW(q);
150  }
151  }
152 
153  //now solve for all right-hand sides simultaneously
154  ::FullMatrix<double> solution (n_dofs, dpc);
155  fine_mass.gauss_jordan();
156  fine_mass.mmult(solution, coarse_rhs_matrix);
157 
158  //and distribute to the fine cell matrices
159  for (unsigned int child_no=0; child_no<nc; ++child_no)
160  for (unsigned int i=0; i<dpc; ++i)
161  for (unsigned int j=0; j<dpc; ++j)
162  {
163  const unsigned int gdi=child_ldi[child_no][i];
164  //remove small entries
165  if (std::fabs(solution(gdi, j)) > 1.e-12)
166  matrices[ref_case-1][child_no](i,j)=solution(gdi, j);
167  }
168  }
169  }
170  }
171  }
172 }
173 
174 
175 template <int dim, int spacedim>
176 FE_Q_Bubbles<dim,spacedim>::FE_Q_Bubbles (const unsigned int q_degree)
177  :
178  FE_Q_Base<TensorProductPolynomialsBubbles<dim>, dim, spacedim> (
179  TensorProductPolynomialsBubbles<dim>(Polynomials::generate_complete_Lagrange_basis(QGaussLobatto<1>(q_degree+1).get_points())),
180  FiniteElementData<dim>(get_dpo_vector(q_degree),
181  1, q_degree+1,
182  FiniteElementData<dim>::H1),
183  get_riaf_vector(q_degree)),
184  n_bubbles((q_degree<=1)?1:dim)
185 {
186  Assert (q_degree > 0,
187  ExcMessage ("This element can only be used for polynomial degrees "
188  "greater than zero"));
189 
190  this->initialize(QGaussLobatto<1>(q_degree+1).get_points());
191 
192  // adjust unit support point for discontinuous node
193  Point<dim> point;
194  for (unsigned int d=0; d<dim; ++d)
195  point[d] = 0.5;
196  for (unsigned int i=0; i<n_bubbles; ++i)
197  this->unit_support_points.push_back(point);
199 
201  if (dim == spacedim)
202  {
203  internal::FE_Q_Bubbles::compute_embedding_matrices
204  (*this, this->prolongation, false);
205  // Fill restriction matrices with L2-projection
207  }
208 
209 }
210 
211 
212 
213 template <int dim, int spacedim>
215  :
216  FE_Q_Base<TensorProductPolynomialsBubbles<dim>, dim, spacedim> (
217  TensorProductPolynomialsBubbles<dim>(Polynomials::generate_complete_Lagrange_basis(points.get_points())),
218  FiniteElementData<dim>(get_dpo_vector(points.size()-1),
219  1, points.size(),
220  FiniteElementData<dim>::H1),
221  get_riaf_vector(points.size()-1)),
222  n_bubbles((points.size()-1<=1)?1:dim)
223 {
224  Assert (points.size() > 1,
225  ExcMessage ("This element can only be used for polynomial degrees "
226  "at least one"));
227 
228  this->initialize(points.get_points());
229 
230  // adjust unit support point for discontinuous node
231  Point<dim> point;
232  for (unsigned int d=0; d<dim; ++d)
233  point[d] = 0.5;
234  for (unsigned int i=0; i< n_bubbles; ++i)
235  this->unit_support_points.push_back(point);
237 
239  if (dim == spacedim)
240  {
241  internal::FE_Q_Bubbles::compute_embedding_matrices
242  (*this, this->prolongation, false);
243  // Fill restriction matrices with L2-projection
245  }
246 }
247 
248 
249 
250 template <int dim, int spacedim>
251 std::string
253 {
254  // note that the FETools::get_fe_by_name function depends on the
255  // particular format of the string this function returns, so they have to be
256  // kept in synch
257 
258  std::ostringstream namebuf;
259  bool type = true;
260  const unsigned int n_points = this->degree;
261  std::vector<double> points(n_points);
262  const unsigned int dofs_per_cell = this->dofs_per_cell;
263  const std::vector<Point<dim> > &unit_support_points = this->unit_support_points;
264  unsigned int index = 0;
265 
266  // Decode the support points in one coordinate direction.
267  for (unsigned int j=0; j<dofs_per_cell; j++)
268  {
269  if ((dim>1) ? (unit_support_points[j](1)==0 &&
270  ((dim>2) ? unit_support_points[j](2)==0: true)) : true)
271  {
272  if (index == 0)
273  points[index] = unit_support_points[j](0);
274  else if (index == 1)
275  points[n_points-1] = unit_support_points[j](0);
276  else
277  points[index-1] = unit_support_points[j](0);
278 
279  index++;
280  }
281  }
282  // Do not consider the discontinuous node for dimension 1
283  Assert (index == n_points || (dim==1 && index == n_points+n_bubbles),
284  ExcMessage ("Could not decode support points in one coordinate direction."));
285 
286  // Check whether the support points are equidistant.
287  for (unsigned int j=0; j<n_points; j++)
288  if (std::fabs(points[j] - (double)j/(this->degree-1)) > 1e-15)
289  {
290  type = false;
291  break;
292  }
293 
294  if (type == true)
295  {
296  if (this->degree > 3)
297  namebuf << "FE_Q_Bubbles<"
298  << Utilities::dim_string(dim,spacedim)
299  << ">(QIterated(QTrapez()," << this->degree-1 << "))";
300  else
301  namebuf << "FE_Q_Bubbles<"
302  << Utilities::dim_string(dim,spacedim)
303  << ">(" << this->degree-1 << ")";
304  }
305  else
306  {
307 
308  // Check whether the support points come from QGaussLobatto.
309  const QGaussLobatto<1> points_gl(n_points);
310  type = true;
311  for (unsigned int j=0; j<n_points; j++)
312  if (points[j] != points_gl.point(j)(0))
313  {
314  type = false;
315  break;
316  }
317  if (type == true)
318  namebuf << "FE_Q_Bubbles<" << dim << ">(" << this->degree-1 << ")";
319  else
320  namebuf << "FE_Q_Bubbles<" << dim << ">(QUnknownNodes(" << this->degree << "))";
321  }
322  return namebuf.str();
323 }
324 
325 
326 
327 template <int dim, int spacedim>
328 std::unique_ptr<FiniteElement<dim,spacedim> >
330 {
331  return std_cxx14::make_unique<FE_Q_Bubbles<dim,spacedim>>(*this);
332 }
333 
334 
335 
336 template <int dim, int spacedim>
337 void
340  std::vector<double> &nodal_values) const
341 {
342  Assert (support_point_values.size() == this->unit_support_points.size(),
343  ExcDimensionMismatch(support_point_values.size(),
344  this->unit_support_points.size()));
345  Assert (nodal_values.size() == this->dofs_per_cell,
346  ExcDimensionMismatch(nodal_values.size(),this->dofs_per_cell));
347  Assert (support_point_values[0].size() == this->n_components(),
348  ExcDimensionMismatch(support_point_values[0].size(),this->n_components()));
349 
350  for (unsigned int i=0; i<this->dofs_per_cell-1; ++i)
351  {
352  const std::pair<unsigned int, unsigned int> index
353  = this->system_to_component_index(i);
354  nodal_values[i] = support_point_values[i](index.first);
355  }
356 
357  // We don't use the bubble functions for local interpolation
358  for (unsigned int i = 0; i<n_bubbles; ++i)
359  nodal_values[nodal_values.size()-i-1] = 0.;
360 }
361 
362 
363 
364 template <int dim, int spacedim>
365 void
368  FullMatrix<double> &interpolation_matrix) const
369 {
370  // We don't know how to do this properly, yet.
371  // However, for SolutionTransfer to work we need to provide an implementation
372  // for the case that the x_source_fe is identical to this FE
373  typedef FE_Q_Bubbles<dim,spacedim> FEQBUBBLES;
374 
375  AssertThrow ((x_source_fe.get_name().find ("FE_Q_Bubbles<") == 0)
376  ||
377  (dynamic_cast<const FEQBUBBLES *>(&x_source_fe) != nullptr)
378  ,
379  (typename FiniteElement<dim,spacedim>::
380  ExcInterpolationNotImplemented()));
381  Assert (interpolation_matrix.m() == this->dofs_per_cell,
382  ExcDimensionMismatch (interpolation_matrix.m(),
383  this->dofs_per_cell));
384  Assert (interpolation_matrix.n() == x_source_fe.dofs_per_cell,
385  ExcDimensionMismatch (interpolation_matrix.m(),
386  x_source_fe.dofs_per_cell));
387 
388  //Provide a short cut in case we are just inquiring the identity
389  auto casted_fe = dynamic_cast<const FEQBUBBLES *>(&x_source_fe);
390  if (casted_fe != nullptr && casted_fe->degree == this->degree)
391  for (unsigned int i=0; i<interpolation_matrix.m(); ++i)
392  interpolation_matrix.set(i,i,1.);
393  //else we need to do more...
394  else
395  Assert(false, (typename FiniteElement<dim,spacedim>::
396  ExcInterpolationNotImplemented()));
397 }
398 
399 
400 
401 template <int dim, int spacedim>
402 std::vector<bool>
404 {
405  unsigned int n_cont_dofs = Utilities::fixed_power<dim>(q_deg+1);
406  const unsigned int n_bubbles = (q_deg<=1?1:dim);
407  return std::vector<bool> (n_cont_dofs+n_bubbles,true);
408 }
409 
410 
411 
412 template <int dim, int spacedim>
413 std::vector<unsigned int>
415 {
416  std::vector<unsigned int> dpo(dim+1, 1U);
417  for (unsigned int i=1; i<dpo.size(); ++i)
418  dpo[i]=dpo[i-1]*(q_deg-1);
419 
420  dpo[dim]+=(q_deg<=1?1:dim);//all the bubble functions are discontinuous
421  return dpo;
422 }
423 
424 
425 
426 template <int dim, int spacedim>
427 bool
429 (const unsigned int shape_index,
430  const unsigned int face_index) const
431 {
432  // discontinuous functions have no support on faces
433  if (shape_index >= this->n_dofs_per_cell()-n_bubbles)
434  return false;
435  else
436  return FE_Q_Base<TensorProductPolynomialsBubbles<dim>,dim,spacedim>::has_support_on_face(shape_index, face_index);
437 }
438 
439 
440 
441 template <int dim, int spacedim>
442 const FullMatrix<double> &
444 (const unsigned int child,
445  const RefinementCase<dim> &refinement_case) const
446 {
449  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
450  ExcMessage("Prolongation matrices are only available for refined cells!"));
451  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
452  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
453 
454  Assert (this->prolongation[refinement_case-1][child].n() != 0,
455  ExcMessage("This prolongation matrix has not been computed yet!"));
456  // finally return the matrix
457  return this->prolongation[refinement_case-1][child];
458 }
459 
460 
461 
462 template <int dim, int spacedim>
463 const FullMatrix<double> &
465 (const unsigned int child,
466  const RefinementCase<dim> &refinement_case) const
467 {
470  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
471  ExcMessage("Restriction matrices are only available for refined cells!"));
472  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
473  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
474 
475  Assert(this->restriction[refinement_case-1][child].n() != 0,
476  ExcMessage("This restriction matrix has not been computed yet!"));
477 
478  //finally return the matrix
479  return this->restriction[refinement_case-1][child];
480 }
481 
482 // explicit instantiations
483 #include "fe_q_bubbles.inst"
484 
485 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
const unsigned int n_bubbles
Definition: fe_q_bubbles.h:168
Shape function values.
size_type m() const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
const std::vector< Point< dim > > & get_points() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
const Point< dim > & point(const unsigned int i) const
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10508
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)
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
Definition: point.h:104
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2413
void set(const size_type i, const size_type j, const number value)
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
virtual void execute_coarsening_and_refinement()
Definition: tria.cc:11739
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
#define Assert(cond, exc)
Definition: exceptions.h:1142
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::string get_name() const =0
virtual std::string get_name() const
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:176
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
void initialize(const std::vector< Point< 1 > > &support_points_1d)
Definition: fe_q_base.cc:434
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static std::vector< bool > get_riaf_vector(const unsigned int degree)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
Definition: fe.h:33
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
FE_Q_Bubbles(const unsigned int p)
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
static ::ExceptionBase & ExcInternalError()