Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_BlockUpperTriInverseOp.cpp
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 #include "Teko_BlockUpperTriInverseOp.hpp"
11 
12 #include "Teuchos_Utils.hpp"
13 
14 namespace Teko {
15 
16 using Teuchos::RCP;
17 
28 BlockUpperTriInverseOp::BlockUpperTriInverseOp(BlockedLinearOp& U,
29  const std::vector<LinearOp>& invDiag)
30  : U_(U) {
31  invDiag_ = invDiag;
32 
33  // sanity check
34  int blocks = blockRowCount(U_);
35  TEUCHOS_ASSERT(blocks > 0);
36  TEUCHOS_ASSERT(blocks == blockColCount(U_));
37  TEUCHOS_ASSERT(blocks == (int)invDiag_.size());
38 
39  // create the range and product space
41 
42  // just flip flop them!
43  productRange_ = U->productDomain();
44  productDomain_ = U->productRange();
45 }
46 
47 void BlockUpperTriInverseOp::implicitApply(const BlockedMultiVector& src, BlockedMultiVector& dst,
48  const double alpha, const double beta) const {
49  // call the no tranpose version
50  implicitApply(Thyra::NOTRANS, src, dst, alpha, beta);
51 }
52 
65 void BlockUpperTriInverseOp::implicitApply(const Thyra::EOpTransp M_trans,
66  const BlockedMultiVector& src, BlockedMultiVector& dst,
67  const double alpha, const double beta) const {
68  int blocks = blockCount(src);
69 
70  TEUCHOS_ASSERT(blocks == blockRowCount(U_));
71  TEUCHOS_ASSERT(blocks == blockCount(dst));
72 
73  if (!allocated) {
74  srcScrap_ = deepcopy(src);
75  dstScrap_ = deepcopy(dst);
76  allocated = true;
77  }
78 
79  // build a scrap vector for storing work
80  Thyra::assign<double>(srcScrap_.ptr(), *src);
81  BlockedMultiVector dstCopy;
82  if (beta != 0.0) {
83  Thyra::assign<double>(dstScrap_.ptr(), *dst);
84  dstCopy = dstScrap_;
85  } else
86  dstCopy = dst; // shallow copy
87 
88  // extract the blocks components from
89  // the source and destination vectors
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_));
95  }
96 
97  // run back-substituion: run over each column
98  // From Heath pg. 66
99  if (M_trans == Thyra::NOTRANS) {
100  for (int b = blocks - 1; b >= 0; b--) {
101  applyOp(invDiag_[b], scrapVec[b], dstVec[b]);
102 
103  // loop over each row
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);
108  }
109  }
110  }
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]);
114 
115  // loop over each row
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);
120  }
121  }
122  }
123  } else {
124  TEUCHOS_TEST_FOR_EXCEPT(true);
125  }
126 
127  // scale result by alpha
128  if (beta != 0)
129  update(alpha, dstCopy, beta, dst); // dst = alpha * dstCopy + beta * dst
130  else if (alpha != 1.0)
131  scale(alpha, dst); // dst = alpha * dst
132 }
133 
134 void BlockUpperTriInverseOp::describe(Teuchos::FancyOStream& out_arg,
135  const Teuchos::EVerbosityLevel verbLevel) const {
136  using Teuchos::OSTab;
137 
138  RCP<Teuchos::FancyOStream> out = rcp(&out_arg, false);
139  OSTab tab(out);
140  switch (verbLevel) {
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";
149  {
150  OSTab tab2(out);
151  *out << "[U Operator] = ";
152  *out << Teuchos::describe(*U_, verbLevel);
153  }
154  {
155  OSTab tab2(out);
156  *out << "[invDiag Operators]:\n";
157  tab.incrTab();
158  for (int i = 0; i < blockRowCount(U_); i++) {
159  *out << "[invD(" << i << ")] = ";
160  *out << Teuchos::describe(*invDiag_[i], verbLevel);
161  }
162  }
163  break;
164  }
165  default: TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
166  }
167 }
168 
169 } // end namespace Teko
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.