KASKADE 7 development version
|
An adaptor presenting a Dune::LinearOperator <domain_type,range_type> interface to a contiguous sub-blockmatrix of an assembled heterogeneous Galerkin operator. More...
#include <istlinterface.hh>
An adaptor presenting a Dune::LinearOperator <domain_type,range_type> interface to a contiguous sub-blockmatrix of an assembled heterogeneous Galerkin operator.
Assembler | the VariationalFunctionalAssembler containing the Galerkin matrix |
firstRow | the first block row to be included |
lastRow | one behind the final row to be included (half-open range) |
firstCol | the first block column to be included |
lastCol | one behind the final column to be included (half-open range) |
BlockFilter | a boost::fusion predicate working on the assembler's matrix blocks, defining which blocks to include. The default is to take all blocks in the range given by [firstRow,lastRow[ x [firstCol,lastCol[. |
symmetric | whether the resulting operator is symmetric. The default is a conservative check whether the operator is symmetric or not. |
Due to data sharing, the GalerkinOperator is still valid even if the referenced assembler is deleted or otherwise invalidated. However, reassembling without invalidation of the assembler's data structures also modifies the values in the GalerkinOperator.
Definition at line 756 of file istlinterface.hh.
Classes | |
struct | ApplyScaleAdd |
struct | TransposedApplyScaleAdd |
Public Types | |
typedef IstlInterfaceDetail::Base< Assembler_, firstRow, lastRow, firstCol, lastCol, symmetric >::type | Base |
typedef Assembler_ | Assembler |
using | Domain = typename Assembler::AnsatzVariableSetDescription::template CoefficientVector< firstCol, lastCol > |
using | Range = typename Assembler::TestVariableSetDescription::template CoefficientVector< firstRow, lastRow > |
typedef Assembler::Scalar | Scalar |
typedef IstlInterfaceDetail::BlockInfo< firstRow, lastRow, firstCol, lastCol > | BlockInfo |
Public Member Functions | |
AssembledGalerkinOperator (Assembler const &assembler_, bool onlyLowerTriangle_=false, int nThreads_=0) | |
Constructor. More... | |
virtual | ~AssembledGalerkinOperator () |
void | update () |
update operator if grid has changed or assemble(...) has been called. More... | |
virtual MatrixAsTriplet< Scalar > | getTriplet () const __attribute__((deprecated)) |
DEPRECATED. Use get<MatrixAsTriplet<Scalar>>() instead. More... | |
template<class Matrix > | |
Matrix | get (int rowFirst=0, int rowLast=lastRow-firstRow, int colFirst=0, int colLast=lastCol-firstCol) const |
Access matrix (or a submatrix). More... | |
template<class Matrix > | |
std::unique_ptr< Matrix > | getPointer () const |
Access matrix via unique_ptr. Use if Matrix does not support move-semantics. More... | |
MatrixAsTriplet< Scalar > | getmat () const __attribute__((deprecated)) |
DEPRECATED. Use get<MatrixAsTriplet<Scalar>>() instead. More... | |
template<class Vector > | |
void | rangeToVector (Range const &y, Vector &coefficients) const |
Get coefficient vector \(coefficients\in\mathbb{K}^m\) from \(y\in Y\). 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... | |
Assembler const & | getAssembler () const |
Provides access to the underlying assembler. More... | |
template<int i> | |
auto const & | ansatzSpace () const |
Returns a reference to the FE space of the i-th variable in this GOP. More... | |
bool | getOnlyLowerTriangle () const |
virtual Dune::SolverCategory::Category | category () const override |
returns the category of the operator More... | |
Application of the operator | |
virtual void | apply (Domain const &x, Range &b) const |
compute \( b \leftarrow Ax \) More... | |
void | applyTransposed (Range const &x, Domain &b) const |
compute \( b \leftarrow A^T x\) 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 \leftarrow b+\alpha Ax\). More... | |
virtual void | applyscaleaddTransposed (Scalar alpha, Range const &x, Domain &b) const |
Compute \(b \leftarrow b+\alpha A^T x\). More... | |
Protected Types | |
typedef boost::fusion::result_of::as_vector< typenameIstlInterfaceDetail::AllBlocks< Assembler >::result >::type | allblocks |
typedef boost::fusion::result_of::as_vector< typenameboost::fusion::result_of::filter_if< allblocks, BlockFilter >::type >::type | filtered |
Protected Attributes | |
bool | onlyLowerTriangle |
Assembler const & | assembler |
boost::fusion::result_of::as_vector< typenameboost::fusion::result_of::transform< filtered, IstlInterfaceDetail::Translate< firstRow, firstCol > >::type >::type | blocks |
bool | isSymmetric |
int | nThreads |
Related Functions | |
(Note that these are not member functions.) | |
template<class Assembler , int firstRow = 0, int lastRow = Assembler::TestVariableSetDescription::noOfVariables, int firstCol = 0, int lastCol = Assembler::AnsatzVariableSetDescription::noOfVariables> | |
AssembledGalerkinOperator< Assembler, firstRow, lastRow, firstCol, lastCol > | makeAssembledGalerkinOperator (Assembler const &assembler) |
Convenience routine for creating assembled operators. More... | |
|
protected |
Definition at line 994 of file istlinterface.hh.
typedef Assembler_ Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::Assembler |
Definition at line 761 of file istlinterface.hh.
typedef IstlInterfaceDetail::Base<Assembler_,firstRow,lastRow,firstCol,lastCol,symmetric>::type Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::Base |
Definition at line 760 of file istlinterface.hh.
typedef IstlInterfaceDetail::BlockInfo<firstRow,lastRow,firstCol,lastCol> Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::BlockInfo |
Definition at line 766 of file istlinterface.hh.
using Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::Domain = typename Assembler::AnsatzVariableSetDescription::template CoefficientVector<firstCol,lastCol> |
Definition at line 762 of file istlinterface.hh.
|
protected |
Definition at line 998 of file istlinterface.hh.
using Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::Range = typename Assembler::TestVariableSetDescription::template CoefficientVector<firstRow,lastRow> |
Definition at line 763 of file istlinterface.hh.
typedef Assembler::Scalar Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::Scalar |
Definition at line 764 of file istlinterface.hh.
|
inlineexplicit |
Constructor.
op | the variational assembler. |
onlyLowerTriangle | if true, only the lower triangle is returned if a matrix is requested |
nThreads | use given number of threads (or a machine-dependent default if nThreads<=0). |
Definition at line 775 of file istlinterface.hh.
|
inlinevirtual |
Definition at line 784 of file istlinterface.hh.
|
inline |
Returns a reference to the FE space of the i-th variable in this GOP.
Definition at line 971 of file istlinterface.hh.
|
inlinevirtual |
compute \( b \leftarrow Ax \)
Definition at line 899 of file istlinterface.hh.
Referenced by Kaskade::ddxpy_template(), Kaskade::GoalOrientedErrorEstimator< TemplateFunctional, OriginalVariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace >::operator()(), and Kaskade::MartinsErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSystemSolver_H, LinearSystemSolver_L, LinearSolverA_H, LinearSolverA_L, RefinementStrategy >::operator()().
|
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 919 of file istlinterface.hh.
|
inlinevirtual |
Compute \(b \leftarrow b+\alpha Ax\).
Note that x and b must not refer to the same memory locations (in case Domain==Range).
Definition at line 938 of file istlinterface.hh.
Referenced by Kaskade::YetAnotherHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components, ReferenceSolution, ReferenceOperator >::operator()(), Kaskade::YetAnotherHBErrorEstimator_Elasticity< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, RefinementStrategy, lump, components >::operator()(), Kaskade::GoalOrientedErrorEstimator< TemplateFunctional, OriginalVariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace >::operator()(), Kaskade::HierarchicalBasisErrorEstimator2< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverLA, LinearSolverHA, LinearSolverLM, LinearSolverHM, lumpM, RefinementStrategy >::operator()(), Kaskade::StateEquationHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverHA, RefinementStrategy >::operator()(), Kaskade::VariationalEquationHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverHA, RefinementStrategy >::operator()(), Kaskade::AdjointEquationHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverHA, RefinementStrategy >::operator()(), Kaskade::AdjointEquationLinearPropagationHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverLA, LinearSolverHA, LinearSolverLU, RefinementStrategy >::operator()(), Kaskade::AnotherHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverLA, LinearSolverHA, LinearSolverHU, LinearSolverLU, RefinementStrategy, lump >::operator()(), Kaskade::StupidHBErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverHA, RefinementStrategy >::operator()(), Kaskade::MartinsErrorEstimator< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSystemSolver_H, LinearSystemSolver_L, LinearSolverA_H, LinearSolverA_L, RefinementStrategy >::operator()(), and Kaskade::HierarchicalBasisErrorEstimator3< Functional, VariableSetDescription, ExtensionVariableSetDescription, ExtensionSpace, NormFunctional, LinearSolverLA, LinearSolverHA, LinearSolverLM, LinearSolverHM, lumpM, RefinementStrategy >::operator()().
|
inlinevirtual |
Compute \(b \leftarrow b+\alpha A^T x\).
Note that x and b must not refer to the same memory locations (in case Domain==Range).
Definition at line 949 of file istlinterface.hh.
|
inline |
compute \( b \leftarrow A^T x\)
Definition at line 908 of file istlinterface.hh.
|
inlineoverridevirtual |
returns the category of the operator
From the Dune doxygen documentation it is unclear what this is supposed to mean. We return a dummy here.
Definition at line 986 of file istlinterface.hh.
|
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:
Definition at line 858 of file istlinterface.hh.
|
inlinevirtual |
Definition at line 927 of file istlinterface.hh.
Referenced by Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyDp().
|
inline |
Access matrix (or a submatrix).
Matrix | a matrix type for which a partial specialization of AssemblyDetail::Fill is defined. |
rowFirst | index of first row (starting at default 0) |
rowLast | index of one behind last row (half-open range, defaults to maximum range) |
colFirst | index of first column (starting at default 0) |
colLast | index of one behind last column (half-open range, defaults to maximum range) |
Use this if Matrix supports move semantics. (NumaBCRSMatrix does.)
Definition at line 814 of file istlinterface.hh.
|
inline |
Provides access to the underlying assembler.
Definition at line 962 of file istlinterface.hh.
MatrixAsTriplet< Scalar > Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::getmat | ( | ) | const |
DEPRECATED. Use get<MatrixAsTriplet<Scalar>>() instead.
|
inline |
Definition at line 979 of file istlinterface.hh.
|
inline |
Access matrix via unique_ptr. Use if Matrix does not support move-semantics.
Definition at line 821 of file istlinterface.hh.
|
inlinevirtual |
DEPRECATED. Use get<MatrixAsTriplet<Scalar>>() instead.
Definition at line 797 of file istlinterface.hh.
|
inline |
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:
Definition at line 842 of file istlinterface.hh.
|
inline |
update operator if grid has changed or assemble(...) has been called.
Definition at line 787 of file istlinterface.hh.
|
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:
Definition at line 873 of file istlinterface.hh.
|
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:
Definition at line 888 of file istlinterface.hh.
|
related |
Convenience routine for creating assembled operators.
Definition at line 1111 of file istlinterface.hh.
|
protected |
Definition at line 1001 of file istlinterface.hh.
Referenced by Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::ansatzSpace(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::apply(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyDp(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyscaleadd(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyscaleaddTransposed(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyTransposed(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::get(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::getAssembler(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::getPointer(), and Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::update().
|
protected |
Definition at line 1006 of file istlinterface.hh.
Referenced by Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::apply(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyDp(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyscaleadd(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyscaleaddTransposed(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::applyTransposed(), and Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::update().
|
protected |
Definition at line 1007 of file istlinterface.hh.
|
protected |
Definition at line 1100 of file istlinterface.hh.
Referenced by Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::update().
|
protected |
Definition at line 1000 of file istlinterface.hh.
Referenced by Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::get(), Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::getOnlyLowerTriangle(), and Kaskade::AssembledGalerkinOperator< Assembler_, firstRow, lastRow, firstCol, lastCol, BlockFilter, symmetric >::getPointer().