AMDiS 2.10
The Adaptive Multi-Dimensional Simulation Toolbox
DOFVector< GB, T, Traits > Class Template Reference

The basic container that stores a base vector and a corresponding basis. More...

#include <DOFVector.hpp>

Inherits VectorFacade< double, Traits::template Vector< GB >::template Impl >, Observer< event::preAdapt >, Observer< event::adapt >, and Observer< event::postAdapt >.

Public Types

using GlobalBasis = GB
 A global basis associated to the coefficients.
 
using Coefficients = VectorFacade< T, Traits::template Vector< GB >::template Impl >
 The internal coefficient vector storage.
 
using size_type = typename GlobalBasis::size_type
 The index/size - type.
 
using value_type = T
 The type of the elements of the DOFVector.
 
- Public Types inherited from VectorFacade< double, Traits::template Vector< GB >::template Impl >
using size_type = typename Impl::size_type
 
using value_type = typename Impl::value_type
 
using HasLocalSize = decltype(std::declval< V >().localSize())
 
using HasGlobalSize = decltype(std::declval< V >().globalSize())
 

Public Member Functions

template<class DataTransferTag = tag::interpolation_datatransfer>
 DOFVector (std::shared_ptr< GB > const &basis, DataTransferTag={})
 
template<class GB_ , class DataTransferTag = tag::interpolation_datatransfer, REQUIRES(Concepts::Similar< GB_, GB >) >
 DOFVector (GB_ &&basis, DataTransferTag tag={})
 
template<class GV , class PBF , class DataTransferTag = tag::interpolation_datatransfer, REQUIRES(Concepts::PreBasisFactory< PBF, GV >) >
 DOFVector (GV const &gridView, PBF const &preBasisFactory, DataTransferTag tag={})
 
GlobalBasis const & basis () const
 Return the global basis.
 
Coefficients const & coefficients () const
 
Coefficientscoefficients ()
 
void backup (std::string const &filename)
 Write DOFVector to file.
 
void restore (std::string const &filename)
 Read backup data from file.
 
void resize ()
 
void resizeZero ()
 
auto const & dataTransfer () const
 Return the associated DataTransfer object.
 
auto & dataTransfer ()
 Return the associated DataTransfer object.
 
void setDataTransfer (DataTransfer< GB, Coefficients > dataTransfer)
 Assign the DataTransfer object.
 
template<class DataTransferTag >
void setDataTransfer (DataTransferTag tag)
 Create a new DataTransfer object based on the operation type.
 
template<class Type >
auto makeDerivative (Type type) const
 
template<class Basis >
void resize (Basis const &basis)
 Resize the vector to the size of the basis.
 
template<class Basis >
void resizeZero (Basis const &basis)
 Resize the vector to the size of the basis and set to zero.
 
- Public Member Functions inherited from VectorFacade< double, Traits::template Vector< GB >::template Impl >
 VectorFacade (GlobalBasis const &basis)
 
Impl const & impl () const
 Return the underlying linear algebra backend.
 
Implimpl ()
 
std::size_t localSize () const
 Return the number of entries in the local part of the vector.
 
std::size_t globalSize () const
 Return the number of entries in the global vector.
 
void resize (Basis const &basis)
 Resize the vector to the size of the basis.
 
void resizeZero (Basis const &basis)
 Resize the vector to the size of the basis and set to zero.
 
void init (Basis const &basis, bool clear)
 
void finish ()
 Finish the insertion of values, see init()
 
void insert (Index const &idx, typename Impl::value_type const &value, Assign assign={})
 Insert a single value into the matrix (add or overwrite to existing value) More...
 
void set (Index const &idx, typename Impl::value_type const &value)
 See insert for assignment operation Assigner::assign.
 
void add (Index const &idx, typename Impl::value_type const &value)
 See insert for assignment operation Assigner::plus_assign.
 
Impl::value_type get (Index const &idx) const
 Return the value of the vector at the local index idx.
 
void gather (LocalView const &localView, Node const &node, Buffer &buffer) const
 Extract values from the vector referring to the given local indices and store it into a buffer. More...
 
void gather (LocalView const &localView, Buffer &buffer) const
 
void gather (std::vector< Index > const &indices, Buffer &buffer)
 Call gather the values associated to the indices into the buffer.
 
void scatter (LocalView const &localView, Node const &node, NodeVector const &localVector, MaskRange const &mask, Assign assign)
 Insert a block of values into the vector (add or overwrite to existing values) More...
 
void scatter (LocalView const &localView, Node const &node, NodeVector const &localVector, Assign assign)
 Call scatter with MaskRange set to FakeContainer.
 
void scatter (LocalView const &localView, LocalVector const &localVector, Assign assign)
 Call scatter with Node given by the tree of the localView.
 
void scatter (std::vector< Index > const &indices, Buffer const &values)
 Call scatter the values associated to the indices into the vector.
 
void forEach (LocalInd const &localInd, Func &&func)
 Apply func to each value at given indices localInd. More...
 
void forEach (LocalInd const &localInd, Func &&func) const
 Apply func to each value at given indices localInd. More...
 

Static Public Member Functions

template<class Self >
static auto & coefficients (Self &self)
 

Protected Member Functions

void updateImpl (event::preAdapt e) override
 
void updateImpl (event::adapt e) override
 
void updateImpl (event::postAdapt) override
 

Friends

template<class , class , class , class >
class DiscreteFunction
 
template<class Range = void, class... Indices>
auto valueOf (DOFVector &dofVec, Indices... ii)
 A Generator to transform a DOFVector into a DiscreteFunction.
 
template<class Range = void, class... Indices>
auto valueOf (DOFVector const &dofVec, Indices... ii)
 A Generator to transform a const DOFVector into a DiscreteFunction.
 
template<class Type >
auto derivativeOf (DOFVector const &dofVec, Type type)
 
auto gradientOf (DOFVector const &dofVec)
 Return a derivative grid-function representing the gradient (jacobian) of the DOFVector.
 
auto divergenceOf (DOFVector const &dofVec)
 Return a derivative grid-function representing the divergence of the DOFVector.
 
auto partialDerivativeOf (DOFVector const &dofVec, std::size_t comp)
 Return a derivative grid-function representing the partial derivative wrt. the components comp of the DOFVector.
 

Detailed Description

template<class GB, class T = double, class Traits = BackendTraits>
class AMDiS::DOFVector< GB, T, Traits >

The basic container that stores a base vector and a corresponding basis.

Template Parameters
GBBasis of the vector
TType of the coefficients
TraitsCollection of parameter for the linear-algebra backend

Constructor & Destructor Documentation

◆ DOFVector() [1/3]

DOFVector ( std::shared_ptr< GB > const &  basis,
DataTransferTag  = {} 
)
inline

(1) Constructor. Stores the shared_ptr of the basis and creates a new DataTransfer.

◆ DOFVector() [2/3]

DOFVector ( GB_ &&  basis,
DataTransferTag  tag = {} 
)
inline

(2) Constructor. Forwards to (1) by wrapping reference into non-destroying shared_ptr, see Dune::wrap_or_move.

◆ DOFVector() [3/3]

DOFVector ( GV const &  gridView,
PBF const &  preBasisFactory,
DataTransferTag  tag = {} 
)
inline

(3) Constructor. Forwards to (1) by creating a new basis from a dune-functions pre-basis factory.

Member Function Documentation

◆ updateImpl() [1/3]

void updateImpl ( event::adapt  e)
inlineoverrideprotected

Override of Observer update(event::adapt) method. Redirects to adapt method of the DataTransfer object.

◆ updateImpl() [2/3]

void updateImpl ( event::postAdapt  )
inlineoverrideprotected

Override of Observer update(event::postAdapt) method. Redirects to postAdapt method of the DataTransfer object.

◆ updateImpl() [3/3]

void updateImpl ( event::preAdapt  e)
inlineoverrideprotected

Override of Observer update(event::preAdapt) method. Redirects to preAdapt method of the DataTransfer object.

Friends And Related Function Documentation

◆ derivativeOf

auto derivativeOf ( DOFVector< GB, T, Traits > const &  dofVec,
Type  type 
)
friend

Return a derivative grid-function of the given type, e.g., tag::jacobian, tag::divergence, or tag::partialDerivative


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