Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_StridedTpetraOperator.cpp
1 /*
2 // @HEADER
3 //
4 // ***********************************************************************
5 //
6 // Teko: A package for block and physics based preconditioning
7 // Copyright 2010 Sandia Corporation
8 //
9 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
10 // the U.S. Government retains certain rights in this software.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
40 //
41 // ***********************************************************************
42 //
43 // @HEADER
44 
45 */
46 
47 #include "Teko_StridedTpetraOperator.hpp"
48 #include "Teko_TpetraStridedMappingStrategy.hpp"
49 #include "Teko_TpetraReorderedMappingStrategy.hpp"
50 
51 #include "Teuchos_VerboseObject.hpp"
52 
53 #include "Thyra_LinearOpBase.hpp"
54 #include "Thyra_TpetraLinearOp.hpp"
55 #include "Thyra_TpetraThyraWrappers.hpp"
56 #include "Thyra_DefaultProductMultiVector.hpp"
57 #include "Thyra_DefaultProductVectorSpace.hpp"
58 #include "Thyra_DefaultBlockedLinearOp.hpp"
59 
60 #include "Tpetra_Vector.hpp"
61 #include "MatrixMarket_Tpetra.hpp"
62 
63 #include "Teko_Utilities.hpp"
64 
65 namespace Teko {
66 namespace TpetraHelpers {
67 
68 using Teuchos::RCP;
69 using Teuchos::rcp;
70 using Teuchos::rcp_dynamic_cast;
71 
72 StridedTpetraOperator::StridedTpetraOperator(
73  int numVars, const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> >& content,
74  const std::string& label)
75  : Teko::TpetraHelpers::TpetraOperatorWrapper(), label_(label) {
76  std::vector<int> vars;
77 
78  // build vector describing the sub maps
79  for (int i = 0; i < numVars; i++) vars.push_back(1);
80 
81  SetContent(vars, content);
82 }
83 
84 StridedTpetraOperator::StridedTpetraOperator(
85  const std::vector<int>& vars,
86  const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> >& content, const std::string& label)
87  : Teko::TpetraHelpers::TpetraOperatorWrapper(), label_(label) {
88  SetContent(vars, content);
89 }
90 
91 void StridedTpetraOperator::SetContent(
92  const std::vector<int>& vars,
93  const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> >& content) {
94  fullContent_ = content;
95  stridedMapping_ = rcp(new TpetraStridedMappingStrategy(vars, fullContent_->getDomainMap(),
96  *fullContent_->getDomainMap()->getComm()));
97  SetMapStrategy(stridedMapping_);
98 
99  // build thyra operator
100  BuildBlockedOperator();
101 }
102 
103 void StridedTpetraOperator::BuildBlockedOperator() {
104  TEUCHOS_ASSERT(stridedMapping_ != Teuchos::null);
105 
106  // get a CRS matrix
107  const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > crsContent =
108  rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT> >(fullContent_);
109 
110  // ask the strategy to build the Thyra operator for you
111  if (stridedOperator_ == Teuchos::null) {
112  stridedOperator_ = stridedMapping_->buildBlockedThyraOp(crsContent, label_);
113  } else {
114  const RCP<Thyra::BlockedLinearOpBase<ST> > blkOp =
115  rcp_dynamic_cast<Thyra::BlockedLinearOpBase<ST> >(stridedOperator_, true);
116  stridedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
117  }
118 
119  // set whatever is returned
120  SetOperator(stridedOperator_, false);
121 
122  // reorder if neccessary
123  if (reorderManager_ != Teuchos::null) Reorder(*reorderManager_);
124 }
125 
126 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > StridedTpetraOperator::GetBlock(
127  int i, int j) const {
128  const RCP<const Thyra::BlockedLinearOpBase<ST> > blkOp =
129  Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<ST> >(getThyraOp());
130 
131  RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT> > tOp =
132  rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT> >(blkOp->getBlock(i, j), true);
133  return tOp->getConstTpetraOperator();
134 }
135 
139 void StridedTpetraOperator::Reorder(const BlockReorderManager& brm) {
140  reorderManager_ = rcp(new BlockReorderManager(brm));
141 
142  // build reordered objects
143  RCP<const MappingStrategy> reorderMapping =
144  rcp(new TpetraReorderedMappingStrategy(*reorderManager_, stridedMapping_));
145  RCP<const Thyra::BlockedLinearOpBase<ST> > blockOp =
146  rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<ST> >(stridedOperator_);
147 
148  RCP<const Thyra::LinearOpBase<ST> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
149 
150  // set them as working values
151  SetMapStrategy(reorderMapping);
152  SetOperator(A, false);
153 }
154 
156 void StridedTpetraOperator::RemoveReording() {
157  SetMapStrategy(stridedMapping_);
158  SetOperator(stridedOperator_, false);
159  reorderManager_ = Teuchos::null;
160 }
161 
164 void StridedTpetraOperator::WriteBlocks(const std::string& prefix) const {
165  RCP<Thyra::PhysicallyBlockedLinearOpBase<ST> > blockOp =
166  rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<ST> >(stridedOperator_);
167 
168  // get size of strided block operator
169  int rows = Teko::blockRowCount(blockOp);
170 
171  for (int i = 0; i < rows; i++) {
172  for (int j = 0; j < rows; j++) {
173  // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered)
174  RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT> > tOp =
175  rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT> >(blockOp->getBlock(i, j));
176  RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > mat =
177  Teuchos::rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT> >(
178  tOp->getConstTpetraOperator());
179 
180  // write to file
181  Tpetra::MatrixMarket::Writer<Tpetra::CrsMatrix<ST, LO, GO, NT> >::writeSparseFile(
182  formatBlockName(prefix, i, j, rows).c_str(), mat);
183  }
184  }
185 }
186 
194 std::string StridedTpetraOperator::PrintNorm(const eNormType& nrmType, const char newline) {
195  BlockedLinearOp blockOp = toBlockedLinearOp(stridedOperator_);
196 
197  // get size of strided block operator
198  int rowCount = Teko::blockRowCount(blockOp);
199  int colCount = Teko::blockRowCount(blockOp);
200 
201  std::stringstream ss;
202  ss.precision(4);
203  ss << std::scientific;
204  for (int row = 0; row < rowCount; row++) {
205  for (int col = 0; col < colCount; col++) {
206  // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered)
207  RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT> > Aij =
208  rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT> >(
209  blockOp->getBlock(row, col), true);
210  RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > mat =
211  Teuchos::rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT> >(
212  Aij->getConstTpetraOperator(), true);
213 
214  // compute the norm
215  ST norm = 0.0;
216  switch (nrmType) {
217  // case Inf:
218  // norm = mat->normInf();
219  // break;
220  // case One:
221  // norm = mat->normOne();
222  // break;
223  case Frobenius: norm = mat->getFrobeniusNorm(); break;
224  default:
225  TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error,
226  "StridedTpetraOperator::NormType incorrectly specified. Only "
227  "Frobenius norm implemented for Tpetra matrices.");
228  }
229 
230  ss << norm << " ";
231  }
232  ss << newline;
233  }
234 
235  return ss.str();
236 }
237 
238 bool StridedTpetraOperator::testAgainstFullOperator(int count, ST tol) const {
239  Tpetra::Vector<ST, LO, GO, NT> xf(getRangeMap());
240  Tpetra::Vector<ST, LO, GO, NT> xs(getRangeMap());
241  Tpetra::Vector<ST, LO, GO, NT> y(getDomainMap());
242 
243  // test operator many times
244  bool result = true;
245  ST diffNorm = 0.0, trueNorm = 0.0;
246  for (int i = 0; i < count; i++) {
247  xf.putScalar(0.0);
248  xs.putScalar(0.0);
249  y.randomize();
250 
251  // apply operator
252  apply(y, xs); // xs = A*y
253  fullContent_->apply(y, xf); // xf = A*y
254 
255  // compute norms
256  xs.update(-1.0, xf, 1.0);
257  diffNorm = xs.norm2();
258  trueNorm = xf.norm2();
259 
260  // check result
261  result &= (diffNorm / trueNorm < tol);
262  }
263 
264  return result;
265 }
266 
267 } // end namespace TpetraHelpers
268 } // end namespace Teko