Reference documentation for deal.II version 8.5.1
fe_dgp_nonparametric.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 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/grid/tria.h>
19 #include <deal.II/grid/tria_iterator.h>
20 #include <deal.II/dofs/dof_accessor.h>
21 #include <deal.II/fe/fe.h>
22 #include <deal.II/fe/mapping.h>
23 #include <deal.II/fe/fe_dgp_nonparametric.h>
24 #include <deal.II/fe/fe_values.h>
25 
26 #include <sstream>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 template <int dim, int spacedim>
32  :
33  FiniteElement<dim,spacedim> (
34  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree,
35  FiniteElementData<dim>::L2),
36  std::vector<bool>(
37  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree).dofs_per_cell,true),
38  std::vector<ComponentMask>(
39  FiniteElementData<dim>(get_dpo_vector(degree),1, degree).dofs_per_cell,
40  std::vector<bool>(1,true))),
41  polynomial_space (Polynomials::Legendre::generate_complete_basis(degree))
42 {
43  const unsigned int n_dofs = this->dofs_per_cell;
44  for (unsigned int ref_case = RefinementCase<dim>::cut_x;
45  ref_case<RefinementCase<dim>::isotropic_refinement+1; ++ref_case)
46  {
47  if (dim!=2 && ref_case!=RefinementCase<dim>::isotropic_refinement)
48  // do nothing, as anisotropic
49  // refinement is not
50  // implemented so far
51  continue;
52 
53  const unsigned int nc = GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
54  for (unsigned int i=0; i<nc; ++i)
55  {
56  this->prolongation[ref_case-1][i].reinit (n_dofs, n_dofs);
57  // Fill prolongation matrices with
58  // embedding operators
59  for (unsigned int j=0; j<n_dofs; ++j)
60  this->prolongation[ref_case-1][i](j,j) = 1.;
61  }
62  }
63 
64  // restriction can be defined
65  // through projection for
66  // discontinuous elements, but is
67  // presently not implemented for DGPNonparametric
68  // elements.
69  //
70  // if it were, then the following
71  // snippet would be the right code
72 // if ((degree < Matrices::n_projection_matrices) &&
73 // (Matrices::projection_matrices[degree] != 0))
74 // {
75 // restriction[0].fill (Matrices::projection_matrices[degree]);
76 // }
77 // else
78 // // matrix undefined, set size to zero
79 // for (unsigned int i=0;i<GeometryInfo<dim>::max_children_per_cell;++i)
80 // restriction[i].reinit(0, 0);
81  // since not implemented, set to
82  // "empty". however, that is done in the
83  // default constructor already, so do nothing
84 // for (unsigned int i=0;i<GeometryInfo<dim>::max_children_per_cell;++i)
85 // this->restriction[i].reinit(0, 0);
86 
87  // note further, that these
88  // elements have neither support
89  // nor face-support points, so
90  // leave these fields empty
91 }
92 
93 
94 
95 template <int dim, int spacedim>
96 std::string
98 {
99  // note that the
100  // FETools::get_fe_by_name
101  // function depends on the
102  // particular format of the string
103  // this function returns, so they
104  // have to be kept in synch
105 
106  std::ostringstream namebuf;
107  namebuf << "FE_DGPNonparametric<"
108  << Utilities::dim_string(dim,spacedim)
109  << ">(" << this->degree << ")";
110 
111  return namebuf.str();
112 }
113 
114 
115 
116 template <int dim, int spacedim>
119 {
120  return new FE_DGPNonparametric<dim,spacedim>(*this);
121 }
122 
123 
124 
125 template <int dim, int spacedim>
126 double
128  const Point<dim> &p) const
129 {
130  (void)i;
131  (void)p;
132  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
134  return 0;
135 }
136 
137 
138 
139 template <int dim, int spacedim>
140 double
142  const Point<dim> &p,
143  const unsigned int component) const
144 {
145  (void)i;
146  (void)p;
147  (void)component;
148  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
149  Assert (component == 0, ExcIndexRange (component, 0, 1));
151  return 0;
152 }
153 
154 
155 
156 template <int dim, int spacedim>
159  const Point<dim> &p) const
160 {
161  (void)i;
162  (void)p;
163  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
165  return Tensor<1,dim>();
166 }
167 
168 
169 template <int dim, int spacedim>
172  const Point<dim> &p,
173  const unsigned int component) const
174 {
175  (void)i;
176  (void)p;
177  (void)component;
178  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
179  Assert (component == 0, ExcIndexRange (component, 0, 1));
181  return Tensor<1,dim>();
182 }
183 
184 
185 
186 template <int dim, int spacedim>
189  const Point<dim> &p) const
190 {
191  (void)i;
192  (void)p;
193  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
195  return Tensor<2,dim>();
196 }
197 
198 
199 
200 template <int dim, int spacedim>
203  const Point<dim> &p,
204  const unsigned int component) const
205 {
206  (void)i;
207  (void)p;
208  (void)component;
209  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
210  Assert (component == 0, ExcIndexRange (component, 0, 1));
212  return Tensor<2,dim>();
213 }
214 
215 
216 //---------------------------------------------------------------------------
217 // Auxiliary functions
218 //---------------------------------------------------------------------------
219 
220 
221 template <int dim, int spacedim>
222 std::vector<unsigned int>
224 {
225  std::vector<unsigned int> dpo(dim+1, static_cast<unsigned int>(0));
226  dpo[dim] = deg+1;
227  for (unsigned int i=1; i<dim; ++i)
228  {
229  dpo[dim] *= deg+1+i;
230  dpo[dim] /= i+1;
231  }
232  return dpo;
233 }
234 
235 
236 
237 template <int dim, int spacedim>
240 {
241  UpdateFlags out = flags;
242 
244  out |= update_quadrature_points ;
245 
246  return out;
247 }
248 
249 
250 
251 //---------------------------------------------------------------------------
252 // Data field initialization
253 //---------------------------------------------------------------------------
254 
255 template <int dim, int spacedim>
258 get_data (const UpdateFlags update_flags,
259  const Mapping<dim,spacedim> &,
260  const Quadrature<dim> &,
262 {
263  // generate a new data object
266  data->update_each = requires_update_flags(update_flags);
267 
268  // other than that, there is nothing we can add here as discussed
269  // in the general documentation of this class
270 
271  return data;
272 }
273 
274 
275 
276 //---------------------------------------------------------------------------
277 // Fill data of FEValues
278 //---------------------------------------------------------------------------
279 
280 template <int dim, int spacedim>
281 void
285  const Quadrature<dim> &,
286  const Mapping<dim,spacedim> &,
288  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
289  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
291 {
293 
294  const unsigned int n_q_points = mapping_data.quadrature_points.size();
295 
296  std::vector<double> values((fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
297  std::vector<Tensor<1,dim> > grads((fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
298  std::vector<Tensor<2,dim> > grad_grads((fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
299  std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
300  std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
301 
302  if (fe_internal.update_each & (update_values | update_gradients))
303  for (unsigned int i=0; i<n_q_points; ++i)
304  {
305  polynomial_space.compute(mapping_data.quadrature_points[i],
306  values, grads, grad_grads,
307  empty_vector_of_3rd_order_tensors,
308  empty_vector_of_4th_order_tensors);
309 
310  if (fe_internal.update_each & update_values)
311  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
312  output_data.shape_values[k][i] = values[k];
313 
314  if (fe_internal.update_each & update_gradients)
315  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
316  output_data.shape_gradients[k][i] = grads[k];
317 
318  if (fe_internal.update_each & update_hessians)
319  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
320  output_data.shape_hessians[k][i] = grad_grads[k];
321  }
322 }
323 
324 
325 
326 template <int dim, int spacedim>
327 void
330  const unsigned int ,
331  const Quadrature<dim-1> &,
332  const Mapping<dim,spacedim> &,
334  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
335  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
337 {
339 
340  const unsigned int n_q_points = mapping_data.quadrature_points.size();
341 
342  std::vector<double> values((fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
343  std::vector<Tensor<1,dim> > grads((fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
344  std::vector<Tensor<2,dim> > grad_grads((fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
345  std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
346  std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
347 
348  if (fe_internal.update_each & (update_values | update_gradients))
349  for (unsigned int i=0; i<n_q_points; ++i)
350  {
351  polynomial_space.compute(mapping_data.quadrature_points[i],
352  values, grads, grad_grads,
353  empty_vector_of_3rd_order_tensors,
354  empty_vector_of_4th_order_tensors);
355 
356  if (fe_internal.update_each & update_values)
357  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
358  output_data.shape_values[k][i] = values[k];
359 
360  if (fe_internal.update_each & update_gradients)
361  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
362  output_data.shape_gradients[k][i] = grads[k];
363 
364  if (fe_internal.update_each & update_hessians)
365  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
366  output_data.shape_hessians[k][i] = grad_grads[k];
367  }
368 }
369 
370 
371 
372 template <int dim, int spacedim>
373 void
376  const unsigned int ,
377  const unsigned int ,
378  const Quadrature<dim-1> &,
379  const Mapping<dim,spacedim> &,
381  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
382  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
384 {
386 
387  const unsigned int n_q_points = mapping_data.quadrature_points.size();
388 
389  std::vector<double> values((fe_internal.update_each & update_values) ? this->dofs_per_cell : 0);
390  std::vector<Tensor<1,dim> > grads((fe_internal.update_each & update_gradients) ? this->dofs_per_cell : 0);
391  std::vector<Tensor<2,dim> > grad_grads((fe_internal.update_each & update_hessians) ? this->dofs_per_cell : 0);
392  std::vector<Tensor<3,dim> > empty_vector_of_3rd_order_tensors;
393  std::vector<Tensor<4,dim> > empty_vector_of_4th_order_tensors;
394 
395  if (fe_internal.update_each & (update_values | update_gradients))
396  for (unsigned int i=0; i<n_q_points; ++i)
397  {
398  polynomial_space.compute(mapping_data.quadrature_points[i],
399  values, grads, grad_grads,
400  empty_vector_of_3rd_order_tensors,
401  empty_vector_of_4th_order_tensors);
402 
403  if (fe_internal.update_each & update_values)
404  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
405  output_data.shape_values[k][i] = values[k];
406 
407  if (fe_internal.update_each & update_gradients)
408  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
409  output_data.shape_gradients[k][i] = grads[k];
410 
411  if (fe_internal.update_each & update_hessians)
412  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
413  output_data.shape_hessians[k][i] = grad_grads[k];
414  }
415 }
416 
417 
418 
419 template <int dim, int spacedim>
420 void
423  FullMatrix<double> &interpolation_matrix) const
424 {
425  // this is only implemented, if the source
426  // FE is also a DGPNonparametric element. in that case,
427  // both elements have no dofs on their
428  // faces and the face interpolation matrix
429  // is necessarily empty -- i.e. there isn't
430  // much we need to do here.
431  (void)interpolation_matrix;
432  typedef FiniteElement<dim,spacedim> FEE;
433  AssertThrow ((x_source_fe.get_name().find ("FE_DGPNonparametric<") == 0)
434  ||
435  (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&x_source_fe) != 0),
436  typename FEE::
437  ExcInterpolationNotImplemented());
438 
439  Assert (interpolation_matrix.m() == 0,
440  ExcDimensionMismatch (interpolation_matrix.m(),
441  0));
442  Assert (interpolation_matrix.n() == 0,
443  ExcDimensionMismatch (interpolation_matrix.n(),
444  0));
445 }
446 
447 
448 
449 template <int dim, int spacedim>
450 void
453  const unsigned int ,
454  FullMatrix<double> &interpolation_matrix) const
455 {
456  // this is only implemented, if the source
457  // FE is also a DGPNonparametric element. in that case,
458  // both elements have no dofs on their
459  // faces and the face interpolation matrix
460  // is necessarily empty -- i.e. there isn't
461  // much we need to do here.
462  (void)interpolation_matrix;
463  typedef FiniteElement<dim,spacedim> FEE;
464  AssertThrow ((x_source_fe.get_name().find ("FE_DGPNonparametric<") == 0)
465  ||
466  (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&x_source_fe) != 0),
467  typename FEE::
468  ExcInterpolationNotImplemented());
469 
470  Assert (interpolation_matrix.m() == 0,
471  ExcDimensionMismatch (interpolation_matrix.m(),
472  0));
473  Assert (interpolation_matrix.n() == 0,
474  ExcDimensionMismatch (interpolation_matrix.n(),
475  0));
476 }
477 
478 
479 
480 template <int dim, int spacedim>
481 bool
483 {
484  return true;
485 }
486 
487 
488 
489 template <int dim, int spacedim>
490 std::vector<std::pair<unsigned int, unsigned int> >
493 {
494  // there are no such constraints for DGPNonparametric
495  // elements at all
496  if (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&fe_other) != 0)
497  return
498  std::vector<std::pair<unsigned int, unsigned int> > ();
499  else
500  {
501  Assert (false, ExcNotImplemented());
502  return std::vector<std::pair<unsigned int, unsigned int> > ();
503  }
504 }
505 
506 
507 
508 template <int dim, int spacedim>
509 std::vector<std::pair<unsigned int, unsigned int> >
512 {
513  // there are no such constraints for DGPNonparametric
514  // elements at all
515  if (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&fe_other) != 0)
516  return
517  std::vector<std::pair<unsigned int, unsigned int> > ();
518  else
519  {
520  Assert (false, ExcNotImplemented());
521  return std::vector<std::pair<unsigned int, unsigned int> > ();
522  }
523 }
524 
525 
526 
527 template <int dim, int spacedim>
528 std::vector<std::pair<unsigned int, unsigned int> >
531 {
532  // there are no such constraints for DGPNonparametric
533  // elements at all
534  if (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&fe_other) != 0)
535  return
536  std::vector<std::pair<unsigned int, unsigned int> > ();
537  else
538  {
539  Assert (false, ExcNotImplemented());
540  return std::vector<std::pair<unsigned int, unsigned int> > ();
541  }
542 }
543 
544 
545 
546 template <int dim, int spacedim>
550 {
551  // check whether both are discontinuous
552  // elements, see the description of
553  // FiniteElementDomination::Domination
554  if (dynamic_cast<const FE_DGPNonparametric<dim,spacedim>*>(&fe_other) != 0)
556 
557  Assert (false, ExcNotImplemented());
559 }
560 
561 
562 
563 template <int dim, int spacedim>
564 bool
566  const unsigned int) const
567 {
568  return true;
569 }
570 
571 
572 
573 template <int dim, int spacedim>
574 std::size_t
576 {
577  Assert (false, ExcNotImplemented ());
578  return 0;
579 }
580 
581 
582 
583 template <int dim, int spacedim>
584 unsigned int
586 {
587  return this->degree;
588 }
589 
590 
591 
592 // explicit instantiations
593 #include "fe_dgp_nonparametric.inst"
594 
595 
596 DEAL_II_NAMESPACE_CLOSE
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Shape function values.
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
size_type m() const
friend class FE_DGPNonparametric
virtual std::size_t memory_consumption() const
STL namespace.
Transformed quadrature points.
virtual bool hp_constraints_are_implemented() const
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual FiniteElement< dim, spacedim > * clone() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
size_type n() const
unsigned int get_degree() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2199
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual std::string get_name() const =0
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:161
const unsigned int dofs_per_cell
Definition: fe_base.h:295
Shape function gradients.
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
static ::ExceptionBase & ExcNotImplemented()
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
UpdateFlags update_each
Definition: fe.h:638
static ::ExceptionBase & ExcInternalError()
virtual std::string get_name() const