Reference documentation for deal.II version GIT 7b2de2f2f9 2023-09-24 11:00:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_dgp_nonparametric.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2022 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.md at
12 // the top level directory of deal.II.
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 
23 
24 #include <deal.II/fe/fe.h>
25 #include <deal.II/fe/mapping.h>
26 
28 
270 template <int dim, int spacedim = dim>
271 class FE_DGPNonparametric : public FiniteElement<dim, spacedim>
272 {
273 public:
277  FE_DGPNonparametric(const unsigned int k);
278 
284  virtual std::string
285  get_name() const override;
286 
287  virtual std::unique_ptr<FiniteElement<dim, spacedim>>
288  clone() const override;
289 
290  // for documentation, see the FiniteElement base class
291  virtual UpdateFlags
292  requires_update_flags(const UpdateFlags update_flags) const override;
293 
304  virtual double
305  shape_value(const unsigned int i, const Point<dim> &p) const override;
306 
317  virtual double
318  shape_value_component(const unsigned int i,
319  const Point<dim> &p,
320  const unsigned int component) const override;
321 
332  virtual Tensor<1, dim>
333  shape_grad(const unsigned int i, const Point<dim> &p) const override;
334 
345  virtual Tensor<1, dim>
346  shape_grad_component(const unsigned int i,
347  const Point<dim> &p,
348  const unsigned int component) const override;
349 
360  virtual Tensor<2, dim>
361  shape_grad_grad(const unsigned int i, const Point<dim> &p) const override;
362 
373  virtual Tensor<2, dim>
374  shape_grad_grad_component(const unsigned int i,
375  const Point<dim> &p,
376  const unsigned int component) const override;
377 
382  unsigned int
383  get_degree() const;
384 
396  virtual void
399  const unsigned int face_no = 0) const override;
400 
412  virtual void
414  const FiniteElement<dim, spacedim> &source,
415  const unsigned int subface,
417  const unsigned int face_no = 0) const override;
418 
442  virtual std::vector<std::pair<unsigned int, unsigned int>>
444  const FiniteElement<dim, spacedim> &fe_other) const override;
445 
453  virtual std::vector<std::pair<unsigned int, unsigned int>>
455  const FiniteElement<dim, spacedim> &fe_other) const override;
456 
464  virtual std::vector<std::pair<unsigned int, unsigned int>>
466  const unsigned int face_no = 0) const override;
467 
476  virtual bool
477  hp_constraints_are_implemented() const override;
478 
484  const unsigned int codim = 0) const override final;
485 
494  virtual bool
495  has_support_on_face(const unsigned int shape_index,
496  const unsigned int face_index) const override;
497 
506  virtual std::size_t
507  memory_consumption() const override;
508 
509 protected:
514  virtual std::unique_ptr<
516  get_data(
517  const UpdateFlags update_flags,
518  const Mapping<dim, spacedim> &mapping,
519  const Quadrature<dim> &quadrature,
521  spacedim>
522  &output_data) const override;
523 
524  virtual void
526  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
527  const CellSimilarity::Similarity cell_similarity,
528  const Quadrature<dim> &quadrature,
529  const Mapping<dim, spacedim> &mapping,
530  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
532  &mapping_data,
533  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
535  spacedim>
536  &output_data) const override;
537 
539 
540  virtual void
542  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
543  const unsigned int face_no,
544  const hp::QCollection<dim - 1> &quadrature,
545  const Mapping<dim, spacedim> &mapping,
546  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
548  &mapping_data,
549  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
551  spacedim>
552  &output_data) const override;
553 
554  virtual void
556  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
557  const unsigned int face_no,
558  const unsigned int sub_no,
559  const Quadrature<dim - 1> &quadrature,
560  const Mapping<dim, spacedim> &mapping,
561  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
563  &mapping_data,
564  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
566  spacedim>
567  &output_data) const override;
568 
569 private:
576  static std::vector<unsigned int>
577  get_dpo_vector(const unsigned int degree);
578 
583 
584  // Allow access from other dimensions.
585  template <int, int>
586  friend class FE_DGPNonparametric;
587 };
588 
592 
593 #endif
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
unsigned int get_degree() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
const PolynomialSpace< dim > polynomial_space
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual bool hp_constraints_are_implemented() const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const override
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual std::string get_name() const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
friend class FE_DGPNonparametric
virtual std::size_t memory_consumption() const override
const unsigned int degree
Definition: fe_data.h:453
Definition: point.h:112
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
UpdateFlags
@ matrix
Contents is actually a matrix.