DUNE PDELab (2.8)

Dune::PDELab::WeightedVectorAccumulationView< C > Class Template Reference

An accumulate-only view on a local vector that automatically takes into account an accumulation weight. More...

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

Public Types

typedef C Container
 The type of the underlying LocalVector.
 
typedef Container::BaseContainer BaseContainer
 The type of the storage container underlying the LocalVector.
 
typedef Container::value_type value_type
 The value type of the entries.
 
typedef Container::weight_type weight_type
 The type of the weight applied when accumulating contributions.
 
typedef WeightedVectorAccumulationView WeightedAccumulationView
 Export this type for uniform handling of the containers themselves and their views.
 
typedef Container::size_type size_type
 The size_type of the underlying container.
 

Public Member Functions

WeightedAccumulationView weightedAccumulationView (weight_type weight)
 Returns a WeighedAccumulationView with some weight in addition to this view's weight.
 
weight_type weight () const
 Returns the weight associated with this view. More...
 
void setWeight (weight_type weight)
 Resets the weighting coefficient of the view. More...
 
template<typename LFS >
void accumulate (const LFS &lfs, size_type n, value_type v)
 Applies the current weight to v and adds the result to the n-th degree of freedom of the lfs.
 
template<typename LFS >
void rawAccumulate (const LFS &lfs, size_type n, value_type v)
 Adds v to the n-th degree of freedom of the lfs without applying the current weight. More...
 
 WeightedVectorAccumulationView (C &container, weight_type weight)
 Constructor.
 
size_type size () const
 Returns the size of the underlying container.
 
bool modified () const
 Returns whether this view has been written to.
 
void resetModified ()
 Resets the modification state of the view to not modified. More...
 
Containercontainer ()
 Returns the container (of type LocalVector) that this view is based on.
 
const Containercontainer () const
 Returns the container (of type LocalVector) that this view is based on (const version).
 
BaseContainerbase ()
 Returns the storage container of the underlying LocalVector.
 
const BaseContainerbase () const
 Returns the storage container of the underlying LocalVector (const version).
 
auto data ()
 Access underlying container.
 
const auto data () const
 Access underlying container, const version.
 

Detailed Description

template<typename C>
class Dune::PDELab::WeightedVectorAccumulationView< C >

An accumulate-only view on a local vector that automatically takes into account an accumulation weight.

Member Function Documentation

◆ rawAccumulate()

template<typename C >
template<typename LFS >
void Dune::PDELab::WeightedVectorAccumulationView< C >::rawAccumulate ( const LFS &  lfs,
size_type  n,
value_type  v 
)
inline

Adds v to the n-th degree of freedom of the lfs without applying the current weight.

Warning
When using this method, you must take care of applying the weight yourself, otherwise your program may exhibit strange behavior or even calculate wrong results!

◆ resetModified()

template<typename C >
void Dune::PDELab::WeightedVectorAccumulationView< C >::resetModified ( )
inline

Resets the modification state of the view to not modified.

Warning
Never call this method from within a local operator, or your local residual / matrix contributions will be lost!

◆ setWeight()

template<typename C >
void Dune::PDELab::WeightedVectorAccumulationView< C >::setWeight ( weight_type  weight)
inline

Resets the weighting coefficient of the view.

Warning
Only call this method when you know what you are doing! It is especially not meant to be called from inside local operators.

References Dune::PDELab::WeightedVectorAccumulationView< C >::weight().

◆ weight()

template<typename C >
weight_type Dune::PDELab::WeightedVectorAccumulationView< C >::weight ( ) const
inline

Returns the weight associated with this view.

Note
This can be used together with rawAccumulate() to avoid applying the weight at each loop iteration.

Referenced by Dune::PDELab::WeightedVectorAccumulationView< C >::setWeight(), and Dune::PDELab::WeightedVectorAccumulationView< C >::weightedAccumulationView().


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 (Dec 21, 23:30, 2024)