Reference documentation for deal.II version 8.5.1
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 #include <deal.II/fe/mapping.h>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
29 
268 template <int dim, int spacedim=dim>
269 class FE_DGPNonparametric : public FiniteElement<dim,spacedim>
270 {
271 public:
275  FE_DGPNonparametric (const unsigned int k);
276 
282  virtual std::string get_name () const;
283 
284  // for documentation, see the FiniteElement base class
285  virtual
287  requires_update_flags (const UpdateFlags update_flags) const;
288 
299  virtual double shape_value (const unsigned int i,
300  const Point<dim> &p) const;
301 
312  virtual double shape_value_component (const unsigned int i,
313  const Point<dim> &p,
314  const unsigned int component) const;
315 
326  virtual Tensor<1,dim> shape_grad (const unsigned int i,
327  const Point<dim> &p) const;
328 
339  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
340  const Point<dim> &p,
341  const unsigned int component) const;
342 
353  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
354  const Point<dim> &p) const;
355 
366  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
367  const Point<dim> &p,
368  const unsigned int component) const;
369 
374  unsigned int get_degree () const;
375 
387  virtual void
389  FullMatrix<double> &matrix) const;
390 
402  virtual void
404  const unsigned int subface,
405  FullMatrix<double> &matrix) const;
406 
430  virtual
431  std::vector<std::pair<unsigned int, unsigned int> >
433 
441  virtual
442  std::vector<std::pair<unsigned int, unsigned int> >
443  hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
444 
452  virtual
453  std::vector<std::pair<unsigned int, unsigned int> >
454  hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
455 
464  virtual bool hp_constraints_are_implemented () const;
465 
475  virtual
478 
487  virtual bool has_support_on_face (const unsigned int shape_index,
488  const unsigned int face_index) const;
489 
498  virtual std::size_t memory_consumption () const;
499 
500 
501 private:
508  struct Matrices
509  {
515 
521  static const unsigned int n_embedding_matrices;
522 
527 
531  static const unsigned int n_projection_matrices;
532  };
533 
534 protected:
535 
541  virtual FiniteElement<dim,spacedim> *clone() const;
542 
547  virtual
549  get_data (const UpdateFlags update_flags,
550  const Mapping<dim,spacedim> &mapping,
551  const Quadrature<dim> &quadrature,
553 
554  virtual
555  void
556  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
557  const CellSimilarity::Similarity cell_similarity,
558  const Quadrature<dim> &quadrature,
559  const Mapping<dim,spacedim> &mapping,
560  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
561  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
562  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
564 
565  virtual
566  void
567  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
568  const unsigned int face_no,
569  const Quadrature<dim-1> &quadrature,
570  const Mapping<dim,spacedim> &mapping,
571  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
572  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
573  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
575 
576  virtual
577  void
578  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
579  const unsigned int face_no,
580  const unsigned int sub_no,
581  const Quadrature<dim-1> &quadrature,
582  const Mapping<dim,spacedim> &mapping,
583  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
584  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
585  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
587 
588 private:
589 
596  static
597  std::vector<unsigned int>
598  get_dpo_vector (const unsigned int degree);
599 
604 
608  template <int, int> friend class FE_DGPNonparametric;
609 };
610 
613 #ifndef DOXYGEN
614 
615 // declaration of explicit specializations of member variables, if the
616 // compiler allows us to do that (the standard says we must)
617 #ifndef DEAL_II_MEMBER_VAR_SPECIALIZATION_BUG
618 template <>
620 
621 template <>
623 
624 template <>
626 
627 template <>
629 
630 template <>
632 
633 template <>
635 
636 template <>
638 
639 template <>
641 
642 template <>
644 
645 template <>
647 
648 template <>
650 
651 template <>
653 #endif
654 
655 #endif // DOXYGEN
656 
657 DEAL_II_NAMESPACE_CLOSE
658 
659 #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 unsigned int degree
Definition: fe_base.h:311
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:46
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)
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