DUNE PDELab (2.7)

Dune::PDELab::VectorDiscreteGridFunctionDiv< T, X > Class Template Reference

Compute divergence of vector-valued functions. More...

#include <dune/pdelab/gridfunctionspace/gridfunctionspaceutilities.hh>

Public Types

typedef T::GridViewType GridViewType
 Type of the GridView.
 
typedef LeafNodeTag NodeTag
 The type tag that describes a LeafNode.
 

Public Member Functions

 VectorDiscreteGridFunctionDiv (const GFS &gfs, const X &x_)
 Construct a VectorDiscreteGridFunctionDiv. More...
 
 VectorDiscreteGridFunctionDiv (std::shared_ptr< const GFS > gfs, std::shared_ptr< const X > x_)
 Construct a VectorDiscreteGridFunctionDiv. More...
 
const Traits::GridViewType & getGridView () const
 get a reference to the GridView
 
Output::DataSetType dataSetType () const
 Return the data set type of this function.
 
void setDataSetType (Output::DataSetType dataSetType)
 Set the data set type of this function.
 

Static Public Attributes

static const bool isLeaf = true
 Mark this class as a leaf in a dune-typetree.
 
static const bool isPower = false
 Mark this class as a non power in the dune-typetree.
 
static const bool isComposite = false
 Mark this class as a non composite in the dune-typetree.
 
static const std::size_t CHILDREN = 0
 Leafs have no children.
 

Detailed Description

template<typename T, typename X>
class Dune::PDELab::VectorDiscreteGridFunctionDiv< T, X >

Compute divergence of vector-valued functions.

Template Parameters
TType of VectorGridFunctionSpace.
XType of coefficients vector.

The grid function space should be a vector grid function space consisting of scalar-valued function spaces. The divergence will be a single component function.

Constructor & Destructor Documentation

◆ VectorDiscreteGridFunctionDiv() [1/2]

template<typename T , typename X >
Dune::PDELab::VectorDiscreteGridFunctionDiv< T, X >::VectorDiscreteGridFunctionDiv ( const GFS &  gfs,
const X &  x_ 
)
inline

Construct a VectorDiscreteGridFunctionDiv.

Parameters
gfsGridFunctionSpace.
x_Coefficient vector.

◆ VectorDiscreteGridFunctionDiv() [2/2]

template<typename T , typename X >
Dune::PDELab::VectorDiscreteGridFunctionDiv< T, X >::VectorDiscreteGridFunctionDiv ( std::shared_ptr< const GFS >  gfs,
std::shared_ptr< const X >  x_ 
)
inline

Construct a VectorDiscreteGridFunctionDiv.

Parameters
gfsshared pointer to the GridFunctionsSpace
x_shared pointer to the coefficients vector

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)