Teko  Version of the Day
 All Classes Files Functions Variables Pages
Public Member Functions | Protected Attributes | List of all members
Teko::BlockDiagonalInverseOp Class Reference

This linear operator computes the inverse of a diagonal matrix. More...

#include <Teko_BlockDiagonalInverseOp.hpp>

Inheritance diagram for Teko::BlockDiagonalInverseOp:
Inheritance graph
[legend]

Public Member Functions

 BlockDiagonalInverseOp (BlockedLinearOp &A, const std::vector< LinearOp > &invDiag)
 This constructor explicitly takes a diagonal matrix and inverse diagonal operators and builds a back substitution operator. More...
 
Inherited methods from Thyra::LinearOpBase
virtual VectorSpace range () const
 Range space of this operator. More...
 
virtual VectorSpace domain () const
 Domain space of this operator. More...
 
virtual void implicitApply (const BlockedMultiVector &x, BlockedMultiVector &y, const double alpha=1.0, const double beta=0.0) const
 Perform a matrix vector multiply with this operator. More...
 
virtual void implicitApply (const Thyra::EOpTransp M_trans, const BlockedMultiVector &x, BlockedMultiVector &y, const double alpha=1.0, const double beta=0.0) const
 Perform a matrix vector multiply with this implicitly defined blocked operator. More...
 

Protected Attributes

std::vector< LinearOp > invDiag_
 (Approximate) Inverses of the diagonal operators More...
 
Teuchos::RCP< const
Thyra::ProductVectorSpaceBase
< double > > 
productRange_
 Range vector space. More...
 
Teuchos::RCP< const
Thyra::ProductVectorSpaceBase
< double > > 
productDomain_
 Domain vector space. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Teko::BlockImplicitLinearOp
virtual bool opSupportedImpl (const Thyra::EOpTransp M_trans) const
 Functions required by Thyra::LinearOpBase. More...
 

Detailed Description

This linear operator computes the inverse of a diagonal matrix.

This linear operator computes the inverse of a diagonal matrix. This requires the inverse of the operators on the diagonal.

Definition at line 25 of file Teko_BlockDiagonalInverseOp.hpp.

Constructor & Destructor Documentation

Teko::BlockDiagonalInverseOp::BlockDiagonalInverseOp ( BlockedLinearOp &  A,
const std::vector< LinearOp > &  invDiag 
)

This constructor explicitly takes a diagonal matrix and inverse diagonal operators and builds a back substitution operator.

This constructor explicitly takes the inverse diagonal operators and builds a back substitution operator.

Parameters
[in]AMatrix object
[in]invDiagVector containing the inverse of the diagonal blocks

Definition at line 18 of file Teko_BlockDiagonalInverseOp.cpp.

Member Function Documentation

virtual VectorSpace Teko::BlockDiagonalInverseOp::range ( ) const
inlinevirtual

Range space of this operator.

Implements Teko::BlockImplicitLinearOp.

Definition at line 45 of file Teko_BlockDiagonalInverseOp.hpp.

virtual VectorSpace Teko::BlockDiagonalInverseOp::domain ( ) const
inlinevirtual

Domain space of this operator.

Implements Teko::BlockImplicitLinearOp.

Definition at line 48 of file Teko_BlockDiagonalInverseOp.hpp.

void Teko::BlockDiagonalInverseOp::implicitApply ( const BlockedMultiVector &  x,
BlockedMultiVector &  y,
const double  alpha = 1.0,
const double  beta = 0.0 
) const
virtual

Perform a matrix vector multiply with this operator.

The apply function takes one vector as input and applies $ D^{-1} $. The result is returned in $y$. If this operator is represented as $M$ then $ y = \alpha M x + \beta y $.

Parameters
[in]x
[in,out]y
[in]alpha(default=1)
[in]beta(default=0)

Implements Teko::BlockImplicitLinearOp.

Definition at line 35 of file Teko_BlockDiagonalInverseOp.cpp.

void Teko::BlockDiagonalInverseOp::implicitApply ( const Thyra::EOpTransp  M_trans,
const BlockedMultiVector &  x,
BlockedMultiVector &  y,
const double  alpha = 1.0,
const double  beta = 0.0 
) const
virtual

Perform a matrix vector multiply with this implicitly defined blocked operator.

The apply function takes one vector as input and applies a linear operator. The result is returned in $y$. If this operator is represented as $M$ then $ y = \alpha M x + \beta y $

Parameters
[in]x
[in,out]y
[in]alpha(default=1)
[in]beta(default=0)

Reimplemented from Teko::BlockImplicitLinearOp.

Definition at line 41 of file Teko_BlockDiagonalInverseOp.cpp.

Member Data Documentation

std::vector<LinearOp> Teko::BlockDiagonalInverseOp::invDiag_
protected

(Approximate) Inverses of the diagonal operators

Definition at line 88 of file Teko_BlockDiagonalInverseOp.hpp.

Teuchos::RCP<const Thyra::ProductVectorSpaceBase<double> > Teko::BlockDiagonalInverseOp::productRange_
protected

Range vector space.

Definition at line 91 of file Teko_BlockDiagonalInverseOp.hpp.

Teuchos::RCP<const Thyra::ProductVectorSpaceBase<double> > Teko::BlockDiagonalInverseOp::productDomain_
protected

Domain vector space.

Definition at line 93 of file Teko_BlockDiagonalInverseOp.hpp.


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