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