10 #include "Teko_BlockUpperTriInverseOp.hpp"
12 #include "Teuchos_Utils.hpp"
28 BlockUpperTriInverseOp::BlockUpperTriInverseOp(BlockedLinearOp& U,
29 const std::vector<LinearOp>& invDiag)
34 int blocks = blockRowCount(
U_);
35 TEUCHOS_ASSERT(blocks > 0);
36 TEUCHOS_ASSERT(blocks == blockColCount(
U_));
37 TEUCHOS_ASSERT(blocks == (
int)
invDiag_.size());
48 const double alpha,
const double beta)
const {
66 const BlockedMultiVector& src, BlockedMultiVector& dst,
67 const double alpha,
const double beta)
const {
68 int blocks = blockCount(src);
70 TEUCHOS_ASSERT(blocks == blockRowCount(
U_));
71 TEUCHOS_ASSERT(blocks == blockCount(dst));
74 srcScrap_ = deepcopy(src);
75 dstScrap_ = deepcopy(dst);
80 Thyra::assign<double>(srcScrap_.ptr(), *src);
81 BlockedMultiVector dstCopy;
83 Thyra::assign<double>(dstScrap_.ptr(), *dst);
90 std::vector<MultiVector> dstVec;
91 std::vector<MultiVector> scrapVec;
92 for (
int b = 0; b < blocks; b++) {
93 dstVec.push_back(getBlock(b, dstCopy));
94 scrapVec.push_back(getBlock(b, srcScrap_));
99 if (M_trans == Thyra::NOTRANS) {
100 for (
int b = blocks - 1; b >= 0; b--) {
101 applyOp(
invDiag_[b], scrapVec[b], dstVec[b]);
104 for (
int i = 0; i < b; i++) {
105 LinearOp u_ib = getBlock(i, b,
U_);
106 if (u_ib != Teuchos::null) {
107 applyOp(u_ib, dstVec[b], scrapVec[i], -1.0, 1.0);
111 }
else if (M_trans == Thyra::TRANS || M_trans == Thyra::CONJTRANS) {
112 for (
int b = 0; b < blocks; b++) {
113 applyTransposeOp(
invDiag_[b], scrapVec[b], dstVec[b]);
116 for (
int i = b + 1; i < blocks; i++) {
117 LinearOp u_bi = getBlock(b, i,
U_);
118 if (u_bi != Teuchos::null) {
119 applyTransposeOp(u_bi, dstVec[b], scrapVec[i], -1.0, 1.0);
124 TEUCHOS_TEST_FOR_EXCEPT(
true);
129 update(alpha, dstCopy, beta, dst);
130 else if (alpha != 1.0)
134 void BlockUpperTriInverseOp::describe(Teuchos::FancyOStream& out_arg,
135 const Teuchos::EVerbosityLevel verbLevel)
const {
136 using Teuchos::OSTab;
138 RCP<Teuchos::FancyOStream> out = rcp(&out_arg,
false);
141 case Teuchos::VERB_DEFAULT:
142 case Teuchos::VERB_LOW: *out << this->description() << std::endl;
break;
143 case Teuchos::VERB_MEDIUM:
144 case Teuchos::VERB_HIGH:
145 case Teuchos::VERB_EXTREME: {
146 *out << Teuchos::Describable::description() <<
"{"
147 <<
"rangeDim=" << this->
range()->dim() <<
",domainDim=" << this->
domain()->dim()
148 <<
",rows=" << blockRowCount(
U_) <<
",cols=" << blockColCount(
U_) <<
"}\n";
151 *out <<
"[U Operator] = ";
152 *out << Teuchos::describe(*
U_, verbLevel);
156 *out <<
"[invDiag Operators]:\n";
158 for (
int i = 0; i < blockRowCount(
U_); i++) {
159 *out <<
"[invD(" << i <<
")] = ";
160 *out << Teuchos::describe(*
invDiag_[i], verbLevel);
165 default: TEUCHOS_TEST_FOR_EXCEPT(
true);
virtual VectorSpace domain() const
Domain space of this operator.
std::vector< LinearOp > invDiag_
(Approximate) Inverses of the diagonal operators
virtual VectorSpace range() const
Range space of this operator.
Teuchos::RCP< const Thyra::ProductVectorSpaceBase< double > > productRange_
Range 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.
const BlockedLinearOp U_
operator
Teuchos::RCP< const Thyra::ProductVectorSpaceBase< double > > productDomain_
Domain vector space.