DUNE PDELab (git)

Dune::Functions::DefaultLocalView< GB > Class Template Reference

The restriction of a finite element basis to a single element. More...

#include <dune/functions/functionspacebases/defaultlocalview.hh>

Public Types

using GlobalBasis = GB
 The global FE basis that this is a view on.
 
using GridView = typename GlobalBasis::GridView
 The grid view the global FE basis lives on.
 
using Element = typename GridView::template Codim< 0 >::Entity
 Type of the grid element we are bound to.
 
using size_type = std::size_t
 The type used for sizes.
 
using Tree = typename GlobalBasis::PreBasis::Node
 Tree of local finite elements / local shape function sets.
 
using MultiIndex = std::conditional_t<(PreBasis::minMultiIndexSize==PreBasis::maxMultiIndexSize), StaticMultiIndex< size_type, PreBasis::maxMultiIndexSize >, Dune::ReservedVector< size_type, PreBasis::multiIndexBufferSize > >
 Type used for global numbering of the basis vectors.
 

Public Member Functions

 DefaultLocalView (const GlobalBasis &globalBasis)
 Construct local view for a given global finite element basis.
 
void bind (const Element &e)
 Bind the view to a grid element. More...
 
bool bound () const
 Return if the view is bound to a grid element.
 
const Elementelement () const
 Return the grid element that the view is bound to. More...
 
void unbind ()
 Unbind from the current element. More...
 
const Treetree () const
 Return the local ansatz tree associated to the bound entity. More...
 
size_type size () const
 Total number of degrees of freedom on this element.
 
size_type maxSize () const
 Maximum local size for any element on the GridView. More...
 
const MultiIndexindex (size_type i) const
 Maps from subtree index set [0..size-1] to a globally unique multi index in global basis.
 
const GlobalBasisglobalBasis () const
 Return the global basis that we are a view on.
 

Detailed Description

template<class GB>
class Dune::Functions::DefaultLocalView< GB >

The restriction of a finite element basis to a single element.

Member Function Documentation

◆ bind()

template<class GB >
void Dune::Functions::DefaultLocalView< GB >::bind ( const Element e)
inline

Bind the view to a grid element.

Having to bind the view to an element before being able to actually access any of its data members offers to centralize some expensive setup code in the 'bind' method, which can save a lot of run-time.

References Dune::Functions::DefaultLocalView< GB >::size().

◆ element()

template<class GB >
const Element& Dune::Functions::DefaultLocalView< GB >::element ( ) const
inline

Return the grid element that the view is bound to.

Exceptions
Dune::Exceptionif the view is not bound to anything

◆ maxSize()

template<class GB >
size_type Dune::Functions::DefaultLocalView< GB >::maxSize ( ) const
inline

Maximum local size for any element on the GridView.

This is the maximal size needed for local matrices and local vectors, i.e., the result is

◆ tree()

template<class GB >
const Tree& Dune::Functions::DefaultLocalView< GB >::tree ( ) const
inline

Return the local ansatz tree associated to the bound entity.

Returns
Tree // This is tree

◆ unbind()

template<class GB >
void Dune::Functions::DefaultLocalView< GB >::unbind ( )
inline

Unbind from the current element.

Calling this method should only be a hint that the view can be unbound.


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.80.0 (May 16, 22:29, 2024)