DUNE-FEM (unstable)

Dune::Fem::hpDG::HierarchicLegendreDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage > Class Template Reference

Implementation of an \(hp\)-adaptive discrete function space using product Legendre polynomials. More...

#include <dune/fem/space/hpdg/legendre.hh>

Public Types

using KeyType = typename BasisFunctionSetsType::KeyType
 key type identifying a basis function set
 
using BasisFunctionSetType = typename BaseType::BasisFunctionSetType
 basis function set type
 
using BlockMapperType = typename BaseType::BlockMapperType
 block mapper type
 
using InterpolationImplType = DiscontinuousGalerkinLocalL2Projection< GridPartType, BasisFunctionSetType >
 local interpolation type

 
using InterpolationType = LocalInterpolationWrapper< DiscreteFunctionSpaceType >
 local interpolation type

 
typedef DofManager< GridTypeDofManagerType
 type of DoF manager
 
typedef CommunicationManager< DiscreteFunctionSpaceTypeCommunicationManagerType
 type of communication manager
 
typedef Traits::FunctionSpaceType FunctionSpaceType
 type of function space
 
typedef GridPartType::IntersectionType IntersectionType
 type of the intersections
 
typedef AuxiliaryDofsType SlaveDofsType
 deprecated type
 

Public Member Functions

bool continuous (const IntersectionType &intersection) const
 returns true if discrete functions over this space have zero jump over the given intersection. More...
 
void adapt (const Vector &polynomialOrders, const int polOrderShift=0) const
 default implementation of adapt does nothing, its only used in PAdaptiveLagrangeSpace
 
int sequence () const
 get index of the sequence in grid sequences More...
 
const GridTypegrid () const
 get reference to grid this discrete function space belongs to More...
 
GridTypegrid ()
 get reference to grid this discrete function space belongs to More...
 
GridPartType & gridPart () const
 
GridPartType & gridPart ()
 get a reference to the associated grid partition More...
 
const IndexSetTypeindexSet () const
 Get a reference to the associated index set. More...
 
int size () const
 get number of DoFs for this space More...
 
int primarySize () const
 get number of primary DoFs for this space More...
 
int auxiliarySize () const
 get number of auxiliary DoFs for this space More...
 
int maxNumDofs () const
 return the maximal number of dofs on entities
 
IteratorType begin () const
 get iterator pointing to the first entity of the associated grid partition More...
 
IteratorType end () const
 get iterator pointing behind the last entity of the associated grid partition More...
 
void forEach (FunctorType &f) const
 apply a functor to each entity in the associated grid partition More...
 
bool multipleGeometryTypes () const
 returns true if the grid has more than one geometry type More...
 
InterfaceType communicationInterface () const
 return the communication interface appropriate for this space More...
 
CommunicationDirection communicationDirection () const
 return the communication interface appropriate for this space More...
 
const CommunicationManagerTypecommunicator () const
 return reference to communicator (see CommunicationManager) More...
 
void communicate (DiscreteFunction &discreteFunction) const
 communicate data for given discrete function using the space's default communication operation More...
 
void communicate (DiscreteFunction &discreteFunction, const Operation &op) const
 communicate data for given discrete function More...
 
BaseType::template CommDataHandle< DiscreteFunction, Operation >::Type createDataHandle (DiscreteFunction &discreteFunction, const Operation &operation) const
 
const AuxiliaryDofsTypeauxiliaryDofs () const
 get auxiliary dofs
 
void addFunction (DiscreteFunction &df) const
 default implementation of addFunction does nothing at the moment
 
void removeFunction (DiscreteFunction &df) const
 default implementation of removeFunction does nothing at the moment
 
const std::vector< GeometryType > & geomTypes (int codim) const
 returns true if the grid has more than one geometry type More...
 
const AuxiliaryDofsTypeslaveDofs () const
 deprecated method, use auxiliaryDofs
 
Query methods
bool continuous () const
 please doc me
 
bool continuous (const typename BaseType::IntersectionType &intersection) const
 please doc me
 
bool multipleBasisFunctionSets () const
 please doc me
 
Basis function set methods
int order () const
 return polynomial order
 
int order (const EntityType &entity) const
 return polynomial order
 
BasisFunctionSetType basisFunctionSet (const EntityType &entity) const
 return basis function set
 
Interpolation
InterpolationType interpolation () const
 return interpolation More...
 
InterpolationImplType interpolation (const EntityType &entity) const
 return interpolation More...
 
InterpolationImplType localInterpolation (const EntityType &entity) const
 return interpolation More...
 
Block mapper
BlockMapperTypeblockMapper () const
 return block mapper
 
Adaptation
const KeyTypekey (const EntityType &entity) const
 get identifiying basis function set key assigned to given entity More...
 
void mark (const KeyType &key, const EntityType &entity)
 assign new key to given entity More...
 
KeyType getMark (const EntityType &entity) const
 get key to be assigned to an entity after next call to adapt() More...
 
bool adapt ()
 please doc me
 
bool adapt (DataProjection< DiscreteFunctionSpace, Implementation > &projection)
 please doc me
 
Deprecated methods
DFSpaceIdentifier type () const
 
Non-interface methods
const BasisFunctionSetsType & basisFunctionSets () const
 return basis function sets
 

Static Public Attributes

static constexpr std::size_t localBlockSize
 size of local blocks
 

Related Functions

(Note that these are not member functions.)

bool operator== (const DiscreteFunctionSpaceInterface< Traits > &X, const DiscreteFunctionSpaceInterface< Traits > &Y)
 check two spaces for equality More...
 

Detailed Description

template<class FunctionSpace, class GridPart, int order, class Storage>
class Dune::Fem::hpDG::HierarchicLegendreDiscontinuousGalerkinSpace< FunctionSpace, GridPart, order, Storage >

Implementation of an \(hp\)-adaptive discrete function space using product Legendre polynomials.

Template Parameters
FunctionSpacea Dune::Fem::FunctionSpace
GridParta Dune::Fem::GridPart
ordermaximum polynomial order per coordinate
Storagefor certain caching features

Member Function Documentation

◆ auxiliarySize()

int Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::auxiliarySize ( ) const
inlineinherited

get number of auxiliary DoFs for this space

Returns
number of auxiliary DoFs (degrees of freedom that are NOT owned by this process)

◆ begin()

IteratorType Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::begin ( ) const
inlineinherited

get iterator pointing to the first entity of the associated grid partition

Returns
iterator pointing to first entity
Note
The default implementation uses the codim 0 iterators of the associated grid partition.

◆ communicate() [1/2]

void Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::communicate ( DiscreteFunction &  discreteFunction) const
inlineinherited

communicate data for given discrete function using the space's default communication operation

Parameters
discreteFunctiondiscrete function to be communicated

◆ communicate() [2/2]

void Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::communicate ( DiscreteFunction &  discreteFunction,
const Operation &  op 
) const
inlineinherited

communicate data for given discrete function

Parameters
discreteFunctiondiscrete function to be communicated
[in]opcommunication operation to use (see DFCommunicationOperation)

◆ communicationDirection()

CommunicationDirection Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::communicationDirection ( ) const
inlineinherited

return the communication interface appropriate for this space

Returns
communication interface

◆ communicationInterface()

InterfaceType Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::communicationInterface ( ) const
inlineinherited

return the communication interface appropriate for this space

Returns
communication interface

◆ communicator()

const CommunicationManagerType & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::communicator ( ) const
inlineinherited

return reference to communicator (see CommunicationManager)

Returns
reference to communicator

◆ continuous()

bool Dune::Fem::DiscreteFunctionSpaceInterface< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::continuous ( const IntersectionType intersection) const
inlineinherited

returns true if discrete functions over this space have zero jump over the given intersection.

For example, a Lagrange space returns true iff the intersection is conforming while a discontiuous Galerkin space always returns false.

Parameters
intersectionIntersection for which we want to know the continuety
Returns
true if the space contians functions which are continuous over the intersection, false otherwise

◆ createDataHandle()

BaseType::template CommDataHandle< DiscreteFunction, Operation >::Type Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::createDataHandle ( DiscreteFunction &  discreteFunction,
const Operation &  operation 
) const
inlineinherited

Note
The default implementation is
return CommDataHandle< DiscreteFunction, Operation > :: Type( discreteFunction );

◆ end()

IteratorType Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::end ( ) const
inlineinherited

get iterator pointing behind the last entity of the associated grid partition

Returns
iterator pointing behind last entity
Note
The default implementation uses the codim 0 iterators of the associated grid partition.

◆ forEach()

void Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::forEach ( FunctorType &  f) const
inlineinherited

apply a functor to each entity in the associated grid partition

The functor must provide an the following operator

template< class EntityType >
void operator() ( const EntityType & );
Parameters
[in]ffunctor to apply
Note
The default implementation simply does the following:
const IteratorType end = asImp().end();
for( IteratorType it = asImp().begin(); it != end; ++it )
f( *it );
IteratorType end() const
get iterator pointing behind the last entity of the associated grid partition
Definition: discretefunctionspace.hh:817
IteratorType begin() const
get iterator pointing to the first entity of the associated grid partition
Definition: discretefunctionspace.hh:807
GridPartType::template Codim< Traits::codimension >::IteratorType IteratorType
type of iterator for grid traversal
Definition: discretefunctionspace.hh:222

◆ geomTypes()

const std::vector< GeometryType > & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::geomTypes ( int  codim) const
inlineinherited

returns true if the grid has more than one geometry type

Returns
true if the underlying grid has more than one geometry type (hybrid grid), false otherwise

◆ getMark()

KeyType Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::getMark ( const EntityType entity) const
inlineinherited

get key to be assigned to an entity after next call to adapt()

Parameters
[in]entitygrid part entity
Returns
key

◆ grid() [1/2]

GridType & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::grid ( )
inlineinherited

get reference to grid this discrete function space belongs to

Returns
reference to grid

◆ grid() [2/2]

const GridType & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::grid ( ) const
inlineinherited

get reference to grid this discrete function space belongs to

Returns
constant reference to grid

◆ gridPart() [1/2]

GridPartType & Dune::Fem::DiscreteFunctionSpaceInterface< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::gridPart ( )
inlineinherited

get a reference to the associated grid partition

Returns
reference to the grid partition

◆ gridPart() [2/2]

GridPartType & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::gridPart ( ) const
inlineinherited

◆ indexSet()

const IndexSetType & Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::indexSet ( ) const
inlineinherited

Get a reference to the associated index set.

Returns
const reference to index set

◆ interpolation() [1/2]

InterpolationType Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::interpolation ( ) const
inlineinherited

return interpolation

Parameters
[in]entitya grid part entity
Returns
local interpolation

◆ interpolation() [2/2]

InterpolationImplType Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::interpolation ( const EntityType entity) const
inlineinherited

return interpolation

Parameters
[in]entitya grid part entity
Returns
local interpolation

◆ key()

const KeyType & Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::key ( const EntityType entity) const
inlineinherited

get identifiying basis function set key assigned to given entity

Parameters
[in]entitygrid part entity
Returns
key

◆ localInterpolation()

InterpolationImplType Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::localInterpolation ( const EntityType entity) const
inlineinherited

return interpolation

Parameters
[in]entitya grid part entity
Returns
local interpolation

◆ mark()

void Dune::Fem::hpDG::DiscontinuousGalerkinSpace< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::mark ( const KeyType key,
const EntityType entity 
)
inlineinherited

assign new key to given entity

Parameters
[in]keykey identifying basis function set
[in]entitygrid part entity

◆ multipleGeometryTypes()

bool Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::multipleGeometryTypes ( ) const
inlineinherited

returns true if the grid has more than one geometry type

Returns
true if the underlying grid has more than one geometry type (hybrid grid), false otherwise

◆ primarySize()

int Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::primarySize ( ) const
inlineinherited

get number of primary DoFs for this space

Returns
number of primary DoFs (degrees of freedom that are owned by this process )

◆ sequence()

int Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::sequence ( ) const
inlineinherited

get index of the sequence in grid sequences

Returns
number of current sequence

◆ size()

int Dune::Fem::DiscreteFunctionSpaceDefault< LegendreDiscontinuousGalerkinSpaceTraits< FunctionSpace, GridPart, order, true, Storage > >::size ( ) const
inlineinherited

get number of DoFs for this space

Returns
number of DoFs (degrees of freedom)

Friends And Related Function Documentation

◆ operator==()

bool operator== ( const DiscreteFunctionSpaceInterface< Traits > &  X,
const DiscreteFunctionSpaceInterface< Traits > &  Y 
)
related

check two spaces for equality

This is a default implemented equality operator for discrete function spaces. It assumes the mapper to be a singleton and then compares the addresses of the two mappers.

Note that this method can be specialized by implementing another version that uses the exact traits of the discrete function space.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 25, 23:30, 2024)