Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
fe_dgp_nonparametric.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2018 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 
21 #include <deal.II/base/polynomial.h>
22 #include <deal.II/base/polynomial_space.h>
23 
24 #include <deal.II/fe/fe.h>
25 #include <deal.II/fe/mapping.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
31 
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
398  FullMatrix<double> &matrix) const override;
399 
411  virtual void
413  const unsigned int subface,
414  FullMatrix<double> &matrix) const override;
415 
439  virtual std::vector<std::pair<unsigned int, unsigned int>>
441  const FiniteElement<dim, spacedim> &fe_other) const override;
442 
450  virtual std::vector<std::pair<unsigned int, unsigned int>>
452  const FiniteElement<dim, spacedim> &fe_other) const override;
453 
461  virtual std::vector<std::pair<unsigned int, unsigned int>>
463  const FiniteElement<dim, spacedim> &fe_other) const override;
464 
473  virtual bool
474  hp_constraints_are_implemented() const override;
475 
481  const unsigned int codim = 0) const override final;
482 
491  virtual bool
492  has_support_on_face(const unsigned int shape_index,
493  const unsigned int face_index) const override;
494 
503  virtual std::size_t
504  memory_consumption() const override;
505 
506 protected:
511  virtual std::unique_ptr<
513  get_data(
514  const UpdateFlags update_flags,
515  const Mapping<dim, spacedim> &mapping,
516  const Quadrature<dim> & quadrature,
518  spacedim>
519  &output_data) const override;
520 
521  virtual void
522  fill_fe_values(
523  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
524  const CellSimilarity::Similarity cell_similarity,
525  const Quadrature<dim> & quadrature,
526  const Mapping<dim, spacedim> & mapping,
527  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
528  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
529  spacedim>
530  & mapping_data,
531  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
533  spacedim>
534  &output_data) const override;
535 
536  virtual void
537  fill_fe_face_values(
538  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
539  const unsigned int face_no,
540  const Quadrature<dim - 1> & quadrature,
541  const Mapping<dim, spacedim> & mapping,
542  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
543  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
544  spacedim>
545  & mapping_data,
546  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
548  spacedim>
549  &output_data) const override;
550 
551  virtual void
552  fill_fe_subface_values(
553  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
554  const unsigned int face_no,
555  const unsigned int sub_no,
556  const Quadrature<dim - 1> & quadrature,
557  const Mapping<dim, spacedim> & mapping,
558  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
559  const ::internal::FEValuesImplementation::MappingRelatedData<dim,
560  spacedim>
561  & mapping_data,
562  const typename FiniteElement<dim, spacedim>::InternalDataBase &fe_internal,
564  spacedim>
565  &output_data) const override;
566 
567 private:
574  static std::vector<unsigned int>
575  get_dpo_vector(const unsigned int degree);
576 
581 
585  template <int, int>
586  friend class FE_DGPNonparametric;
587 };
588 
591 DEAL_II_NAMESPACE_CLOSE
592 
593 #endif
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
friend class FE_DGPNonparametric
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual std::size_t memory_consumption() const override
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
const unsigned int degree
Definition: fe_base.h:298
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
const PolynomialSpace< dim > polynomial_space
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual bool hp_constraints_are_implemented() const override
unsigned int get_degree() const
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const override final
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
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
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const override
virtual std::string get_name() 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
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const override