Reference documentation for deal.II version 9.0.0
fe_dgp_monomial.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 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/fe/fe_dgp_monomial.h>
18 #include <deal.II/fe/fe_tools.h>
19 
20 #include <sstream>
21 #include <deal.II/base/std_cxx14/memory.h>
22 
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 namespace internal
27 {
28  namespace FE_DGPMonomial
29  {
30  namespace
31  {
32  // storage of hand-chosen support
33  // points
34  //
35  // For dim=2, dofs_per_cell of
36  // FE_DGPMonomial(k) is given by
37  // 0.5(k+1)(k+2), i.e.
38  //
39  // k 0 1 2 3 4 5 6 7
40  // dofs 1 3 6 10 15 21 28 36
41  //
42  // indirect access of unit points:
43  // the points for degree k are
44  // located at
45  //
46  // points[start_index[k]..start_index[k+1]-1]
47  const unsigned int start_index2d[6]= {0,1,4,10,20,35};
48  const double points2d[35][2]=
49  {
50  {0,0},
51  {0,0},{1,0},{0,1},
52  {0,0},{1,0},{0,1},{1,1},{0.5,0},{0,0.5},
53  {0,0},{1,0},{0,1},{1,1},{1./3.,0},{2./3.,0},{0,1./3.},{0,2./3.},{0.5,1},{1,0.5},
54  {0,0},{1,0},{0,1},{1,1},{0.25,0},{0.5,0},{0.75,0},{0,0.25},{0,0.5},{0,0.75},{1./3.,1},{2./3.,1},{1,1./3.},{1,2./3.},{0.5,0.5}
55  };
56 
57  // For dim=3, dofs_per_cell of
58  // FE_DGPMonomial(k) is given by
59  // 1./6.(k+1)(k+2)(k+3), i.e.
60  //
61  // k 0 1 2 3 4 5 6 7
62  // dofs 1 4 10 20 35 56 84 120
63  const unsigned int start_index3d[6]= {0,1,5,15/*,35*/};
64  const double points3d[35][3]=
65  {
66  {0,0,0},
67  {0,0,0},{1,0,0},{0,1,0},{0,0,1},
68  {0,0,0},{1,0,0},{0,1,0},{0,0,1},{0.5,0,0},{0,0.5,0},{0,0,0.5},{1,1,0},{1,0,1},{0,1,1}
69  };
70 
71 
72  template <int dim>
73  void generate_unit_points (const unsigned int,
74  std::vector<Point<dim> > &);
75 
76  template <>
77  void generate_unit_points (const unsigned int k,
78  std::vector<Point<1> > &p)
79  {
80  Assert(p.size()==k+1, ExcDimensionMismatch(p.size(), k+1));
81  const double h = 1./k;
82  for (unsigned int i=0; i<p.size(); ++i)
83  p[i](0)=i*h;
84  }
85 
86  template <>
87  void generate_unit_points (const unsigned int k,
88  std::vector<Point<2> > &p)
89  {
90  Assert(k<=4, ExcNotImplemented());
91  Assert(p.size()==start_index2d[k+1]-start_index2d[k], ExcInternalError());
92  for (unsigned int i=0; i<p.size(); ++i)
93  {
94  p[i](0)=points2d[start_index2d[k]+i][0];
95  p[i](1)=points2d[start_index2d[k]+i][1];
96  }
97  }
98 
99  template <>
100  void generate_unit_points (const unsigned int k,
101  std::vector<Point<3> > &p)
102  {
103  Assert(k<=2, ExcNotImplemented());
104  Assert(p.size()==start_index3d[k+1]-start_index3d[k], ExcInternalError());
105  for (unsigned int i=0; i<p.size(); ++i)
106  {
107  p[i](0)=points3d[start_index3d[k]+i][0];
108  p[i](1)=points3d[start_index3d[k]+i][1];
109  p[i](2)=points3d[start_index3d[k]+i][2];
110  }
111  }
112  }
113  }
114 }
115 
116 
117 
118 template <int dim>
119 FE_DGPMonomial<dim>::FE_DGPMonomial (const unsigned int degree)
120  :
121  FE_Poly<PolynomialsP<dim>, dim> (
122  PolynomialsP<dim>(degree),
123  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
124  std::vector<bool>(FiniteElementData<dim>(get_dpo_vector(degree), 1, degree).dofs_per_cell,true),
125  std::vector<ComponentMask>(FiniteElementData<dim>(
126  get_dpo_vector(degree), 1, degree).dofs_per_cell, std::vector<bool>(1,true)))
127 {
128  Assert(this->poly_space.n()==this->dofs_per_cell, ExcInternalError());
129  Assert(this->poly_space.degree()==this->degree, ExcInternalError());
130 
131  // DG doesn't have constraints, so
132  // leave them empty
133 
134  // Reinit the vectors of
135  // restriction and prolongation
136  // matrices to the right sizes
138  // Fill prolongation matrices with embedding operators
140  // Fill restriction matrices with L2-projection
142 }
143 
144 
145 
146 template <int dim>
147 std::string
149 {
150  // note that the
151  // FETools::get_fe_by_name
152  // function depends on the
153  // particular format of the string
154  // this function returns, so they
155  // have to be kept in synch
156 
157  std::ostringstream namebuf;
158  namebuf << "FE_DGPMonomial<" << dim << ">(" << this->degree << ")";
159 
160  return namebuf.str();
161 }
162 
163 
164 
165 template <int dim>
166 std::unique_ptr<FiniteElement<dim,dim> >
168 {
169  return std_cxx14::make_unique<FE_DGPMonomial<dim>>(*this);
170 }
171 
172 
173 
174 //TODO: Remove this function and use the one in FETools, if needed
175 template <int dim>
176 void
179  FullMatrix<double> &interpolation_matrix) const
180 {
181  const FE_DGPMonomial<dim> *source_dgp_monomial
182  = dynamic_cast<const FE_DGPMonomial<dim> *>(&source_fe);
183 
184  if (source_dgp_monomial)
185  {
186  // ok, source_fe is a DGP_Monomial
187  // element. Then, the interpolation
188  // matrix is simple
189  const unsigned int m=interpolation_matrix.m();
190  const unsigned int n=interpolation_matrix.n();
191  (void)m;
192  (void)n;
193  Assert (m == this->dofs_per_cell, ExcDimensionMismatch (m, this->dofs_per_cell));
194  Assert (n == source_dgp_monomial->dofs_per_cell,
195  ExcDimensionMismatch (n, source_dgp_monomial->dofs_per_cell));
196 
197  const unsigned int min_mn=
198  interpolation_matrix.m()<interpolation_matrix.n() ?
199  interpolation_matrix.m() : interpolation_matrix.n();
200 
201  for (unsigned int i=0; i<min_mn; ++i)
202  interpolation_matrix(i,i)=1.;
203  }
204  else
205  {
206  std::vector<Point<dim> > unit_points(this->dofs_per_cell);
207  internal::FE_DGPMonomial::generate_unit_points(this->degree, unit_points);
208 
209  FullMatrix<double> source_fe_matrix(unit_points.size(), source_fe.dofs_per_cell);
210  for (unsigned int j=0; j<source_fe.dofs_per_cell; ++j)
211  for (unsigned int k=0; k<unit_points.size(); ++k)
212  source_fe_matrix(k,j)=source_fe.shape_value(j, unit_points[k]);
213 
214  FullMatrix<double> this_matrix(this->dofs_per_cell, this->dofs_per_cell);
215  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
216  for (unsigned int k=0; k<unit_points.size(); ++k)
217  this_matrix(k,j)=this->poly_space.compute_value (j, unit_points[k]);
218 
219  this_matrix.gauss_jordan();
220 
221  this_matrix.mmult(interpolation_matrix, source_fe_matrix);
222  }
223 }
224 
225 
226 
227 template <int dim>
228 void
230 {
231  Assert(false, ExcNotImplemented());
232 }
233 
234 
235 //---------------------------------------------------------------------------
236 // Auxiliary functions
237 //---------------------------------------------------------------------------
238 
239 
240 template <int dim>
241 std::vector<unsigned int>
242 FE_DGPMonomial<dim>::get_dpo_vector (const unsigned int deg)
243 {
244  std::vector<unsigned int> dpo(dim+1, 0U);
245  dpo[dim] = deg+1;
246  for (unsigned int i=1; i<dim; ++i)
247  {
248  dpo[dim] *= deg+1+i;
249  dpo[dim] /= i+1;
250  }
251  return dpo;
252 }
253 
254 
255 template <int dim>
256 void
259  FullMatrix<double> &interpolation_matrix) const
260 {
261  // this is only implemented, if the source
262  // FE is also a DGPMonomial element. in that case,
263  // both elements have no dofs on their
264  // faces and the face interpolation matrix
265  // is necessarily empty -- i.e. there isn't
266  // much we need to do here.
267  (void)interpolation_matrix;
268  AssertThrow ((x_source_fe.get_name().find ("FE_DGPMonomial<") == 0)
269  ||
270  (dynamic_cast<const FE_DGPMonomial<dim>*>(&x_source_fe) != nullptr),
271  typename FiniteElement<dim>::
272  ExcInterpolationNotImplemented());
273 
274  Assert (interpolation_matrix.m() == 0,
275  ExcDimensionMismatch (interpolation_matrix.m(),
276  0));
277  Assert (interpolation_matrix.n() == 0,
278  ExcDimensionMismatch (interpolation_matrix.n(),
279  0));
280 }
281 
282 
283 
284 template <int dim>
285 void
288  const unsigned int,
289  FullMatrix<double> &interpolation_matrix) const
290 {
291  // this is only implemented, if the source
292  // FE is also a DGPMonomial element. in that case,
293  // both elements have no dofs on their
294  // faces and the face interpolation matrix
295  // is necessarily empty -- i.e. there isn't
296  // much we need to do here.
297  (void)interpolation_matrix;
298  AssertThrow ((x_source_fe.get_name().find ("FE_DGPMonomial<") == 0)
299  ||
300  (dynamic_cast<const FE_DGPMonomial<dim>*>(&x_source_fe) != nullptr),
301  typename FiniteElement<dim>::
302  ExcInterpolationNotImplemented());
303 
304  Assert (interpolation_matrix.m() == 0,
305  ExcDimensionMismatch (interpolation_matrix.m(),
306  0));
307  Assert (interpolation_matrix.n() == 0,
308  ExcDimensionMismatch (interpolation_matrix.n(),
309  0));
310 }
311 
312 
313 
314 template <int dim>
315 bool
317 {
318  return true;
319 }
320 
321 
322 
323 template <int dim>
324 std::vector<std::pair<unsigned int, unsigned int> >
327 {
328  // there are no such constraints for DGPMonomial
329  // elements at all
330  if (dynamic_cast<const FE_DGPMonomial<dim>*>(&fe_other) != nullptr)
331  return
332  std::vector<std::pair<unsigned int, unsigned int> > ();
333  else
334  {
335  Assert (false, ExcNotImplemented());
336  return std::vector<std::pair<unsigned int, unsigned int> > ();
337  }
338 }
339 
340 
341 
342 template <int dim>
343 std::vector<std::pair<unsigned int, unsigned int> >
346 {
347  // there are no such constraints for DGPMonomial
348  // elements at all
349  if (dynamic_cast<const FE_DGPMonomial<dim>*>(&fe_other) != nullptr)
350  return
351  std::vector<std::pair<unsigned int, unsigned int> > ();
352  else
353  {
354  Assert (false, ExcNotImplemented());
355  return std::vector<std::pair<unsigned int, unsigned int> > ();
356  }
357 }
358 
359 
360 
361 template <int dim>
362 std::vector<std::pair<unsigned int, unsigned int> >
365 {
366  // there are no such constraints for DGPMonomial
367  // elements at all
368  if (dynamic_cast<const FE_DGPMonomial<dim>*>(&fe_other) != nullptr)
369  return
370  std::vector<std::pair<unsigned int, unsigned int> > ();
371  else
372  {
373  Assert (false, ExcNotImplemented());
374  return std::vector<std::pair<unsigned int, unsigned int> > ();
375  }
376 }
377 
378 
379 
380 template <int dim>
384 {
385  // check whether both are discontinuous
386  // elements, see
387  // the description of
388  // FiniteElementDomination::Domination
389  if (dynamic_cast<const FE_DGPMonomial<dim>*>(&fe_other) != nullptr)
391 
392  Assert (false, ExcNotImplemented());
394 }
395 
396 
397 
398 template <>
399 bool
401  const unsigned int face_index) const
402 {
403  return face_index==1 || (face_index==0 && this->degree==0);
404 }
405 
406 
407 
408 template <>
409 bool
410 FE_DGPMonomial<2>::has_support_on_face (const unsigned int shape_index,
411  const unsigned int face_index) const
412 {
413  bool support_on_face=false;
414  if (face_index==1 || face_index==2)
415  support_on_face=true;
416  else
417  {
418  unsigned int degrees[2];
419  this->poly_space.directional_degrees(shape_index, degrees);
420  if ((face_index==0 && degrees[1]==0) ||
421  (face_index==3 && degrees[0]==0))
422  support_on_face=true;
423  }
424  return support_on_face;
425 }
426 
427 
428 
429 template <>
430 bool
431 FE_DGPMonomial<3>::has_support_on_face (const unsigned int shape_index,
432  const unsigned int face_index) const
433 {
434  bool support_on_face=false;
435  if (face_index==1 || face_index==3 || face_index==4)
436  support_on_face=true;
437  else
438  {
439  unsigned int degrees[3];
440  this->poly_space.directional_degrees(shape_index, degrees);
441  if ((face_index==0 && degrees[1]==0) ||
442  (face_index==2 && degrees[2]==0) ||
443  (face_index==5 && degrees[0]==0))
444  support_on_face=true;
445  }
446  return support_on_face;
447 }
448 
449 
450 
451 template <int dim>
452 std::size_t
454 {
455  Assert (false, ExcNotImplemented ());
456  return 0;
457 }
458 
459 
460 
461 // explicit instantiations
462 #include "fe_dgp_monomial.inst"
463 
464 
465 DEAL_II_NAMESPACE_CLOSE
size_type m() const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
virtual void get_subface_interpolation_matrix(const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
FE_DGPMonomial(const unsigned int p)
virtual void get_face_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix) const
void gauss_jordan()
const unsigned int degree
Definition: fe_base.h:311
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)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim > &fe_other) const
virtual bool hp_constraints_are_implemented() const
void initialize_restriction()
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim > &fe_other) const
size_type n() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim > &fe_other) const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
virtual void get_interpolation_matrix(const FiniteElement< dim > &source, FullMatrix< double > &matrix) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::unique_ptr< FiniteElement< dim, dim > > clone() const
virtual std::string get_name() const =0
const unsigned int dofs_per_cell
Definition: fe_base.h:295
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:160
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim > &fe_other) const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::size_t memory_consumption() const
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
static ::ExceptionBase & ExcNotImplemented()
virtual std::string get_name() const
PolynomialsP< dim > poly_space
Definition: fe_poly.h:445
static ::ExceptionBase & ExcInternalError()