Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_BlockLowerTriInverseOp.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Teko: A package for block and physics based preconditioning
4 //
5 // Copyright 2010 NTESS and the Teko contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef __Teko_BlockLowerTriInverseOp_hpp__
11 #define __Teko_BlockLowerTriInverseOp_hpp__
12 
13 #include "Teko_Utilities.hpp"
14 #include "Teko_BlockImplicitLinearOp.hpp"
15 
16 namespace Teko {
17 
27  public:
29 
39  BlockLowerTriInverseOp(BlockedLinearOp &L, const std::vector<LinearOp> &invDiag);
40 
42 
43 
45  virtual VectorSpace range() const { return productRange_; }
46 
48  virtual VectorSpace domain() const { return productDomain_; }
49 
62  virtual void implicitApply(const BlockedMultiVector &x, BlockedMultiVector &y,
63  const double alpha = 1.0, const double beta = 0.0) const;
65 
66  virtual void describe(Teuchos::FancyOStream &out_arg,
67  const Teuchos::EVerbosityLevel verbLevel) const;
68 
69  protected:
70  // fundamental operators to use
71  const BlockedLinearOp L_;
72  std::vector<LinearOp> invDiag_;
73 
74  Teuchos::RCP<const Thyra::ProductVectorSpaceBase<double> >
76  Teuchos::RCP<const Thyra::ProductVectorSpaceBase<double> >
78 
79  // scratch space...so we don't have to reallocate
80  mutable BlockedMultiVector srcScrap_;
81  mutable BlockedMultiVector dstScrap_;
82  mutable bool allocated = false;
83 
84  private:
85  // hide me!
88 };
89 
90 inline LinearOp createBlockLowerTriInverseOp(BlockedLinearOp &U,
91  const std::vector<LinearOp> &invDiag) {
92  return Teuchos::rcp(new BlockLowerTriInverseOp(U, invDiag));
93 }
94 
95 inline LinearOp createBlockLowerTriInverseOp(BlockedLinearOp &L,
96  const std::vector<LinearOp> &invDiag,
97  const std::string &str) {
98  Teuchos::RCP<Thyra::LinearOpBase<double> > result =
99  Teuchos::rcp(new BlockLowerTriInverseOp(L, invDiag));
100  result->setObjectLabel(str);
101 
102  return result;
103 }
104 
105 } // end namespace Teko
106 
107 #endif
virtual void implicitApply(const BlockedMultiVector &x, BlockedMultiVector &y, const double alpha=1.0, const double beta=0.0) const =0
Perform a matrix vector multiply with this implicitly defined blocked operator.
virtual VectorSpace domain() const
Domain space of this operator.
const BlockedLinearOp L_
operator
Teuchos::RCP< const Thyra::ProductVectorSpaceBase< double > > productDomain_
Domain vector space.
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.
std::vector< LinearOp > invDiag_
(Approximate) Inverses of the diagonal operators
A virtual class that simplifies the construction of custom operators.
This linear operator computes the inverse of a lower triangular matrix.
Teuchos::RCP< const Thyra::ProductVectorSpaceBase< double > > productRange_
Range vector space.
virtual VectorSpace range() const
Range space of this operator.