Reference documentation for deal.II version 8.4.2
fe_dgp_nonparametric.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2016 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 #ifndef dealii__fe_dgp_nonparametric_h
17 #define dealii__fe_dgp_nonparametric_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/polynomial.h>
21 #include <deal.II/base/polynomial_space.h>
22 #include <deal.II/fe/fe.h>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 template <int dim> class PolynomialSpace;
27 template <int dim, int spacedim> class MappingQ;
28 
29 
32 
271 template <int dim, int spacedim=dim>
272 class FE_DGPNonparametric : public FiniteElement<dim,spacedim>
273 {
274 public:
278  FE_DGPNonparametric (const unsigned int k);
279 
285  virtual std::string get_name () const;
286 
287  // for documentation, see the FiniteElement base class
288  virtual
290  requires_update_flags (const UpdateFlags update_flags) const;
291 
302  virtual double shape_value (const unsigned int i,
303  const Point<dim> &p) const;
304 
315  virtual double shape_value_component (const unsigned int i,
316  const Point<dim> &p,
317  const unsigned int component) const;
318 
329  virtual Tensor<1,dim> shape_grad (const unsigned int i,
330  const Point<dim> &p) const;
331 
342  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
343  const Point<dim> &p,
344  const unsigned int component) const;
345 
356  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
357  const Point<dim> &p) const;
358 
369  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
370  const Point<dim> &p,
371  const unsigned int component) const;
372 
377  unsigned int get_degree () const;
378 
390  virtual void
392  FullMatrix<double> &matrix) const;
393 
405  virtual void
407  const unsigned int subface,
408  FullMatrix<double> &matrix) const;
409 
433  virtual
434  std::vector<std::pair<unsigned int, unsigned int> >
436 
444  virtual
445  std::vector<std::pair<unsigned int, unsigned int> >
446  hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
447 
455  virtual
456  std::vector<std::pair<unsigned int, unsigned int> >
457  hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
458 
467  virtual bool hp_constraints_are_implemented () const;
468 
478  virtual
481 
490  virtual bool has_support_on_face (const unsigned int shape_index,
491  const unsigned int face_index) const;
492 
501  virtual std::size_t memory_consumption () const;
502 
503 
504 private:
511  struct Matrices
512  {
518 
524  static const unsigned int n_embedding_matrices;
525 
530 
534  static const unsigned int n_projection_matrices;
535  };
536 
537 protected:
538 
544  virtual FiniteElement<dim,spacedim> *clone() const;
545 
550  virtual
552  get_data (const UpdateFlags update_flags,
553  const Mapping<dim,spacedim> &mapping,
554  const Quadrature<dim> &quadrature,
556 
557  virtual
558  void
559  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
560  const CellSimilarity::Similarity cell_similarity,
561  const Quadrature<dim> &quadrature,
562  const Mapping<dim,spacedim> &mapping,
563  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
564  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
565  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
567 
568  virtual
569  void
570  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
571  const unsigned int face_no,
572  const Quadrature<dim-1> &quadrature,
573  const Mapping<dim,spacedim> &mapping,
574  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
575  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
576  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
578 
579  virtual
580  void
581  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
582  const unsigned int face_no,
583  const unsigned int sub_no,
584  const Quadrature<dim-1> &quadrature,
585  const Mapping<dim,spacedim> &mapping,
586  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
587  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
588  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
590 
591 private:
592 
599  static
600  std::vector<unsigned int>
601  get_dpo_vector (const unsigned int degree);
602 
606  const unsigned int degree;
607 
612 
616  template <int, int> friend class FE_DGPNonparametric;
617 };
618 
621 #ifndef DOXYGEN
622 
623 // declaration of explicit specializations of member variables, if the
624 // compiler allows us to do that (the standard says we must)
625 #ifndef DEAL_II_MEMBER_VAR_SPECIALIZATION_BUG
626 template <>
628 
629 template <>
631 
632 template <>
634 
635 template <>
637 
638 template <>
640 
641 template <>
643 
644 template <>
646 
647 template <>
649 
650 template <>
652 
653 template <>
655 
656 template <>
658 
659 template <>
661 #endif
662 
663 #endif // DOXYGEN
664 
665 DEAL_II_NAMESPACE_CLOSE
666 
667 #endif
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
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
friend class FE_DGPNonparametric
virtual std::size_t memory_consumption() const
const PolynomialSpace< dim > polynomial_space
virtual bool hp_constraints_are_implemented() const
virtual FiniteElement< dim, spacedim > * clone() const
static const unsigned int n_embedding_matrices
static const double *const embedding[][GeometryInfo< dim >::max_children_per_cell]
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 const double *const projection_matrices[][GeometryInfo< dim >::max_children_per_cell]
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
unsigned int get_degree() const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:52
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
static const unsigned int n_projection_matrices
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
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)
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
virtual std::string get_name() const