template<class NormalStepAssembler, class TangentialStepAssembler, int stateId = 1, int controlId = 0, int adjointId = 2>
class Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >
- Todo:
- docme
Definition at line 1124 of file istlinterface.hh.
|
| LagrangeOperator (NormalStepAssembler const &normalStepAssembler_, TangentialStepAssembler const &tangentialStepAssembler_, bool onlyLowerTriangle_=false, int nThreads_=0) |
|
virtual | ~LagrangeOperator () |
|
void | update () |
| update operator if grid has changed or assemble(...) has been called. More...
|
|
virtual MatrixAsTriplet< Scalar > | getTriplet () const |
|
template<class Vector > |
void | rangeToVector (Range const &y, Vector &coefficients) const |
| returns a reference to the matrix More...
|
|
template<class Vector > |
void | domainToVector (Domain const &x, Vector &coefficients) const |
| Get coefficients vector \(coefficients\in\mathbb{K}^n\) from \(x\in X\). More...
|
|
template<class Vector > |
void | vectorToDomain (Vector const &coefficients, Domain &x) const |
| Get \(x\in X\) from coefficients vector \(coefficients\in\mathbb{K}^n\). More...
|
|
template<class Vector > |
void | vectorToRange (Vector const &coefficients, Range &x) const |
| Get \( y\in Y \) from coefficient vector \(coefficients\in\mathbb{K}^m\) Apply \(S^{-1}_Y\) to \(coefficients\): \(x=S^{-1}_Y(y)\). More...
|
|
virtual void | apply (Domain const &x, Range &b) const |
| compute \( b \leftarrow Ax \) More...
|
|
virtual Scalar | applyDp (Domain const &x, Range &b) const |
| Computes \( b \leftarrow Ax \) and, in case \( A \) is symmetric, also \( \langle x, b \rangle \). More...
|
|
virtual Scalar | dp (Domain const &x, Range const &y) const |
|
virtual void | applyscaleadd (Scalar alpha, Domain const &x, Range &b) const |
| Compute \(b=b+\alpha Ax\) Note that x and b must not refer to the same memory locations (in case Domain==Range). More...
|
|
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
virtual Scalar Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::applyDp |
( |
Domain const & |
x, |
|
|
Range & |
b |
|
) |
| const |
|
inlinevirtual |
Computes \( b \leftarrow Ax \) and, in case \( A \) is symmetric, also \( \langle x, b \rangle \).
If \( A \) is not symmetric, zero is returned.
Definition at line 1277 of file istlinterface.hh.
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
virtual void Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::applyscaleadd |
( |
Scalar |
alpha, |
|
|
Domain const & |
x, |
|
|
Range & |
b |
|
) |
| const |
|
inlinevirtual |
Compute \(b=b+\alpha Ax\) Note that x and b must not refer to the same memory locations (in case Domain==Range).
Definition at line 1294 of file istlinterface.hh.
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
template<class Vector >
void Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::domainToVector |
( |
Domain const & |
x, |
|
|
Vector & |
coefficients |
|
) |
| const |
|
inline |
Get coefficients vector \(coefficients\in\mathbb{K}^n\) from \(x\in X\).
Apply \(S_X\) to \(x\in X\): \(coefficients=S_X(x)\).
The used vector type Vector must provide:
- its iterator type via typename Vector::iterator
- (possibly overloads of) the free functions:
- typename Vector::iterator std::begin(Vector&);
Definition at line 1229 of file istlinterface.hh.
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
template<class Vector >
void Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::rangeToVector |
( |
Range const & |
y, |
|
|
Vector & |
coefficients |
|
) |
| const |
|
inline |
returns a reference to the matrix
Get coefficient vector \(coefficients\in\mathbb{K}^m\) from \(y\in Y\). Apply \(S_Y\) to \(y\in Y\): \(coefficients=S_Y(y)\).
The used vector type Vector must provide:
- its iterator type via typename Vector::iterator
- (possibly overloads of) the free functions:
- typename Vector::iterator std::begin(Vector&);
Definition at line 1215 of file istlinterface.hh.
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
template<class Vector >
void Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::vectorToDomain |
( |
Vector const & |
coefficients, |
|
|
Domain & |
x |
|
) |
| const |
|
inline |
Get \(x\in X\) from coefficients vector \(coefficients\in\mathbb{K}^n\).
Apply \(S^{-1}_X\) to \(coefficients\): \(x=S^{-1}_X(x)\).
The used vector type Vector must provide:
- its iterator type via typename Vector::const_iterator
- (possibly overloads of) the free functions:
- typename Vector::const_iterator std::begin(Vector const&);
Definition at line 1244 of file istlinterface.hh.
template<class NormalStepAssembler , class TangentialStepAssembler , int stateId = 1, int controlId = 0, int adjointId = 2>
template<class Vector >
void Kaskade::LagrangeOperator< NormalStepAssembler, TangentialStepAssembler, stateId, controlId, adjointId >::vectorToRange |
( |
Vector const & |
coefficients, |
|
|
Range & |
x |
|
) |
| const |
|
inline |
Get \( y\in Y \) from coefficient vector \(coefficients\in\mathbb{K}^m\) Apply \(S^{-1}_Y\) to \(coefficients\): \(x=S^{-1}_Y(y)\).
The used vector type Vector must provide:
- its iterator type via typename Vector::const_iterator
- (possibly overloads of) the free functions:
- typename Vector::const_iterator std::begin(Vector const&);
Definition at line 1259 of file istlinterface.hh.