Dune::P0Function< GV, RT, m > Class Template Reference
[Functions Hierarchy]

#include <p0function.hh>

Inheritance diagram for Dune::P0Function< GV, RT, m >:

Dune::GridFunctionGlobalEvalDefault< GV, RT, m > Dune::FunctionDefault< GV::Grid::ctype, RT, GV::Grid::dimension, m > Dune::L2Function< GV::Grid::ctype, RT, GV::Grid::dimension, m > Dune::FunctionBase< DT, RT, n, m > Dune::FunctionBase< DT, RT, n, m > List of all members.

Detailed Description

template<class GV, class RT, int m = 1>
class Dune::P0Function< GV, RT, m >

class for P0 finite element functions on a grid

This class implements the interface of a DifferentiableGridFunction with piecewise linear elements using a Lagrange basis. It is implemented using the general shape functions, thus it should work for all element types and dimensions.

In addition to the DifferentiableGridFunction interface P0 functions can be initialized from a C0GridFunction via Lagrange interpolation. Dereferencing delivers the coefficient vector.


Public Types

enum  
 export dimension of domain and range
typedef DT DomainFieldType
 export type for domain components
typedef RT RangeFieldType
 export type for range components

Public Member Functions

 P0Function (const GV &gv_)
 create function with undefined degrees of freedom
 ~P0Function ()
 deallocate the vector
virtual RT evallocal (int comp, const Entity &e, const Dune::FieldVector< DT, n > &xi) const
 evaluate single component comp in the entity e at local coordinates xi
virtual void evalalllocal (const Entity &e, const Dune::FieldVector< DT, G::dimension > &xi, Dune::FieldVector< RT, m > &y) const
 evaluate all components in the entity e at local coordinates xi
void interpolate (const C0GridFunction< G, RT, m > &u)
 interpolate nodal values from a grid function
const RepresentationType & operator * () const
 return const reference to coefficient vector
RepresentationType & operator * ()
 return reference to coefficient vector
void preAdapt ()
void postAdapt ()
 Initiate update process.
const MultipleCodimMultipleGeomTypeMapper<
G, typename GV::IndexSet,
P0Layout > & 
mapper () const
 export the mapper for external use
void vtkout (VTKWriter< GV > &vtkwriter, std::string s) const
 VTK output.
virtual RT evallocal (int comp, const Entity &e, const Dune::FieldVector< DT, n > &xi) const =0
 evaluate single component comp in the entity e at local coordinates xi
virtual void evalalllocal (const Entity &e, const Dune::FieldVector< DT, GV::Grid::dimension > &xi, Dune::FieldVector< RT, m > &y) const=0
 evaluate all components in the entity e at local coordinates xi
virtual RT eval (int comp, const Dune::FieldVector< DT, n > &x) const=0
 evaluate single component comp at global point x
virtual void evalall (const Dune::FieldVector< DT, n > &x, Dune::FieldVector< RT, m > &y) const=0
 evaluate all components at point x and store result in y
virtual RT eval (int comp, const Dune::FieldVector< DT, n > &xi) const
 implement global evaluation with local evaluation
virtual void evalall (const Dune::FieldVector< GV::Grid::ctype, n > &x, Dune::FieldVector< RT, m > &y) const
 default implemention for evaluation of all components

Member Function Documentation

template<class GV, class RT, int m = 1>
virtual RT Dune::P0Function< GV, RT, m >::evallocal ( int  comp,
const Entity &  e,
const Dune::FieldVector< DT, n > &  xi 
) const [inline, virtual]

evaluate single component comp in the entity e at local coordinates xi

Evaluate the function in an entity at local coordinates.

Parameters:
[in] comp number of component to be evaluated
[in] e reference to grid entity of codimension 0
[in] xi point in local coordinates of the reference element of e
Returns:
value of the component

template<class GV, class RT, int m = 1>
virtual void Dune::P0Function< GV, RT, m >::evalalllocal ( const Entity &  e,
const Dune::FieldVector< DT, G::dimension > &  xi,
Dune::FieldVector< RT, m > &  y 
) const [inline, virtual]

evaluate all components in the entity e at local coordinates xi

Evaluates all components of a function at once.

Parameters:
[in] e reference to grid entity of codimension 0
[in] xi point in local coordinates of the reference element of e
[out] y vector with values to be filled

template<class GV, class RT, int m = 1>
void Dune::P0Function< GV, RT, m >::interpolate ( const C0GridFunction< G, RT, m > &  u  )  [inline]

interpolate nodal values from a grid function

Lagrange interpolation of a P0 finite element function from given continuous grid function. Evaluation is done by visiting the vertices of each element and storing a bitvector of visited vertices.

Parameters:
[in] u a continuous grid function

template<class GV, class RT, int m = 1>
const RepresentationType& Dune::P0Function< GV, RT, m >::operator * (  )  const [inline]

return const reference to coefficient vector

Dereferencing a finite element function returns the coefficient representation of the finite element function. This is the const version.

template<class GV, class RT, int m = 1>
RepresentationType& Dune::P0Function< GV, RT, m >::operator * (  )  [inline]

return reference to coefficient vector

Dereferencing a finite element function returns the coefficient representation of the finite element function. This is the non-const version.

template<class GV, class RT, int m = 1>
void Dune::P0Function< GV, RT, m >::preAdapt (  )  [inline]

empty method to maintain symmetry For vertex data nothing is required in preAdapt but for other finite element functions this method is necessary.

template<class GV, class RT, int m = 1>
void Dune::P0Function< GV, RT, m >::postAdapt (  )  [inline]

Initiate update process.

Call this method after the grid has been adapted. The representation is now updated to the new grid and the finite element function can be used on the new grid. However the data is not initialized. The old representation (with respect to the old grid) can still be accessed if it has been saved. It is deleted in endUpdate().

virtual RT Dune::GridFunction< GV::Grid , RT , m >::evallocal ( int  comp,
const Entity &  e,
const Dune::FieldVector< DT, n > &  xi 
) const [pure virtual, inherited]

evaluate single component comp in the entity e at local coordinates xi

Evaluate the function in an entity at local coordinates.

Parameters:
[in] comp number of component to be evaluated
[in] e reference to grid entity of codimension 0
[in] xi point in local coordinates of the reference element of e
Returns:
value of the component

virtual void Dune::GridFunction< GV::Grid , RT , m >::evalalllocal ( const Entity &  e,
const Dune::FieldVector< DT, GV::Grid ::dimension > &  xi,
Dune::FieldVector< RT , m > &  y 
) const [pure virtual, inherited]

evaluate all components in the entity e at local coordinates xi

Evaluates all components of a function at once.

Parameters:
[in] e reference to grid entity of codimension 0
[in] xi point in local coordinates of the reference element of e
[out] y vector with values to be filled

template<class DT, class RT, int n, int m>
virtual RT Dune::FunctionBase< DT, RT, n, m >::eval ( int  comp,
const Dune::FieldVector< DT, n > &  x 
) const [pure virtual, inherited]

template<class DT, class RT, int n, int m>
virtual void Dune::FunctionBase< DT, RT, n, m >::evalall ( const Dune::FieldVector< DT, n > &  x,
Dune::FieldVector< RT, m > &  y 
) const [pure virtual, inherited]

evaluate all components at point x and store result in y

Evaluation function for all components at once.

Parameters:
[in] x position to be evaluated
[out] y result vector to be filled

Implemented in Dune::FunctionDefault< DT, RT, n, m >.

template<class GV, class RT, int m>
virtual RT Dune::GridFunctionGlobalEvalDefault< GV, RT, m >::eval ( int  comp,
const Dune::FieldVector< DT, n > &  xi 
) const [inline, virtual, inherited]

implement global evaluation with local evaluation

The local entity is searched via a hierarchic search

Parameters:
[in] comp number of component to be evaluated
[in] xi point in local coordinates of the reference element of e
Returns:
value of the component

virtual void Dune::FunctionDefault< GV::Grid::ctype , RT , n, m >::evalall ( const Dune::FieldVector< GV::Grid::ctype , n > &  x,
Dune::FieldVector< RT , m > &  y 
) const [inline, virtual, inherited]

default implemention for evaluation of all components

Evaluate all components at once using componentwise evaluation function

Parameters:
[in] x position to be evaluated
[out] y result vector to be filled


The documentation for this class was generated from the following file:

Generated on 6 Jan 2009 with Doxygen (ver 1.5.1) [logfile].