Teko Version of the Day
Loading...
Searching...
No Matches
Teko_Utilities.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_Config.h"
11#include "Teko_Utilities.hpp"
12
13// Thyra includes
14#include "Thyra_MultiVectorStdOps.hpp"
15#include "Thyra_ZeroLinearOpBase.hpp"
16#include "Thyra_DefaultDiagonalLinearOp.hpp"
17#include "Thyra_DefaultAddedLinearOp.hpp"
18#include "Thyra_DefaultScaledAdjointLinearOp.hpp"
19#include "Thyra_DefaultMultipliedLinearOp.hpp"
20#include "Thyra_DefaultZeroLinearOp.hpp"
21#include "Thyra_DefaultProductMultiVector.hpp"
22#include "Thyra_DefaultProductVectorSpace.hpp"
23#include "Thyra_MultiVectorStdOps.hpp"
24#include "Thyra_VectorStdOps.hpp"
25#include "Thyra_SpmdVectorBase.hpp"
26#include <utility>
27
28#ifdef TEKO_HAVE_EPETRA
29#include "Thyra_EpetraExtDiagScaledMatProdTransformer.hpp"
30#include "Thyra_EpetraExtDiagScalingTransformer.hpp"
31#include "Thyra_EpetraExtAddTransformer.hpp"
32#include "Thyra_get_Epetra_Operator.hpp"
33#include "Thyra_EpetraThyraWrappers.hpp"
34#include "Thyra_EpetraOperatorWrapper.hpp"
35#include "Thyra_EpetraLinearOp.hpp"
36#endif
37
38// Teuchos includes
39#include "Teuchos_Array.hpp"
40
41// Epetra includes
42#ifdef TEKO_HAVE_EPETRA
43#include "Epetra_Operator.h"
44#include "Epetra_CrsGraph.h"
45#include "Epetra_CrsMatrix.h"
46#include "Epetra_Vector.h"
47#include "Epetra_Map.h"
48
49#include "EpetraExt_Transpose_RowMatrix.h"
50#include "EpetraExt_MatrixMatrix.h"
51#include <EpetraExt_BlockMapOut.h>
52#include <EpetraExt_RowMatrixOut.h>
53
54#include "Teko_EpetraHelpers.hpp"
55#include "Teko_EpetraOperatorWrapper.hpp"
56#endif
57
58// Anasazi includes
59#include "AnasaziBasicEigenproblem.hpp"
60#include "AnasaziThyraAdapter.hpp"
61#include "AnasaziBlockKrylovSchurSolMgr.hpp"
62#include "AnasaziBlockKrylovSchur.hpp"
63#include "AnasaziStatusTestMaxIters.hpp"
64
65// Isorropia includes
66#ifdef Teko_ENABLE_Isorropia
67#include "Isorropia_EpetraProber.hpp"
68#endif
69
70// Teko includes
71#include "Teko_TpetraHelpers.hpp"
72#include "Teko_TpetraOperatorWrapper.hpp"
73
74// Tpetra
75#include "Thyra_TpetraLinearOp.hpp"
76#include "Tpetra_CrsMatrix.hpp"
77#include "Tpetra_Vector.hpp"
78#include "Thyra_TpetraThyraWrappers.hpp"
79#include "TpetraExt_MatrixMatrix.hpp"
80#include "Tpetra_RowMatrixTransposer.hpp"
81#include "MatrixMarket_Tpetra.hpp"
82
83#include <cmath>
84
85namespace Teko {
86
87using Teuchos::rcp;
88using Teuchos::RCP;
89using Teuchos::rcp_dynamic_cast;
90#ifdef Teko_ENABLE_Isorropia
91using Isorropia::Epetra::Prober;
92#endif
93
94const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream() {
95 Teuchos::RCP<Teuchos::FancyOStream> os = Teuchos::VerboseObjectBase::getDefaultOStream();
96
97 // os->setShowProcRank(true);
98 // os->setOutputToRootOnly(-1);
99 return os;
100}
101
102// distance function...not parallel...entirely internal to this cpp file
103inline double dist(int dim, double *coords, int row, int col) {
104 double value = 0.0;
105 for (int i = 0; i < dim; i++)
106 value += std::pow(coords[dim * row + i] - coords[dim * col + i], 2.0);
107
108 // the distance between the two
109 return std::sqrt(value);
110}
111
112// distance function...not parallel...entirely internal to this cpp file
113inline double dist(double *x, double *y, double *z, int stride, int row, int col) {
114 double value = 0.0;
115 if (x != 0) value += std::pow(x[stride * row] - x[stride * col], 2.0);
116 if (y != 0) value += std::pow(y[stride * row] - y[stride * col], 2.0);
117 if (z != 0) value += std::pow(z[stride * row] - z[stride * col], 2.0);
118
119 // the distance between the two
120 return std::sqrt(value);
121}
122
141#ifdef TEKO_HAVE_EPETRA
142RCP<Epetra_CrsMatrix> buildGraphLaplacian(int dim, double *coords,
143 const Epetra_CrsMatrix &stencil) {
144 // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra
145 // storage
146 RCP<Epetra_CrsMatrix> gl = rcp(
147 new Epetra_CrsMatrix(Copy, stencil.RowMap(), stencil.ColMap(), stencil.MaxNumEntries() + 1),
148 true);
149
150 // allocate an additional value for the diagonal, if neccessary
151 std::vector<double> rowData(stencil.GlobalMaxNumEntries() + 1);
152 std::vector<int> rowInd(stencil.GlobalMaxNumEntries() + 1);
153
154 // loop over all the rows
155 for (int j = 0; j < gl->NumMyRows(); j++) {
156 int row = gl->GRID(j);
157 double diagValue = 0.0;
158 int diagInd = -1;
159 int rowSz = 0;
160
161 // extract a copy of this row...put it in rowData, rowIndicies
162 stencil.ExtractGlobalRowCopy(row, stencil.MaxNumEntries(), rowSz, &rowData[0], &rowInd[0]);
163
164 // loop over elements of row
165 for (int i = 0; i < rowSz; i++) {
166 int col = rowInd[i];
167
168 // is this a 0 entry masquerading as some thing else?
169 double value = rowData[i];
170 if (value == 0) continue;
171
172 // for nondiagonal entries
173 if (row != col) {
174 double d = dist(dim, coords, row, col);
175 rowData[i] = -1.0 / d;
176 diagValue += rowData[i];
177 } else
178 diagInd = i;
179 }
180
181 // handle diagonal entry
182 if (diagInd < 0) { // diagonal not in row
183 rowData[rowSz] = -diagValue;
184 rowInd[rowSz] = row;
185 rowSz++;
186 } else { // diagonal in row
187 rowData[diagInd] = -diagValue;
188 rowInd[diagInd] = row;
189 }
190
191 // insert row data into graph Laplacian matrix
192 TEUCHOS_TEST_FOR_EXCEPT(gl->InsertGlobalValues(row, rowSz, &rowData[0], &rowInd[0]));
193 }
194
195 gl->FillComplete();
196
197 return gl;
198}
199#endif
200
201RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> buildGraphLaplacian(
202 int dim, ST *coords, const Tpetra::CrsMatrix<ST, LO, GO, NT> &stencil) {
203 // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra
204 // storage
205 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> gl = rcp(new Tpetra::CrsMatrix<ST, LO, GO, NT>(
206 stencil.getRowMap(), stencil.getColMap(), stencil.getGlobalMaxNumRowEntries() + 1));
207
208 // allocate an additional value for the diagonal, if neccessary
209 auto rowInd = typename Tpetra::CrsMatrix<ST, LO, GO, NT>::nonconst_global_inds_host_view_type(
210 Kokkos::ViewAllocateWithoutInitializing("rowIndices"),
211 stencil.getGlobalMaxNumRowEntries() + 1);
212 auto rowData = typename Tpetra::CrsMatrix<ST, LO, GO, NT>::nonconst_values_host_view_type(
213 Kokkos::ViewAllocateWithoutInitializing("rowIndices"),
214 stencil.getGlobalMaxNumRowEntries() + 1);
215
216 // loop over all the rows
217 for (LO j = 0; j < (LO)gl->getLocalNumRows(); j++) {
218 GO row = gl->getRowMap()->getGlobalElement(j);
219 ST diagValue = 0.0;
220 GO diagInd = -1;
221 size_t rowSz = 0;
222
223 // extract a copy of this row...put it in rowData, rowIndicies
224 stencil.getGlobalRowCopy(row, rowInd, rowData, rowSz);
225
226 // loop over elements of row
227 for (size_t i = 0; i < rowSz; i++) {
228 GO col = rowInd(i);
229
230 // is this a 0 entry masquerading as some thing else?
231 ST value = rowData[i];
232 if (value == 0) continue;
233
234 // for nondiagonal entries
235 if (row != col) {
236 ST d = dist(dim, coords, row, col);
237 rowData[i] = -1.0 / d;
238 diagValue += rowData(i);
239 } else
240 diagInd = i;
241 }
242
243 // handle diagonal entry
244 if (diagInd < 0) { // diagonal not in row
245 rowData(rowSz) = -diagValue;
246 rowInd(rowSz) = row;
247 rowSz++;
248 } else { // diagonal in row
249 rowData(diagInd) = -diagValue;
250 rowInd(diagInd) = row;
251 }
252
253 // insert row data into graph Laplacian matrix
254 gl->replaceGlobalValues(row, rowInd, rowData);
255 }
256
257 gl->fillComplete();
258
259 return gl;
260}
261
284#ifdef TEKO_HAVE_EPETRA
285RCP<Epetra_CrsMatrix> buildGraphLaplacian(double *x, double *y, double *z, int stride,
286 const Epetra_CrsMatrix &stencil) {
287 // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra
288 // storage
289 RCP<Epetra_CrsMatrix> gl = rcp(
290 new Epetra_CrsMatrix(Copy, stencil.RowMap(), stencil.ColMap(), stencil.MaxNumEntries() + 1),
291 true);
292
293 // allocate an additional value for the diagonal, if neccessary
294 std::vector<double> rowData(stencil.GlobalMaxNumEntries() + 1);
295 std::vector<int> rowInd(stencil.GlobalMaxNumEntries() + 1);
296
297 // loop over all the rows
298 for (int j = 0; j < gl->NumMyRows(); j++) {
299 int row = gl->GRID(j);
300 double diagValue = 0.0;
301 int diagInd = -1;
302 int rowSz = 0;
303
304 // extract a copy of this row...put it in rowData, rowIndicies
305 stencil.ExtractGlobalRowCopy(row, stencil.MaxNumEntries(), rowSz, &rowData[0], &rowInd[0]);
306
307 // loop over elements of row
308 for (int i = 0; i < rowSz; i++) {
309 int col = rowInd[i];
310
311 // is this a 0 entry masquerading as some thing else?
312 double value = rowData[i];
313 if (value == 0) continue;
314
315 // for nondiagonal entries
316 if (row != col) {
317 double d = dist(x, y, z, stride, row, col);
318 rowData[i] = -1.0 / d;
319 diagValue += rowData[i];
320 } else
321 diagInd = i;
322 }
323
324 // handle diagonal entry
325 if (diagInd < 0) { // diagonal not in row
326 rowData[rowSz] = -diagValue;
327 rowInd[rowSz] = row;
328 rowSz++;
329 } else { // diagonal in row
330 rowData[diagInd] = -diagValue;
331 rowInd[diagInd] = row;
332 }
333
334 // insert row data into graph Laplacian matrix
335 TEUCHOS_TEST_FOR_EXCEPT(gl->InsertGlobalValues(row, rowSz, &rowData[0], &rowInd[0]));
336 }
337
338 gl->FillComplete();
339
340 return gl;
341}
342#endif
343
344RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> buildGraphLaplacian(
345 ST *x, ST *y, ST *z, GO stride, const Tpetra::CrsMatrix<ST, LO, GO, NT> &stencil) {
346 // allocate a new matrix with storage for the laplacian...in case of diagonals add one extra
347 // storage
348 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> gl =
349 rcp(new Tpetra::CrsMatrix<ST, LO, GO, NT>(stencil.getRowMap(), stencil.getColMap(),
350 stencil.getGlobalMaxNumRowEntries() + 1),
351 true);
352
353 // allocate an additional value for the diagonal, if neccessary
354 auto rowInd = typename Tpetra::CrsMatrix<ST, LO, GO, NT>::nonconst_global_inds_host_view_type(
355 Kokkos::ViewAllocateWithoutInitializing("rowIndices"),
356 stencil.getGlobalMaxNumRowEntries() + 1);
357 auto rowData = typename Tpetra::CrsMatrix<ST, LO, GO, NT>::nonconst_values_host_view_type(
358 Kokkos::ViewAllocateWithoutInitializing("rowIndices"),
359 stencil.getGlobalMaxNumRowEntries() + 1);
360
361 // loop over all the rows
362 for (LO j = 0; j < (LO)gl->getLocalNumRows(); j++) {
363 GO row = gl->getRowMap()->getGlobalElement(j);
364 ST diagValue = 0.0;
365 GO diagInd = -1;
366 size_t rowSz = 0;
367
368 // extract a copy of this row...put it in rowData, rowIndicies
369 stencil.getGlobalRowCopy(row, rowInd, rowData, rowSz);
370
371 // loop over elements of row
372 for (size_t i = 0; i < rowSz; i++) {
373 GO col = rowInd(i);
374
375 // is this a 0 entry masquerading as some thing else?
376 ST value = rowData[i];
377 if (value == 0) continue;
378
379 // for nondiagonal entries
380 if (row != col) {
381 ST d = dist(x, y, z, stride, row, col);
382 rowData[i] = -1.0 / d;
383 diagValue += rowData(i);
384 } else
385 diagInd = i;
386 }
387
388 // handle diagonal entry
389 if (diagInd < 0) { // diagonal not in row
390 rowData(rowSz) = -diagValue;
391 rowInd(rowSz) = row;
392 rowSz++;
393 } else { // diagonal in row
394 rowData(diagInd) = -diagValue;
395 rowInd(diagInd) = row;
396 }
397
398 // insert row data into graph Laplacian matrix
399 gl->replaceGlobalValues(row, rowInd, rowData);
400 }
401
402 gl->fillComplete();
403
404 return gl;
405}
406
422void applyOp(const LinearOp &A, const MultiVector &x, MultiVector &y, double alpha, double beta) {
423 Thyra::apply(*A, Thyra::NOTRANS, *x, y.ptr(), alpha, beta);
424}
425
441void applyTransposeOp(const LinearOp &A, const MultiVector &x, MultiVector &y, double alpha,
442 double beta) {
443 Thyra::apply(*A, Thyra::TRANS, *x, y.ptr(), alpha, beta);
444}
445
448void update(double alpha, const MultiVector &x, double beta, MultiVector &y) {
449 Teuchos::Array<double> scale;
450 Teuchos::Array<Teuchos::Ptr<const Thyra::MultiVectorBase<double>>> vec;
451
452 // build arrays needed for linear combo
453 scale.push_back(alpha);
454 vec.push_back(x.ptr());
455
456 // compute linear combination
457 Thyra::linear_combination<double>(scale, vec, beta, y.ptr());
458}
459
461BlockedLinearOp getUpperTriBlocks(const BlockedLinearOp &blo, bool callEndBlockFill) {
462 int rows = blockRowCount(blo);
463
464 TEUCHOS_ASSERT(rows == blockColCount(blo));
465
466 RCP<const Thyra::ProductVectorSpaceBase<double>> range = blo->productRange();
467 RCP<const Thyra::ProductVectorSpaceBase<double>> domain = blo->productDomain();
468
469 // allocate new operator
470 BlockedLinearOp upper = createBlockedOp();
471
472 // build new operator
473 upper->beginBlockFill(rows, rows);
474
475 for (int i = 0; i < rows; i++) {
476 // put zero operators on the diagonal
477 // this gurantees the vector space of
478 // the new operator are fully defined
479 RCP<const Thyra::LinearOpBase<double>> zed =
480 Thyra::zero<double>(range->getBlock(i), domain->getBlock(i));
481 upper->setBlock(i, i, zed);
482
483 for (int j = i + 1; j < rows; j++) {
484 // get block i,j
485 LinearOp uij = blo->getBlock(i, j);
486
487 // stuff it in U
488 if (uij != Teuchos::null) upper->setBlock(i, j, uij);
489 }
490 }
491 if (callEndBlockFill) upper->endBlockFill();
492
493 return upper;
494}
495
497BlockedLinearOp getLowerTriBlocks(const BlockedLinearOp &blo, bool callEndBlockFill) {
498 int rows = blockRowCount(blo);
499
500 TEUCHOS_ASSERT(rows == blockColCount(blo));
501
502 RCP<const Thyra::ProductVectorSpaceBase<double>> range = blo->productRange();
503 RCP<const Thyra::ProductVectorSpaceBase<double>> domain = blo->productDomain();
504
505 // allocate new operator
506 BlockedLinearOp lower = createBlockedOp();
507
508 // build new operator
509 lower->beginBlockFill(rows, rows);
510
511 for (int i = 0; i < rows; i++) {
512 // put zero operators on the diagonal
513 // this gurantees the vector space of
514 // the new operator are fully defined
515 RCP<const Thyra::LinearOpBase<double>> zed =
516 Thyra::zero<double>(range->getBlock(i), domain->getBlock(i));
517 lower->setBlock(i, i, zed);
518
519 for (int j = 0; j < i; j++) {
520 // get block i,j
521 LinearOp uij = blo->getBlock(i, j);
522
523 // stuff it in U
524 if (uij != Teuchos::null) lower->setBlock(i, j, uij);
525 }
526 }
527 if (callEndBlockFill) lower->endBlockFill();
528
529 return lower;
530}
531
551BlockedLinearOp zeroBlockedOp(const BlockedLinearOp &blo) {
552 int rows = blockRowCount(blo);
553
554 TEUCHOS_ASSERT(rows == blockColCount(blo)); // assert that matrix is square
555
556 RCP<const Thyra::ProductVectorSpaceBase<double>> range = blo->productRange();
557 RCP<const Thyra::ProductVectorSpaceBase<double>> domain = blo->productDomain();
558
559 // allocate new operator
560 BlockedLinearOp zeroOp = createBlockedOp();
561
562 // build new operator
563 zeroOp->beginBlockFill(rows, rows);
564
565 for (int i = 0; i < rows; i++) {
566 // put zero operators on the diagonal
567 // this gurantees the vector space of
568 // the new operator are fully defined
569 RCP<const Thyra::LinearOpBase<double>> zed =
570 Thyra::zero<double>(range->getBlock(i), domain->getBlock(i));
571 zeroOp->setBlock(i, i, zed);
572 }
573
574 return zeroOp;
575}
576
578bool isZeroOp(const LinearOp op) {
579 // if operator is null...then its zero!
580 if (op == Teuchos::null) return true;
581
582 // try to cast it to a zero linear operator
583 LinearOp test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double>>(op);
584
585 // if it works...then its zero...otherwise its null
586 if (test != Teuchos::null) return true;
587
588 // See if the operator is a wrapped zero op
589 ST scalar = 0.0;
590 Thyra::EOpTransp transp = Thyra::NOTRANS;
591 RCP<const Thyra::LinearOpBase<ST>> wrapped_op;
592 Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
593 test = rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<double>>(wrapped_op);
594 return test != Teuchos::null;
595}
596
597std::pair<ModifiableLinearOp, bool> getAbsRowSumMatrixEpetra(const LinearOp &op) {
598#ifndef TEKO_HAVE_EPETRA
599 return std::make_pair(ModifiableLinearOp{}, false);
600#else
601 RCP<const Epetra_CrsMatrix> eCrsOp;
602
603 const auto eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op);
604
605 if (!eOp) {
606 return std::make_pair(ModifiableLinearOp{}, false);
607 }
608
609 eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
610
611 // extract diagonal
612 const auto ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
613 Epetra_Vector &diag = *ptrDiag;
614
615 // compute absolute value row sum
616 diag.PutScalar(0.0);
617 for (int i = 0; i < eCrsOp->NumMyRows(); i++) {
618 double *values = 0;
619 int numEntries;
620 eCrsOp->ExtractMyRowView(i, numEntries, values);
621
622 // build abs value row sum
623 for (int j = 0; j < numEntries; j++) diag[i] += std::abs(values[j]);
624 }
625
626 // build Thyra diagonal operator
627 return std::make_pair(Teko::Epetra::thyraDiagOp(ptrDiag, eCrsOp->RowMap(),
628 "absRowSum( " + op->getObjectLabel() + " )"),
629 true);
630#endif
631}
632
633std::pair<ModifiableLinearOp, bool> getAbsRowSumMatrixTpetra(const LinearOp &op) {
634 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp;
635
636 const auto tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op);
637
638 tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
639 true);
640
641 // extract diagonal
642 const auto ptrDiag = Tpetra::createVector<ST, LO, GO, NT>(tCrsOp->getRowMap());
643 auto &diag = *ptrDiag;
644
645 // compute absolute value row sum
646 diag.putScalar(0.0);
647 for (LO i = 0; i < (LO)tCrsOp->getLocalNumRows(); i++) {
648 auto numEntries = tCrsOp->getNumEntriesInLocalRow(i);
649 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::local_inds_host_view_type indices;
650 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::values_host_view_type values;
651 tCrsOp->getLocalRowView(i, indices, values);
652
653 // build abs value row sum
654 for (size_t j = 0; j < numEntries; j++) diag.sumIntoLocalValue(i, std::abs(values(j)));
655 }
656
657 // build Thyra diagonal operator
658 return std::make_pair(
659 Teko::TpetraHelpers::thyraDiagOp(ptrDiag, *tCrsOp->getRowMap(),
660 "absRowSum( " + op->getObjectLabel() + " ))"),
661 true);
662}
663
672ModifiableLinearOp getAbsRowSumMatrix(const LinearOp &op) {
673 try {
674 auto eResult = getAbsRowSumMatrixEpetra(op);
675 if (eResult.second) {
676 return eResult.first;
677 }
678
679 auto tResult = getAbsRowSumMatrixTpetra(op);
680 if (tResult.second) {
681 return tResult.first;
682 } else {
683 throw std::logic_error("Neither Epetra nor Tpetra");
684 }
685 } catch (std::exception &e) {
686 auto out = Teuchos::VerboseObjectBase::getDefaultOStream();
687
688 *out << "Teko: getAbsRowSumMatrix requires an Epetra_CrsMatrix or a "
689 "Tpetra::CrsMatrix\n";
690 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator "
691 "from a \""
692 << op->description() << std::endl;
693 *out << " OR\n";
694 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or "
695 "a Tpetra_Operator to a Tpetra::CrsMatrix\n";
696 *out << std::endl;
697 *out << "*** THROWN EXCEPTION ***\n";
698 *out << e.what() << std::endl;
699 *out << "************************\n";
700
701 throw e;
702 }
703}
704
713ModifiableLinearOp getAbsRowSumInvMatrix(const LinearOp &op) {
714 // if this is a blocked operator, extract diagonals block by block
715 // FIXME: this does not add in values from off-diagonal blocks
716 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_op =
717 rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(op);
718 if (blocked_op != Teuchos::null) {
719 int numRows = blocked_op->productRange()->numBlocks();
720 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
721 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_diag =
722 Thyra::defaultBlockedLinearOp<double>();
723 blocked_diag->beginBlockFill(numRows, numRows);
724 for (int r = 0; r < numRows; ++r) {
725 for (int c = 0; c < numRows; ++c) {
726 if (r == c)
727 blocked_diag->setNonconstBlock(r, c, getAbsRowSumInvMatrix(blocked_op->getBlock(r, c)));
728 else
729 blocked_diag->setBlock(r, c,
730 Thyra::zero<double>(blocked_op->getBlock(r, c)->range(),
731 blocked_op->getBlock(r, c)->domain()));
732 }
733 }
734 blocked_diag->endBlockFill();
735 return blocked_diag;
736 }
737
738 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
739 ST scalar = 0.0;
740 bool transp = false;
741 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp =
742 Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
743
744 // extract diagonal
745 const RCP<Tpetra::Vector<ST, LO, GO, NT>> ptrDiag =
746 Tpetra::createVector<ST, LO, GO, NT>(tCrsOp->getRowMap());
747 Tpetra::Vector<ST, LO, GO, NT> &diag = *ptrDiag;
748
749 // compute absolute value row sum
750 diag.putScalar(0.0);
751 for (LO i = 0; i < (LO)tCrsOp->getLocalNumRows(); i++) {
752 LO numEntries = tCrsOp->getNumEntriesInLocalRow(i);
753 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::local_inds_host_view_type indices;
754 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::values_host_view_type values;
755 tCrsOp->getLocalRowView(i, indices, values);
756
757 // build abs value row sum
758 for (LO j = 0; j < numEntries; j++) diag.sumIntoLocalValue(i, std::abs(values(j)));
759 }
760 diag.scale(scalar);
761 diag.reciprocal(diag); // invert entries
762
763 // build Thyra diagonal operator
764 return Teko::TpetraHelpers::thyraDiagOp(ptrDiag, *tCrsOp->getRowMap(),
765 "absRowSum( " + op->getObjectLabel() + " ))");
766
767 } else {
768#ifdef TEKO_HAVE_EPETRA
769 RCP<const Thyra::EpetraLinearOp> eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op, true);
770 RCP<const Epetra_CrsMatrix> eCrsOp =
771 rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
772
773 // extract diagonal
774 const RCP<Epetra_Vector> ptrDiag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
775 Epetra_Vector &diag = *ptrDiag;
776
777 // compute absolute value row sum
778 diag.PutScalar(0.0);
779 for (int i = 0; i < eCrsOp->NumMyRows(); i++) {
780 double *values = 0;
781 int numEntries;
782 eCrsOp->ExtractMyRowView(i, numEntries, values);
783
784 // build abs value row sum
785 for (int j = 0; j < numEntries; j++) diag[i] += std::abs(values[j]);
786 }
787 diag.Reciprocal(diag); // invert entries
788
789 // build Thyra diagonal operator
790 return Teko::Epetra::thyraDiagOp(ptrDiag, eCrsOp->RowMap(),
791 "absRowSum( " + op->getObjectLabel() + " )");
792#else
793 throw std::logic_error(
794 "getAbsRowSumInvMatrix is trying to use Epetra "
795 "code, but TEKO_HAVE_EPETRA is disabled!");
796#endif
797 }
798}
799
807ModifiableLinearOp getLumpedMatrix(const LinearOp &op) {
808 RCP<Thyra::VectorBase<ST>> ones = Thyra::createMember(op->domain());
809 RCP<Thyra::VectorBase<ST>> diag = Thyra::createMember(op->range());
810
811 // set to all ones
812 Thyra::assign(ones.ptr(), 1.0);
813
814 // compute lumped diagonal
815 // Thyra::apply(*op,Thyra::NONCONJ_ELE,*ones,&*diag);
816 Thyra::apply(*op, Thyra::NOTRANS, *ones, diag.ptr());
817
818 return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
819}
820
829ModifiableLinearOp getInvLumpedMatrix(const LinearOp &op) {
830 RCP<Thyra::VectorBase<ST>> ones = Thyra::createMember(op->domain());
831 RCP<Thyra::VectorBase<ST>> diag = Thyra::createMember(op->range());
832
833 // set to all ones
834 Thyra::assign(ones.ptr(), 1.0);
835
836 // compute lumped diagonal
837 Thyra::apply(*op, Thyra::NOTRANS, *ones, diag.ptr());
838 Thyra::reciprocal(*diag, diag.ptr());
839
840 return rcp(new Thyra::DefaultDiagonalLinearOp<ST>(diag));
841}
842
843const std::pair<ModifiableLinearOp, bool> getDiagonalOpEpetra(const LinearOp &op) {
844#ifndef TEKO_HAVE_EPETRA
845 return std::make_pair(ModifiableLinearOp{}, false);
846#else
847 RCP<const Epetra_CrsMatrix> eCrsOp;
848
849 const auto eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op);
850 if (!eOp) {
851 return std::make_pair(ModifiableLinearOp{}, false);
852 }
853
854 eCrsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
855
856 // extract diagonal
857 const auto diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
858 TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
859
860 // build Thyra diagonal operator
861 return std::make_pair(Teko::Epetra::thyraDiagOp(diag, eCrsOp->RowMap(),
862 "inv(diag( " + op->getObjectLabel() + " ))"),
863 true);
864#endif
865}
866
867const std::pair<ModifiableLinearOp, bool> getDiagonalOpTpetra(const LinearOp &op) {
868 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp;
869
870 const auto tOp = rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op);
871 if (!tOp) {
872 return std::make_pair(ModifiableLinearOp{}, false);
873 }
874
875 tCrsOp = rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
876 true);
877
878 // extract diagonal
879 const auto diag = Tpetra::createVector<ST, LO, GO, NT>(tCrsOp->getRowMap());
880 tCrsOp->getLocalDiagCopy(*diag);
881
882 // build Thyra diagonal operator
883 return std::make_pair(
884 Teko::TpetraHelpers::thyraDiagOp(diag, *tCrsOp->getRowMap(),
885 "inv(diag( " + op->getObjectLabel() + " ))"),
886 true);
887}
888
900const ModifiableLinearOp getDiagonalOp(const LinearOp &op) {
901 try {
902 // get Epetra or Tpetra Operator
903 const auto eDiagOp = getDiagonalOpEpetra(op);
904
905 if (eDiagOp.second) {
906 return eDiagOp.first;
907 }
908
909 const auto tDiagOp = getDiagonalOpTpetra(op);
910 if (tDiagOp.second) {
911 return tDiagOp.first;
912 } else
913 throw std::logic_error("Neither Epetra nor Tpetra");
914 } catch (std::exception &e) {
915 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
916
917 *out << "Teko: getDiagonalOp requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
918 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \""
919 << op->description() << std::endl;
920 *out << " OR\n";
921 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a "
922 "Tpetra::CrsMatrix\n";
923 *out << std::endl;
924 *out << "*** THROWN EXCEPTION ***\n";
925 *out << e.what() << std::endl;
926 *out << "************************\n";
927
928 throw e;
929 }
930}
931
932const MultiVector getDiagonal(const LinearOp &op) {
933 try {
934 // get Epetra or Tpetra Operator
935 auto diagOp = getDiagonalOpEpetra(op);
936
937 if (!diagOp.second) {
938 diagOp = getDiagonalOpTpetra(op);
939
940 if (!diagOp.second) {
941 throw std::logic_error("Neither Epetra nor Tpetra");
942 }
943 }
944
945 Teuchos::RCP<const Thyra::MultiVectorBase<double>> v =
946 Teuchos::rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double>>(diagOp.first)
947 ->getDiag();
948 return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double>>(v);
949 } catch (std::exception &e) {
950 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
951
952 *out << "Teko: getDiagonal requires an Epetra_CrsMatrix or a Tpetra::CrsMatrix\n";
953 *out << " Could not extract an Epetra_Operator or a Tpetra_Operator from a \""
954 << op->description() << std::endl;
955 *out << " OR\n";
956 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix or a Tpetra_Operator to a "
957 "Tpetra::CrsMatrix\n";
958 *out << std::endl;
959 *out << "*** THROWN EXCEPTION ***\n";
960 *out << e.what() << std::endl;
961 *out << "************************\n";
962
963 throw e;
964 }
965}
966
967const MultiVector getDiagonal(const Teko::LinearOp &A, const DiagonalType &dt) {
968 LinearOp diagOp = Teko::getDiagonalOp(A, dt);
969
970 Teuchos::RCP<const Thyra::MultiVectorBase<double>> v =
971 Teuchos::rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double>>(diagOp)->getDiag();
972 return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double>>(v);
973}
974
986const ModifiableLinearOp getInvDiagonalOp(const LinearOp &op) {
987 // if this is a diagonal linear op already, just take the reciprocal
988 auto diagonal_op = rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double>>(op);
989 if (diagonal_op != Teuchos::null) {
990 auto diag = diagonal_op->getDiag();
991 auto inv_diag = diag->clone_v();
992 Thyra::reciprocal(*diag, inv_diag.ptr());
993 return rcp(new Thyra::DefaultDiagonalLinearOp<double>(inv_diag));
994 }
995
996 // if this is a blocked operator, extract diagonals block by block
997 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_op =
998 rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(op);
999 if (blocked_op != Teuchos::null) {
1000 int numRows = blocked_op->productRange()->numBlocks();
1001 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == numRows);
1002 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_diag =
1003 Thyra::defaultBlockedLinearOp<double>();
1004 blocked_diag->beginBlockFill(numRows, numRows);
1005 for (int r = 0; r < numRows; ++r) {
1006 for (int c = 0; c < numRows; ++c) {
1007 if (r == c)
1008 blocked_diag->setNonconstBlock(r, c, getInvDiagonalOp(blocked_op->getBlock(r, c)));
1009 else
1010 blocked_diag->setBlock(r, c,
1011 Thyra::zero<double>(blocked_op->getBlock(r, c)->range(),
1012 blocked_op->getBlock(r, c)->domain()));
1013 }
1014 }
1015 blocked_diag->endBlockFill();
1016 return blocked_diag;
1017 }
1018
1019 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
1020 ST scalar = 0.0;
1021 bool transp = false;
1022 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp =
1023 Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
1024
1025 // extract diagonal
1026 const RCP<Tpetra::Vector<ST, LO, GO, NT>> diag =
1027 Tpetra::createVector<ST, LO, GO, NT>(tCrsOp->getRowMap());
1028 diag->scale(scalar);
1029 tCrsOp->getLocalDiagCopy(*diag);
1030 diag->reciprocal(*diag);
1031
1032 // build Thyra diagonal operator
1033 return Teko::TpetraHelpers::thyraDiagOp(diag, *tCrsOp->getRowMap(),
1034 "inv(diag( " + op->getObjectLabel() + " ))");
1035
1036 } else {
1037#ifdef TEKO_HAVE_EPETRA
1038 RCP<const Thyra::EpetraLinearOp> eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op, true);
1039 RCP<const Epetra_CrsMatrix> eCrsOp =
1040 rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
1041
1042 // extract diagonal
1043 const RCP<Epetra_Vector> diag = rcp(new Epetra_Vector(eCrsOp->RowMap()));
1044 TEUCHOS_TEST_FOR_EXCEPT(eCrsOp->ExtractDiagonalCopy(*diag));
1045 diag->Reciprocal(*diag);
1046
1047 // build Thyra diagonal operator
1048 return Teko::Epetra::thyraDiagOp(diag, eCrsOp->RowMap(),
1049 "inv(diag( " + op->getObjectLabel() + " ))");
1050#else
1051 throw std::logic_error(
1052 "getInvDiagonalOp is trying to use Epetra "
1053 "code, but TEKO_HAVE_EPETRA is disabled!");
1054#endif
1055 }
1056}
1057
1070const LinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opm, const LinearOp &opr) {
1071 // if this is a blocked operator, multiply block by block
1072 // it is possible that not every factor in the product is blocked and these situations are handled
1073 // separately
1074
1075 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1076 bool isBlockedM = isPhysicallyBlockedLinearOp(opm);
1077 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1078
1079 // all factors blocked
1080 if ((isBlockedL && isBlockedM && isBlockedR)) {
1081 double scalarl = 0.0;
1082 bool transpl = false;
1083 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1084 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1085 double scalarm = 0.0;
1086 bool transpm = false;
1087 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opm =
1088 getPhysicallyBlockedLinearOp(opm, &scalarm, &transpm);
1089 double scalarr = 0.0;
1090 bool transpr = false;
1091 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1092 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1093 double scalar = scalarl * scalarm * scalarr;
1094
1095 int numRows = blocked_opl->productRange()->numBlocks();
1096 int numCols = blocked_opr->productDomain()->numBlocks();
1097 int numMiddle = blocked_opm->productRange()->numBlocks();
1098
1099 // Assume that the middle block is block nxn and that it's diagonal. Otherwise use the two
1100 // argument explicitMultiply twice
1101 TEUCHOS_ASSERT(blocked_opm->productDomain()->numBlocks() == numMiddle);
1102 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1103 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1104
1105 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1106 Thyra::defaultBlockedLinearOp<double>();
1107 blocked_product->beginBlockFill(numRows, numCols);
1108 for (int r = 0; r < numRows; ++r) {
1109 for (int c = 0; c < numCols; ++c) {
1110 LinearOp product_rc = explicitMultiply(
1111 blocked_opl->getBlock(r, 0), blocked_opm->getBlock(0, 0), blocked_opr->getBlock(0, c));
1112 for (int m = 1; m < numMiddle; ++m) {
1113 LinearOp product_m =
1114 explicitMultiply(blocked_opl->getBlock(r, m), blocked_opm->getBlock(m, m),
1115 blocked_opr->getBlock(m, c));
1116 product_rc = explicitAdd(product_rc, product_m);
1117 }
1118 blocked_product->setBlock(r, c, product_rc);
1119 }
1120 }
1121 blocked_product->endBlockFill();
1122 return Thyra::scale<double>(scalar, blocked_product.getConst());
1123 }
1124
1125 // left and right factors blocked
1126 if (isBlockedL && !isBlockedM && isBlockedR) {
1127 double scalarl = 0.0;
1128 bool transpl = false;
1129 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1130 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1131 double scalarr = 0.0;
1132 bool transpr = false;
1133 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1134 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1135 double scalar = scalarl * scalarr;
1136
1137 int numRows = blocked_opl->productRange()->numBlocks();
1138 int numCols = blocked_opr->productDomain()->numBlocks();
1139 int numMiddle = 1;
1140
1141 // Assume that the middle block is 1x1 diagonal. Left must be rx1, right 1xc
1142 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1143 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1144
1145 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1146 Thyra::defaultBlockedLinearOp<double>();
1147 blocked_product->beginBlockFill(numRows, numCols);
1148 for (int r = 0; r < numRows; ++r) {
1149 for (int c = 0; c < numCols; ++c) {
1150 LinearOp product_rc =
1151 explicitMultiply(blocked_opl->getBlock(r, 0), opm, blocked_opr->getBlock(0, c));
1152 blocked_product->setBlock(r, c, product_rc);
1153 }
1154 }
1155 blocked_product->endBlockFill();
1156 return Thyra::scale<double>(scalar, blocked_product.getConst());
1157 }
1158
1159 // only right factor blocked
1160 if (!isBlockedL && !isBlockedM && isBlockedR) {
1161 double scalarr = 0.0;
1162 bool transpr = false;
1163 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1164 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1165 double scalar = scalarr;
1166
1167 int numRows = 1;
1168 int numCols = blocked_opr->productDomain()->numBlocks();
1169 int numMiddle = 1;
1170
1171 // Assume that the middle block is 1x1 diagonal, left is 1x1. Right must be 1xc
1172 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1173
1174 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1175 Thyra::defaultBlockedLinearOp<double>();
1176 blocked_product->beginBlockFill(numRows, numCols);
1177 for (int c = 0; c < numCols; ++c) {
1178 LinearOp product_c = explicitMultiply(opl, opm, blocked_opr->getBlock(0, c));
1179 blocked_product->setBlock(0, c, product_c);
1180 }
1181 blocked_product->endBlockFill();
1182 return Thyra::scale<double>(scalar, blocked_product.getConst());
1183 }
1184
1185 // TODO: three more cases (only non-blocked - blocked - non-blocked not possible)
1186
1187 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1188 bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1189 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1190
1191 if (isTpetral && isTpetram &&
1192 isTpetrar) { // Both operators are Tpetra matrices so explicitly multiply them
1193
1194 // Get left and right Tpetra crs operators
1195 ST scalarl = 0.0;
1196 bool transpl = false;
1197 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1198 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1199 ST scalarm = 0.0;
1200 bool transpm = false;
1201 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpm =
1202 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarm, &transpm);
1203 ST scalarr = 0.0;
1204 bool transpr = false;
1205 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1206 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1207
1208 // Build output operator
1209 RCP<Thyra::LinearOpBase<ST>> explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1210 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
1211 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
1212
1213 // Do explicit matrix-matrix multiply
1214 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOplm =
1215 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1216 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1217 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1218 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOpl, transpl, *tCrsOpm, transpm, *tCrsOplm);
1219 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOplm, false, *tCrsOpr, transpr,
1220 *explicitCrsOp);
1221 explicitCrsOp->resumeFill();
1222 explicitCrsOp->scale(scalarl * scalarm * scalarr);
1223 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1224 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1225 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1226 explicitCrsOp);
1227 return tExplicitOp;
1228
1229 } else if (isTpetral && !isTpetram && isTpetrar) { // Assume that the middle operator is diagonal
1230
1231 // Get left and right Tpetra crs operators
1232 ST scalarl = 0.0;
1233 bool transpl = false;
1234 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1235 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1236 ST scalarr = 0.0;
1237 bool transpr = false;
1238 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1239 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1240
1241 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr;
1242
1243 // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1244 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpm =
1245 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opm);
1246 if (dOpm != Teuchos::null) {
1247 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1248 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpm->getDiag(), true);
1249 diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(),
1250 true);
1251 }
1252 // If it's not diagonal, maybe it's zero
1253 else if (rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST>>(opm) != Teuchos::null) {
1254 diagPtr = rcp(new Tpetra::Vector<ST, LO, GO, NT>(tCrsOpl->getDomainMap()));
1255 } else
1256 TEUCHOS_ASSERT(false);
1257
1258 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOplm =
1259 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1260 tCrsOpl, Tpetra::Import<LO, GO, NT>(tCrsOpl->getRowMap(), tCrsOpl->getRowMap()));
1261
1262 // Do the diagonal scaling
1263 tCrsOplm->rightScale(*diagPtr);
1264
1265 // Build output operator
1266 RCP<Thyra::LinearOpBase<ST>> explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1267 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
1268 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
1269
1270 // Do explicit matrix-matrix multiply
1271 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1272 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1273 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOplm, false, *tCrsOpr, transpr,
1274 *explicitCrsOp);
1275 explicitCrsOp->resumeFill();
1276 explicitCrsOp->scale(scalarl * scalarr);
1277 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1278 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1279 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1280 explicitCrsOp);
1281 return tExplicitOp;
1282
1283 } else { // Assume Epetra and we can use transformers
1284#ifdef TEKO_HAVE_EPETRA
1285 // build implicit multiply
1286 const LinearOp implicitOp = Thyra::multiply(opl, opm, opr);
1287
1288 // build transformer
1289 const RCP<Thyra::LinearOpTransformerBase<double>> prodTrans =
1290 Thyra::epetraExtDiagScaledMatProdTransformer();
1291
1292 // build operator and multiply
1293 const RCP<Thyra::LinearOpBase<double>> explicitOp = prodTrans->createOutputOp();
1294 prodTrans->transform(*implicitOp, explicitOp.ptr());
1295 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " * " +
1296 opm->getObjectLabel() + " * " + opr->getObjectLabel() + " )");
1297
1298 return explicitOp;
1299#else
1300 throw std::logic_error(
1301 "explicitMultiply is trying to use Epetra "
1302 "code, but TEKO_HAVE_EPETRA is disabled!");
1303#endif
1304 }
1305}
1306
1321const ModifiableLinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opm,
1322 const LinearOp &opr, const ModifiableLinearOp &destOp) {
1323 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1324 bool isTpetram = Teko::TpetraHelpers::isTpetraLinearOp(opm);
1325 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1326
1327 if (isTpetral && isTpetram &&
1328 isTpetrar) { // Both operators are Tpetra matrices so explicitly multiply them
1329
1330 // Get left and right Tpetra crs operators
1331 ST scalarl = 0.0;
1332 bool transpl = false;
1333 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1334 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1335 ST scalarm = 0.0;
1336 bool transpm = false;
1337 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpm =
1338 Teko::TpetraHelpers::getTpetraCrsMatrix(opm, &scalarm, &transpm);
1339 ST scalarr = 0.0;
1340 bool transpr = false;
1341 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1342 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1343
1344 // Build output operator
1345 auto tExplicitOp = rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(destOp);
1346 if (tExplicitOp.is_null()) tExplicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1347
1348 // Do explicit matrix-matrix multiply
1349 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOplm =
1350 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1351 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1352 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1353 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOpl, transpl, *tCrsOpm, transpm, *tCrsOplm);
1354 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOplm, false, *tCrsOpr, transpr,
1355 *explicitCrsOp);
1356 explicitCrsOp->resumeFill();
1357 explicitCrsOp->scale(scalarl * scalarm * scalarr);
1358 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1359 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1360 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1361 explicitCrsOp);
1362 return tExplicitOp;
1363
1364 } else if (isTpetral && !isTpetram && isTpetrar) { // Assume that the middle operator is diagonal
1365
1366 // Get left and right Tpetra crs operators
1367 ST scalarl = 0.0;
1368 bool transpl = false;
1369 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1370 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1371 ST scalarr = 0.0;
1372 bool transpr = false;
1373 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1374 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1375
1376 // Cast middle operator as DiagonalLinearOp and extract diagonal as Vector
1377 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpm =
1378 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opm, true);
1379 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1380 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpm->getDiag(), true);
1381 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr =
1382 rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(), true);
1383 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOplm =
1384 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1385 tCrsOpl, Tpetra::Import<LO, GO, NT>(tCrsOpl->getRowMap(), tCrsOpl->getRowMap()));
1386
1387 // Do the diagonal scaling
1388 tCrsOplm->rightScale(*diagPtr);
1389
1390 // Build output operator
1391 RCP<Thyra::LinearOpBase<ST>> explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1392 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
1393 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
1394
1395 // Do explicit matrix-matrix multiply
1396 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1397 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1398 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOplm, false, *tCrsOpr, transpr,
1399 *explicitCrsOp);
1400 explicitCrsOp->resumeFill();
1401 explicitCrsOp->scale(scalarl * scalarr);
1402 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1403 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1404 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1405 explicitCrsOp);
1406 return tExplicitOp;
1407
1408 } else { // Assume Epetra and we can use transformers
1409#ifdef TEKO_HAVE_EPETRA
1410 // build implicit multiply
1411 const LinearOp implicitOp = Thyra::multiply(opl, opm, opr);
1412
1413 // build transformer
1414 const RCP<Thyra::LinearOpTransformerBase<double>> prodTrans =
1415 Thyra::epetraExtDiagScaledMatProdTransformer();
1416
1417 // build operator destination operator
1418 ModifiableLinearOp explicitOp;
1419
1420 // if neccessary build a operator to put the explicit multiply into
1421 if (destOp == Teuchos::null)
1422 explicitOp = prodTrans->createOutputOp();
1423 else
1424 explicitOp = destOp;
1425
1426 // perform multiplication
1427 prodTrans->transform(*implicitOp, explicitOp.ptr());
1428
1429 // label it
1430 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " * " +
1431 opm->getObjectLabel() + " * " + opr->getObjectLabel() + " )");
1432
1433 return explicitOp;
1434#else
1435 throw std::logic_error(
1436 "explicitMultiply is trying to use Epetra "
1437 "code, but TEKO_HAVE_EPETRA is disabled!");
1438#endif
1439 }
1440}
1441
1452const LinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opr) {
1453 // if this is a blocked operator, multiply block by block
1454 // it is possible that not every factor in the product is blocked and these situations are handled
1455 // separately
1456
1457 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1458 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1459
1460 // both factors blocked
1461 if ((isBlockedL && isBlockedR)) {
1462 double scalarl = 0.0;
1463 bool transpl = false;
1464 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1465 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1466 double scalarr = 0.0;
1467 bool transpr = false;
1468 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1469 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1470 double scalar = scalarl * scalarr;
1471
1472 int numRows = blocked_opl->productRange()->numBlocks();
1473 int numCols = blocked_opr->productDomain()->numBlocks();
1474 int numMiddle = blocked_opl->productDomain()->numBlocks();
1475
1476 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1477
1478 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1479 Thyra::defaultBlockedLinearOp<double>();
1480 blocked_product->beginBlockFill(numRows, numCols);
1481 for (int r = 0; r < numRows; ++r) {
1482 for (int c = 0; c < numCols; ++c) {
1483 LinearOp product_rc =
1484 explicitMultiply(blocked_opl->getBlock(r, 0), blocked_opr->getBlock(0, c));
1485 for (int m = 1; m < numMiddle; ++m) {
1486 LinearOp product_m =
1487 explicitMultiply(blocked_opl->getBlock(r, m), blocked_opr->getBlock(m, c));
1488 product_rc = explicitAdd(product_rc, product_m);
1489 }
1490 blocked_product->setBlock(r, c, Thyra::scale(scalar, product_rc));
1491 }
1492 }
1493 blocked_product->endBlockFill();
1494 return blocked_product;
1495 }
1496
1497 // only left factor blocked
1498 if ((isBlockedL && !isBlockedR)) {
1499 double scalarl = 0.0;
1500 bool transpl = false;
1501 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1502 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1503 double scalar = scalarl;
1504
1505 int numRows = blocked_opl->productRange()->numBlocks();
1506 int numCols = 1;
1507 int numMiddle = 1;
1508
1509 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1510
1511 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1512 Thyra::defaultBlockedLinearOp<double>();
1513 blocked_product->beginBlockFill(numRows, numCols);
1514 for (int r = 0; r < numRows; ++r) {
1515 LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r, 0), opr);
1516 blocked_product->setBlock(r, 0, Thyra::scale(scalar, product_r));
1517 }
1518 blocked_product->endBlockFill();
1519 return blocked_product;
1520 }
1521
1522 // only right factor blocked
1523 if ((!isBlockedL && isBlockedR)) {
1524 double scalarr = 0.0;
1525 bool transpr = false;
1526 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1527 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1528 double scalar = scalarr;
1529
1530 int numRows = 1;
1531 int numCols = blocked_opr->productDomain()->numBlocks();
1532 int numMiddle = 1;
1533
1534 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1535
1536 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1537 Thyra::defaultBlockedLinearOp<double>();
1538 blocked_product->beginBlockFill(numRows, numCols);
1539 for (int c = 0; c < numCols; ++c) {
1540 LinearOp product_c = explicitMultiply(opl, blocked_opr->getBlock(0, c));
1541 blocked_product->setBlock(0, c, Thyra::scale(scalar, product_c));
1542 }
1543 blocked_product->endBlockFill();
1544 return blocked_product;
1545 }
1546
1547 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1548 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1549
1550 if (isTpetral && isTpetrar) { // Both operators are Tpetra matrices so explicitly multiply them
1551 // Get left and right Tpetra crs operators
1552 ST scalarl = 0.0;
1553 bool transpl = false;
1554 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1555 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1556 ST scalarr = 0.0;
1557 bool transpr = false;
1558 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1559 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1560
1561 // Build output operator
1562 RCP<Thyra::LinearOpBase<ST>> explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1563 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
1564 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
1565
1566 // Do explicit matrix-matrix multiply
1567 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1568 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1569 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOpl, transpl, *tCrsOpr, transpr,
1570 *explicitCrsOp);
1571 explicitCrsOp->resumeFill();
1572 explicitCrsOp->scale(scalarl * scalarr);
1573 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1574 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1575 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1576 explicitCrsOp);
1577 return tExplicitOp;
1578
1579 } else if (isTpetral && !isTpetrar) { // Assume that the right operator is diagonal
1580
1581 // Get left Tpetra crs operator
1582 ST scalarl = 0.0;
1583 bool transpl = false;
1584 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1585 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1586
1587 // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1588 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpr =
1589 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opr, true);
1590 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1591 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpr->getDiag(), true);
1592 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr =
1593 rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(), true);
1594 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1595 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1596 tCrsOpl, Tpetra::Import<LO, GO, NT>(tCrsOpl->getRowMap(), tCrsOpl->getRowMap()));
1597
1598 explicitCrsOp->rightScale(*diagPtr);
1599 explicitCrsOp->resumeFill();
1600 explicitCrsOp->scale(scalarl);
1601 explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(), tCrsOpl->getRangeMap());
1602
1603 return Thyra::constTpetraLinearOp<ST, LO, GO, NT>(
1604 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1605 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()), explicitCrsOp);
1606
1607 } else if (!isTpetral && isTpetrar) { // Assume that the left operator is diagonal
1608
1609 // Get right Tpetra crs operator
1610 ST scalarr = 0.0;
1611 bool transpr = false;
1612 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1613 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1614
1615 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr;
1616
1617 // Cast left operator as DiagonalLinearOp and extract diagonal as Vector
1618 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpl =
1619 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opl);
1620 if (dOpl != Teuchos::null) {
1621 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1622 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpl->getDiag(), true);
1623 diagPtr = rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(),
1624 true);
1625 }
1626 // If it's not diagonal, maybe it's zero
1627 else if (rcp_dynamic_cast<const Thyra::ZeroLinearOpBase<ST>>(opl) != Teuchos::null) {
1628 diagPtr = rcp(new Tpetra::Vector<ST, LO, GO, NT>(tCrsOpr->getRangeMap()));
1629 } else
1630 TEUCHOS_ASSERT(false);
1631
1632 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1633 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1634 tCrsOpr, Tpetra::Import<LO, GO, NT>(tCrsOpr->getRowMap(), tCrsOpr->getRowMap()));
1635
1636 explicitCrsOp->leftScale(*diagPtr);
1637 explicitCrsOp->resumeFill();
1638 explicitCrsOp->scale(scalarr);
1639 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpr->getRangeMap());
1640
1641 return Thyra::constTpetraLinearOp<ST, LO, GO, NT>(
1642 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1643 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()), explicitCrsOp);
1644
1645 } else { // Assume Epetra and we can use transformers
1646#ifdef TEKO_HAVE_EPETRA
1647 // build implicit multiply
1648 const LinearOp implicitOp = Thyra::multiply(opl, opr);
1649
1650 // build a scaling transformer
1651 RCP<Thyra::LinearOpTransformerBase<double>> prodTrans =
1652 Thyra::epetraExtDiagScalingTransformer();
1653
1654 // check to see if a scaling transformer works: if not use the
1655 // DiagScaledMatrixProduct transformer
1656 if (not prodTrans->isCompatible(*implicitOp))
1657 prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1658
1659 // build operator and multiply
1660 const RCP<Thyra::LinearOpBase<double>> explicitOp = prodTrans->createOutputOp();
1661 prodTrans->transform(*implicitOp, explicitOp.ptr());
1662 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " * " +
1663 opr->getObjectLabel() + " )");
1664
1665 return explicitOp;
1666#else
1667 throw std::logic_error(
1668 "explicitMultiply is trying to use Epetra "
1669 "code, but TEKO_HAVE_EPETRA is disabled!");
1670#endif
1671 }
1672}
1673
1687const ModifiableLinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opr,
1688 const ModifiableLinearOp &destOp) {
1689 // if this is a blocked operator, multiply block by block
1690 // it is possible that not every factor in the product is blocked and these situations are handled
1691 // separately
1692
1693 bool isBlockedL = isPhysicallyBlockedLinearOp(opl);
1694 bool isBlockedR = isPhysicallyBlockedLinearOp(opr);
1695
1696 // both factors blocked
1697 if ((isBlockedL && isBlockedR)) {
1698 double scalarl = 0.0;
1699 bool transpl = false;
1700 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1701 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1702 double scalarr = 0.0;
1703 bool transpr = false;
1704 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1705 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1706 double scalar = scalarl * scalarr;
1707
1708 int numRows = blocked_opl->productRange()->numBlocks();
1709 int numCols = blocked_opr->productDomain()->numBlocks();
1710 int numMiddle = blocked_opl->productDomain()->numBlocks();
1711
1712 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1713
1714 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1715 Thyra::defaultBlockedLinearOp<double>();
1716 blocked_product->beginBlockFill(numRows, numCols);
1717 for (int r = 0; r < numRows; ++r) {
1718 for (int c = 0; c < numCols; ++c) {
1719 LinearOp product_rc =
1720 explicitMultiply(blocked_opl->getBlock(r, 0), blocked_opr->getBlock(0, c));
1721 for (int m = 1; m < numMiddle; ++m) {
1722 LinearOp product_m =
1723 explicitMultiply(blocked_opl->getBlock(r, m), blocked_opr->getBlock(m, c));
1724 product_rc = explicitAdd(product_rc, product_m);
1725 }
1726 blocked_product->setBlock(r, c, Thyra::scale(scalar, product_rc));
1727 }
1728 }
1729 blocked_product->endBlockFill();
1730 return blocked_product;
1731 }
1732
1733 // only left factor blocked
1734 if ((isBlockedL && !isBlockedR)) {
1735 double scalarl = 0.0;
1736 bool transpl = false;
1737 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1738 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
1739 double scalar = scalarl;
1740
1741 int numRows = blocked_opl->productRange()->numBlocks();
1742 int numCols = 1;
1743 int numMiddle = 1;
1744
1745 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == numMiddle);
1746
1747 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1748 Thyra::defaultBlockedLinearOp<double>();
1749 blocked_product->beginBlockFill(numRows, numCols);
1750 for (int r = 0; r < numRows; ++r) {
1751 LinearOp product_r = explicitMultiply(blocked_opl->getBlock(r, 0), opr);
1752 blocked_product->setBlock(r, 0, Thyra::scale(scalar, product_r));
1753 }
1754 blocked_product->endBlockFill();
1755 return blocked_product;
1756 }
1757
1758 // only right factor blocked
1759 if ((!isBlockedL && isBlockedR)) {
1760 double scalarr = 0.0;
1761 bool transpr = false;
1762 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1763 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
1764 double scalar = scalarr;
1765
1766 int numRows = 1;
1767 int numCols = blocked_opr->productDomain()->numBlocks();
1768 int numMiddle = 1;
1769
1770 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numMiddle);
1771
1772 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_product =
1773 Thyra::defaultBlockedLinearOp<double>();
1774 blocked_product->beginBlockFill(numRows, numCols);
1775 for (int c = 0; c < numCols; ++c) {
1776 LinearOp product_c = explicitMultiply(opl, blocked_opr->getBlock(0, c));
1777 blocked_product->setBlock(0, c, Thyra::scale(scalar, product_c));
1778 }
1779 blocked_product->endBlockFill();
1780 return blocked_product;
1781 }
1782
1783 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1784 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1785
1786 if (isTpetral && isTpetrar) { // Both operators are Tpetra matrices so use the explicit Tpetra
1787 // matrix-matrix multiply
1788
1789 // Get left and right Tpetra crs operators
1790 ST scalarl = 0.0;
1791 bool transpl = false;
1792 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1793 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1794 ST scalarr = 0.0;
1795 bool transpr = false;
1796 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1797 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1798
1799 // Build output operator
1800 RCP<Thyra::LinearOpBase<ST>> explicitOp;
1801 if (destOp != Teuchos::null)
1802 explicitOp = destOp;
1803 else
1804 explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
1805 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
1806 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
1807
1808 // Do explicit matrix-matrix multiply
1809 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1810 Tpetra::createCrsMatrix<ST, LO, GO, NT>(tCrsOpl->getRowMap());
1811 Tpetra::MatrixMatrix::Multiply<ST, LO, GO, NT>(*tCrsOpl, transpl, *tCrsOpr, transpr,
1812 *explicitCrsOp);
1813 explicitCrsOp->resumeFill();
1814 explicitCrsOp->scale(scalarl * scalarr);
1815 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpl->getRangeMap());
1816 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1817 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
1818 explicitCrsOp);
1819 return tExplicitOp;
1820
1821 } else if (isTpetral && !isTpetrar) { // Assume that the right operator is diagonal
1822
1823 // Get left Tpetra crs operator
1824 ST scalarl = 0.0;
1825 bool transpl = false;
1826 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
1827 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
1828
1829 // Cast right operator as DiagonalLinearOp and extract diagonal as Vector
1830 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpr =
1831 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opr);
1832 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1833 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpr->getDiag(), true);
1834 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr =
1835 rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(), true);
1836
1837 // Scale by the diagonal operator
1838 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1839 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1840 tCrsOpl, Tpetra::Import<LO, GO, NT>(tCrsOpl->getRowMap(), tCrsOpl->getRowMap()));
1841 explicitCrsOp->rightScale(*diagPtr);
1842 explicitCrsOp->resumeFill();
1843 explicitCrsOp->scale(scalarl);
1844 explicitCrsOp->fillComplete(tCrsOpl->getDomainMap(), tCrsOpl->getRangeMap());
1845 return Thyra::tpetraLinearOp<ST, LO, GO, NT>(
1846 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1847 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()), explicitCrsOp);
1848
1849 } else if (!isTpetral && isTpetrar) { // Assume that the left operator is diagonal
1850
1851 // Get right Tpetra crs operator
1852 ST scalarr = 0.0;
1853 bool transpr = false;
1854 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
1855 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
1856
1857 // Cast leftt operator as DiagonalLinearOp and extract diagonal as Vector
1858 RCP<const Thyra::DiagonalLinearOpBase<ST>> dOpl =
1859 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<ST>>(opl, true);
1860 RCP<const Thyra::TpetraVector<ST, LO, GO, NT>> tPtr =
1861 rcp_dynamic_cast<const Thyra::TpetraVector<ST, LO, GO, NT>>(dOpl->getDiag(), true);
1862 RCP<const Tpetra::Vector<ST, LO, GO, NT>> diagPtr =
1863 rcp_dynamic_cast<const Tpetra::Vector<ST, LO, GO, NT>>(tPtr->getConstTpetraVector(), true);
1864
1865 // Scale by the diagonal operator
1866 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
1867 Tpetra::importAndFillCompleteCrsMatrix<Tpetra::CrsMatrix<ST, LO, GO, NT>>(
1868 tCrsOpr, Tpetra::Import<LO, GO, NT>(tCrsOpr->getRowMap(), tCrsOpr->getRowMap()));
1869 explicitCrsOp->leftScale(*diagPtr);
1870 explicitCrsOp->resumeFill();
1871 explicitCrsOp->scale(scalarr);
1872 explicitCrsOp->fillComplete(tCrsOpr->getDomainMap(), tCrsOpr->getRangeMap());
1873 return Thyra::tpetraLinearOp<ST, LO, GO, NT>(
1874 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
1875 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()), explicitCrsOp);
1876
1877 } else { // Assume Epetra and we can use transformers
1878#ifdef TEKO_HAVE_EPETRA
1879 // build implicit multiply
1880 const LinearOp implicitOp = Thyra::multiply(opl, opr);
1881
1882 // build a scaling transformer
1883
1884 RCP<Thyra::LinearOpTransformerBase<double>> prodTrans =
1885 Thyra::epetraExtDiagScalingTransformer();
1886
1887 // check to see if a scaling transformer works: if not use the
1888 // DiagScaledMatrixProduct transformer
1889 if (not prodTrans->isCompatible(*implicitOp))
1890 prodTrans = Thyra::epetraExtDiagScaledMatProdTransformer();
1891
1892 // build operator destination operator
1893 ModifiableLinearOp explicitOp;
1894
1895 // if neccessary build a operator to put the explicit multiply into
1896 if (destOp == Teuchos::null)
1897 explicitOp = prodTrans->createOutputOp();
1898 else
1899 explicitOp = destOp;
1900
1901 // perform multiplication
1902 prodTrans->transform(*implicitOp, explicitOp.ptr());
1903
1904 // label it
1905 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " * " +
1906 opr->getObjectLabel() + " )");
1907
1908 return explicitOp;
1909#else
1910 throw std::logic_error(
1911 "explicitMultiply is trying to use Epetra "
1912 "code, but TEKO_HAVE_EPETRA is disabled!");
1913#endif
1914 }
1915}
1916
1927const LinearOp explicitAdd(const LinearOp &opl_in, const LinearOp &opr_in) {
1928 // if both blocked, add block by block
1929 if (isPhysicallyBlockedLinearOp(opl_in) && isPhysicallyBlockedLinearOp(opr_in)) {
1930 double scalarl = 0.0;
1931 bool transpl = false;
1932 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1933 getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1934
1935 double scalarr = 0.0;
1936 bool transpr = false;
1937 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1938 getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1939
1940 int numRows = blocked_opl->productRange()->numBlocks();
1941 int numCols = blocked_opl->productDomain()->numBlocks();
1942 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
1943 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
1944
1945 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_sum =
1946 Thyra::defaultBlockedLinearOp<double>();
1947 blocked_sum->beginBlockFill(numRows, numCols);
1948 for (int r = 0; r < numRows; ++r)
1949 for (int c = 0; c < numCols; ++c)
1950 blocked_sum->setBlock(r, c,
1951 explicitAdd(Thyra::scale(scalarl, blocked_opl->getBlock(r, c)),
1952 Thyra::scale(scalarr, blocked_opr->getBlock(r, c))));
1953 blocked_sum->endBlockFill();
1954 return blocked_sum;
1955 }
1956
1957 // if only one is blocked, it must be 1x1
1958 LinearOp opl = opl_in;
1959 LinearOp opr = opr_in;
1960 if (isPhysicallyBlockedLinearOp(opl_in)) {
1961 double scalarl = 0.0;
1962 bool transpl = false;
1963 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
1964 getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
1965 TEUCHOS_ASSERT(blocked_opl->productRange()->numBlocks() == 1);
1966 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == 1);
1967 opl = Thyra::scale(scalarl, blocked_opl->getBlock(0, 0));
1968 }
1969 if (isPhysicallyBlockedLinearOp(opr_in)) {
1970 double scalarr = 0.0;
1971 bool transpr = false;
1972 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
1973 getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
1974 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == 1);
1975 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == 1);
1976 opr = Thyra::scale(scalarr, blocked_opr->getBlock(0, 0));
1977 }
1978
1979 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
1980 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
1981
1982 // if one of the operators in the sum is a thyra zero op
1983 if (isZeroOp(opl)) {
1984 if (isZeroOp(opr)) return opr; // return a zero op if both are zero
1985 if (isTpetrar) { // if other op is tpetra, replace this with a zero crs matrix
1986 ST scalar = 0.0;
1987 bool transp = false;
1988 auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
1989 auto zero_crs = Tpetra::createCrsMatrix<ST, LO, GO, NT>(crs_op->getRowMap());
1990 zero_crs->fillComplete();
1991 opl = Thyra::constTpetraLinearOp<ST, LO, GO, NT>(
1992 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crs_op->getRangeMap()),
1993 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crs_op->getDomainMap()), zero_crs);
1994 isTpetral = true;
1995 } else
1996 return opr->clone();
1997 }
1998 if (isZeroOp(opr)) {
1999 if (isTpetral) { // if other op is tpetra, replace this with a zero crs matrix
2000 ST scalar = 0.0;
2001 bool transp = false;
2002 auto crs_op = Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalar, &transp);
2003 auto zero_crs = Tpetra::createCrsMatrix<ST, LO, GO, NT>(crs_op->getRowMap());
2004 zero_crs->fillComplete();
2005 opr = Thyra::constTpetraLinearOp<ST, LO, GO, NT>(
2006 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crs_op->getRangeMap()),
2007 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crs_op->getDomainMap()), zero_crs);
2008 isTpetrar = true;
2009 } else
2010 return opl->clone();
2011 }
2012
2013 if (isTpetral && isTpetrar) { // Both operators are Tpetra matrices so use the explicit Tpetra
2014 // matrix-matrix add
2015
2016 // Get left and right Tpetra crs operators
2017 ST scalarl = 0.0;
2018 bool transpl = false;
2019 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
2020 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
2021 ST scalarr = 0.0;
2022 bool transpr = false;
2023 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
2024 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
2025
2026 // Build output operator
2027 RCP<Thyra::LinearOpBase<ST>> explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
2028 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
2029 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
2030
2031 // Do explicit matrix-matrix add
2032 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp =
2033 Tpetra::MatrixMatrix::add<ST, LO, GO, NT>(scalarl, transpl, *tCrsOpl, scalarr, transpr,
2034 *tCrsOpr);
2035 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
2036 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
2037 explicitCrsOp);
2038 return tExplicitOp;
2039
2040 } else { // Assume Epetra
2041#ifdef TEKO_HAVE_EPETRA
2042 // build implicit add
2043 const LinearOp implicitOp = Thyra::add(opl, opr);
2044
2045 // build transformer
2046 const RCP<Thyra::LinearOpTransformerBase<double>> prodTrans = Thyra::epetraExtAddTransformer();
2047
2048 // build operator and add
2049 const RCP<Thyra::LinearOpBase<double>> explicitOp = prodTrans->createOutputOp();
2050 prodTrans->transform(*implicitOp, explicitOp.ptr());
2051 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " + " +
2052 opr->getObjectLabel() + " )");
2053
2054 return explicitOp;
2055#else
2056 throw std::logic_error(
2057 "explicitAdd is trying to use Epetra "
2058 "code, but TEKO_HAVE_EPETRA is disabled!");
2059#endif
2060 }
2061}
2062
2075const ModifiableLinearOp explicitAdd(const LinearOp &opl_in, const LinearOp &opr_in,
2076 const ModifiableLinearOp &destOp) {
2077 // if blocked, add block by block
2078 if (isPhysicallyBlockedLinearOp(opl_in) && isPhysicallyBlockedLinearOp(opr_in)) {
2079 double scalarl = 0.0;
2080 bool transpl = false;
2081 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
2082 getPhysicallyBlockedLinearOp(opl_in, &scalarl, &transpl);
2083
2084 double scalarr = 0.0;
2085 bool transpr = false;
2086 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
2087 getPhysicallyBlockedLinearOp(opr_in, &scalarr, &transpr);
2088
2089 int numRows = blocked_opl->productRange()->numBlocks();
2090 int numCols = blocked_opl->productDomain()->numBlocks();
2091 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == numRows);
2092 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == numCols);
2093
2094 RCP<Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_sum =
2095 Teuchos::rcp_dynamic_cast<Thyra::DefaultBlockedLinearOp<double>>(destOp);
2096 if (blocked_sum.is_null()) {
2097 // take care of the null case, this means we must alllocate memory
2098 blocked_sum = Thyra::defaultBlockedLinearOp<double>();
2099
2100 blocked_sum->beginBlockFill(numRows, numCols);
2101 for (int r = 0; r < numRows; ++r) {
2102 for (int c = 0; c < numCols; ++c) {
2103 auto block =
2104 explicitAdd(Thyra::scale(scalarl, blocked_opl->getBlock(r, c)),
2105 Thyra::scale(scalarr, blocked_opr->getBlock(r, c)), Teuchos::null);
2106 blocked_sum->setNonconstBlock(r, c, block);
2107 }
2108 }
2109 blocked_sum->endBlockFill();
2110
2111 } else {
2112 // in this case memory can be reused
2113 for (int r = 0; r < numRows; ++r)
2114 for (int c = 0; c < numCols; ++c)
2115 explicitAdd(Thyra::scale(scalarl, blocked_opl->getBlock(r, c)),
2116 Thyra::scale(scalarr, blocked_opr->getBlock(r, c)),
2117 blocked_sum->getNonconstBlock(r, c));
2118 }
2119
2120 return blocked_sum;
2121 }
2122
2123 LinearOp opl = opl_in;
2124 LinearOp opr = opr_in;
2125 // if only one is blocked, it must be 1x1
2126 if (isPhysicallyBlockedLinearOp(opl)) {
2127 double scalarl = 0.0;
2128 bool transpl = false;
2129 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opl =
2130 getPhysicallyBlockedLinearOp(opl, &scalarl, &transpl);
2131 TEUCHOS_ASSERT(blocked_opl->productRange()->numBlocks() == 1);
2132 TEUCHOS_ASSERT(blocked_opl->productDomain()->numBlocks() == 1);
2133 opl = Thyra::scale(scalarl, blocked_opl->getBlock(0, 0));
2134 }
2135 if (isPhysicallyBlockedLinearOp(opr)) {
2136 double scalarr = 0.0;
2137 bool transpr = false;
2138 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_opr =
2139 getPhysicallyBlockedLinearOp(opr, &scalarr, &transpr);
2140 TEUCHOS_ASSERT(blocked_opr->productRange()->numBlocks() == 1);
2141 TEUCHOS_ASSERT(blocked_opr->productDomain()->numBlocks() == 1);
2142 opr = Thyra::scale(scalarr, blocked_opr->getBlock(0, 0));
2143 }
2144
2145 bool isTpetral = Teko::TpetraHelpers::isTpetraLinearOp(opl);
2146 bool isTpetrar = Teko::TpetraHelpers::isTpetraLinearOp(opr);
2147
2148 if (isTpetral && isTpetrar) { // Both operators are Tpetra matrices so use the explicit
2149 // Tpetra matrix-matrix add
2150
2151 // Get left and right Tpetra crs operators
2152 ST scalarl = 0.0;
2153 bool transpl = false;
2154 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpl =
2155 Teko::TpetraHelpers::getTpetraCrsMatrix(opl, &scalarl, &transpl);
2156 ST scalarr = 0.0;
2157 bool transpr = false;
2158 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOpr =
2159 Teko::TpetraHelpers::getTpetraCrsMatrix(opr, &scalarr, &transpr);
2160
2161 // Build output operator
2162 RCP<Thyra::LinearOpBase<ST>> explicitOp;
2163 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> explicitCrsOp;
2164 if (!destOp.is_null()) {
2165 explicitOp = destOp;
2166 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tOp =
2167 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(destOp);
2168 if (!tOp.is_null())
2169 explicitCrsOp =
2170 rcp_dynamic_cast<Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getTpetraOperator());
2171 bool needNewTpetraMatrix =
2172 (explicitCrsOp.is_null()) || (tCrsOpl == explicitCrsOp) || (tCrsOpr == explicitCrsOp);
2173 if (!needNewTpetraMatrix) {
2174 try {
2175 // try to reuse matrix sparsity with Add. If it fails, build new operator with add
2176 Tpetra::MatrixMatrix::Add<ST, LO, GO, NT>(*tCrsOpl, transpl, scalarl, *tCrsOpr, transpr,
2177 scalarr, explicitCrsOp);
2178 } catch (std::logic_error &e) {
2179 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2180 *out << "*** THROWN EXCEPTION ***\n";
2181 *out << e.what() << std::endl;
2182 *out << "************************\n";
2183 *out << "Teko: explicitAdd unable to reuse existing operator. Creating new operator.\n"
2184 << std::endl;
2185 needNewTpetraMatrix = true;
2186 }
2187 }
2188 if (needNewTpetraMatrix)
2189 // Do explicit matrix-matrix add
2190 explicitCrsOp = Tpetra::MatrixMatrix::add<ST, LO, GO, NT>(scalarl, transpl, *tCrsOpl,
2191 scalarr, transpr, *tCrsOpr);
2192 } else {
2193 explicitOp = rcp(new Thyra::TpetraLinearOp<ST, LO, GO, NT>());
2194 // Do explicit matrix-matrix add
2195 explicitCrsOp = Tpetra::MatrixMatrix::add<ST, LO, GO, NT>(scalarl, transpl, *tCrsOpl, scalarr,
2196 transpr, *tCrsOpr);
2197 }
2198 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT>> tExplicitOp =
2199 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT>>(explicitOp);
2200
2201 tExplicitOp->initialize(Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getRangeMap()),
2202 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(explicitCrsOp->getDomainMap()),
2203 explicitCrsOp);
2204 return tExplicitOp;
2205
2206 } else { // Assume Epetra
2207#ifdef TEKO_HAVE_EPETRA
2208 // build implicit add
2209 const LinearOp implicitOp = Thyra::add(opl, opr);
2210
2211 // build transformer
2212 const RCP<Thyra::LinearOpTransformerBase<double>> prodTrans = Thyra::epetraExtAddTransformer();
2213
2214 // build or reuse destination operator
2215 RCP<Thyra::LinearOpBase<double>> explicitOp;
2216 if (destOp != Teuchos::null)
2217 explicitOp = destOp;
2218 else
2219 explicitOp = prodTrans->createOutputOp();
2220
2221 // perform add
2222 prodTrans->transform(*implicitOp, explicitOp.ptr());
2223 explicitOp->setObjectLabel("explicit( " + opl->getObjectLabel() + " + " +
2224 opr->getObjectLabel() + " )");
2225
2226 return explicitOp;
2227#else
2228 throw std::logic_error(
2229 "explicitAdd is trying to use Epetra "
2230 "code, but TEKO_HAVE_EPETRA is disabled!");
2231#endif
2232 }
2233}
2234
2239const ModifiableLinearOp explicitSum(const LinearOp &op, const ModifiableLinearOp &destOp) {
2240#ifdef TEKO_HAVE_EPETRA
2241 // convert operators to Epetra_CrsMatrix
2242 const RCP<const Epetra_CrsMatrix> epetraOp =
2243 rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*op), true);
2244
2245 if (destOp == Teuchos::null) {
2246 Teuchos::RCP<Epetra_Operator> epetraDest = Teuchos::rcp(new Epetra_CrsMatrix(*epetraOp));
2247
2248 return Thyra::nonconstEpetraLinearOp(epetraDest);
2249 }
2250
2251 const RCP<Epetra_CrsMatrix> epetraDest =
2252 rcp_dynamic_cast<Epetra_CrsMatrix>(get_Epetra_Operator(*destOp), true);
2253
2254 EpetraExt::MatrixMatrix::Add(*epetraOp, false, 1.0, *epetraDest, 1.0);
2255
2256 return destOp;
2257#else
2258 throw std::logic_error(
2259 "explicitSum is trying to use Epetra "
2260 "code, but TEKO_HAVE_EPETRA is disabled!");
2261#endif
2262}
2263
2264const LinearOp explicitTranspose(const LinearOp &op) {
2265 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2266 RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT>> tOp =
2267 rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op, true);
2268 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp =
2269 rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
2270 true);
2271
2272 Tpetra::RowMatrixTransposer<ST, LO, GO, NT> transposer(tCrsOp);
2273 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT>> transOp = transposer.createTranspose();
2274
2275 return Thyra::tpetraLinearOp<ST, LO, GO, NT>(
2276 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(transOp->getRangeMap()),
2277 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(transOp->getDomainMap()), transOp);
2278
2279 } else {
2280#ifdef TEKO_HAVE_EPETRA
2281 RCP<const Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2282 TEUCHOS_TEST_FOR_EXCEPTION(eOp == Teuchos::null, std::logic_error,
2283 "Teko::explicitTranspose Not an Epetra_Operator");
2284 RCP<const Epetra_RowMatrix> eRowMatrixOp =
2285 Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(eOp);
2286 TEUCHOS_TEST_FOR_EXCEPTION(eRowMatrixOp == Teuchos::null, std::logic_error,
2287 "Teko::explicitTranspose Not an Epetra_RowMatrix");
2288
2289 // we now have a delete transpose operator
2290 EpetraExt::RowMatrix_Transpose tranposeOp;
2291 Epetra_RowMatrix &eMat = tranposeOp(const_cast<Epetra_RowMatrix &>(*eRowMatrixOp));
2292
2293 // this copy is because of a poor implementation of the EpetraExt::Transform
2294 // implementation
2295 Teuchos::RCP<Epetra_CrsMatrix> crsMat =
2296 Teuchos::rcp(new Epetra_CrsMatrix(dynamic_cast<Epetra_CrsMatrix &>(eMat)));
2297
2298 return Thyra::epetraLinearOp(crsMat);
2299#else
2300 throw std::logic_error(
2301 "explicitTranspose is trying to use Epetra "
2302 "code, but TEKO_HAVE_EPETRA is disabled!");
2303#endif
2304 }
2305}
2306
2307const LinearOp explicitScale(double scalar, const LinearOp &op) {
2308 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2309 RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT>> tOp =
2310 rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op, true);
2311 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp =
2312 rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
2313 true);
2314 auto crsOpNew = rcp(new Tpetra::CrsMatrix<ST, LO, GO, NT>(*tCrsOp, Teuchos::Copy));
2315 crsOpNew->scale(scalar);
2316 return Thyra::tpetraLinearOp<ST, LO, GO, NT>(
2317 Thyra::createVectorSpace<ST, LO, GO, NT>(crsOpNew->getRangeMap()),
2318 Thyra::createVectorSpace<ST, LO, GO, NT>(crsOpNew->getDomainMap()), crsOpNew);
2319 } else {
2320#ifdef TEKO_HAVE_EPETRA
2321 RCP<const Thyra::EpetraLinearOp> eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op, true);
2322 RCP<const Epetra_CrsMatrix> eCrsOp =
2323 rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
2324 Teuchos::RCP<Epetra_CrsMatrix> crsMat = Teuchos::rcp(new Epetra_CrsMatrix(*eCrsOp));
2325
2326 crsMat->Scale(scalar);
2327
2328 return Thyra::epetraLinearOp(crsMat);
2329#else
2330 throw std::logic_error(
2331 "explicitScale is trying to use Epetra "
2332 "code, but TEKO_HAVE_EPETRA is disabled!");
2333#endif
2334 }
2335}
2336
2337double frobeniusNorm(const LinearOp &op_in) {
2338 LinearOp op;
2339 double scalar = 1.0;
2340
2341 // if blocked, must be 1x1
2342 if (isPhysicallyBlockedLinearOp(op_in)) {
2343 bool transp = false;
2344 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_op =
2345 getPhysicallyBlockedLinearOp(op_in, &scalar, &transp);
2346 TEUCHOS_ASSERT(blocked_op->productRange()->numBlocks() == 1);
2347 TEUCHOS_ASSERT(blocked_op->productDomain()->numBlocks() == 1);
2348 op = blocked_op->getBlock(0, 0);
2349 } else
2350 op = op_in;
2351
2352 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2353 const RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT>> tOp =
2354 rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op);
2355 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> crsOp =
2356 rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
2357 true);
2358 return crsOp->getFrobeniusNorm();
2359 } else {
2360#ifdef TEKO_HAVE_EPETRA
2361 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2362 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp, true);
2363 return crsOp->NormFrobenius();
2364#else
2365 throw std::logic_error(
2366 "frobeniusNorm is trying to use Epetra "
2367 "code, but TEKO_HAVE_EPETRA is disabled!");
2368#endif
2369 }
2370}
2371
2372double oneNorm(const LinearOp &op) {
2373 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2374 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
2375 "One norm not currently implemented for Tpetra matrices");
2376
2377 } else {
2378#ifdef TEKO_HAVE_EPETRA
2379 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2380 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp, true);
2381 return crsOp->NormOne();
2382#else
2383 throw std::logic_error(
2384 "oneNorm is trying to use Epetra "
2385 "code, but TEKO_HAVE_EPETRA is disabled!");
2386#endif
2387 }
2388}
2389
2390double infNorm(const LinearOp &op) {
2391 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2392 ST scalar = 0.0;
2393 bool transp = false;
2394 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> tCrsOp =
2395 Teko::TpetraHelpers::getTpetraCrsMatrix(op, &scalar, &transp);
2396
2397 // extract diagonal
2398 const RCP<Tpetra::Vector<ST, LO, GO, NT>> ptrDiag =
2399 Tpetra::createVector<ST, LO, GO, NT>(tCrsOp->getRowMap());
2400 Tpetra::Vector<ST, LO, GO, NT> &diag = *ptrDiag;
2401
2402 // compute absolute value row sum
2403 diag.putScalar(0.0);
2404 for (LO i = 0; i < (LO)tCrsOp->getLocalNumRows(); i++) {
2405 LO numEntries = tCrsOp->getNumEntriesInLocalRow(i);
2406 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::local_inds_host_view_type indices;
2407 typename Tpetra::CrsMatrix<ST, LO, GO, NT>::values_host_view_type values;
2408 tCrsOp->getLocalRowView(i, indices, values);
2409
2410 // build abs value row sum
2411 for (LO j = 0; j < numEntries; j++) diag.sumIntoLocalValue(i, std::abs(values(j)));
2412 }
2413 return diag.normInf() * scalar;
2414
2415 } else {
2416#ifdef TEKO_HAVE_EPETRA
2417 const RCP<const Epetra_Operator> epOp = Thyra::get_Epetra_Operator(*op);
2418 const RCP<const Epetra_CrsMatrix> crsOp = rcp_dynamic_cast<const Epetra_CrsMatrix>(epOp, true);
2419 return crsOp->NormInf();
2420#else
2421 throw std::logic_error(
2422 "infNorm is trying to use Epetra "
2423 "code, but TEKO_HAVE_EPETRA is disabled!");
2424#endif
2425 }
2426}
2427
2428const LinearOp buildDiagonal(const MultiVector &src, const std::string &lbl) {
2429 RCP<Thyra::VectorBase<double>> dst = Thyra::createMember(src->range());
2430 Thyra::copy(*src->col(0), dst.ptr());
2431
2432 return Thyra::diagonal<double>(dst, lbl);
2433}
2434
2435const LinearOp buildInvDiagonal(const MultiVector &src, const std::string &lbl) {
2436 const RCP<const Thyra::VectorBase<double>> srcV = src->col(0);
2437 RCP<Thyra::VectorBase<double>> dst = Thyra::createMember(srcV->range());
2438 Thyra::reciprocal<double>(*srcV, dst.ptr());
2439
2440 return Thyra::diagonal<double>(dst, lbl);
2441}
2442
2444BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> &mvv) {
2445 Teuchos::Array<MultiVector> mvA;
2446 Teuchos::Array<VectorSpace> vsA;
2447
2448 // build arrays of multi vectors and vector spaces
2449 std::vector<MultiVector>::const_iterator itr;
2450 for (itr = mvv.begin(); itr != mvv.end(); ++itr) {
2451 mvA.push_back(*itr);
2452 vsA.push_back((*itr)->range());
2453 }
2454
2455 // construct the product vector space
2456 const RCP<const Thyra::DefaultProductVectorSpace<double>> vs =
2457 Thyra::productVectorSpace<double>(vsA);
2458
2459 return Thyra::defaultProductMultiVector<double>(vs, mvA);
2460}
2461
2472Teuchos::RCP<Thyra::VectorBase<double>> indicatorVector(const std::vector<int> &indices,
2473 const VectorSpace &vs, double onValue,
2474 double offValue)
2475
2476{
2477 using Teuchos::RCP;
2478
2479 // create a new vector
2480 RCP<Thyra::VectorBase<double>> v = Thyra::createMember<double>(vs);
2481 Thyra::put_scalar<double>(offValue, v.ptr()); // fill it with "off" values
2482
2483 // set on values
2484 for (std::size_t i = 0; i < indices.size(); i++)
2485 Thyra::set_ele<double>(indices[i], onValue, v.ptr());
2486
2487 return v;
2488}
2489
2514double computeSpectralRad(const RCP<const Thyra::LinearOpBase<double>> &A, double tol,
2515 bool isHermitian, int numBlocks, int restart, int verbosity) {
2516 typedef Thyra::LinearOpBase<double> OP;
2517 typedef Thyra::MultiVectorBase<double> MV;
2518
2519 int startVectors = 1;
2520
2521 // construct an initial guess
2522 const RCP<MV> ivec = Thyra::createMember(A->domain());
2523 Thyra::randomize(-1.0, 1.0, ivec.ptr());
2524
2525 RCP<Anasazi::BasicEigenproblem<double, MV, OP>> eigProb =
2526 rcp(new Anasazi::BasicEigenproblem<double, MV, OP>(A, ivec));
2527 eigProb->setNEV(1);
2528 eigProb->setHermitian(isHermitian);
2529
2530 // set the problem up
2531 if (not eigProb->setProblem()) {
2532 // big time failure!
2533 return Teuchos::ScalarTraits<double>::nan();
2534 }
2535
2536 // we want largert magnitude eigenvalue
2537 std::string which("LM"); // largest magnitude
2538
2539 // Create the parameter list for the eigensolver
2540 // verbosity+=Anasazi::TimingDetails;
2541 Teuchos::ParameterList MyPL;
2542 MyPL.set("Verbosity", verbosity);
2543 MyPL.set("Which", which);
2544 MyPL.set("Block Size", startVectors);
2545 MyPL.set("Num Blocks", numBlocks);
2546 MyPL.set("Maximum Restarts", restart);
2547 MyPL.set("Convergence Tolerance", tol);
2548
2549 // build status test manager
2550 // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2551 // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(numBlocks*(restart+1)));
2552
2553 // Create the Block Krylov Schur solver
2554 // This takes as inputs the eigenvalue problem and the solver parameters
2555 Anasazi::BlockKrylovSchurSolMgr<double, MV, OP> MyBlockKrylovSchur(eigProb, MyPL);
2556
2557 // Solve the eigenvalue problem, and save the return code
2558 Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2559
2560 if (solverreturn == Anasazi::Unconverged) {
2561 double real = MyBlockKrylovSchur.getRitzValues().begin()->realpart;
2562 double comp = MyBlockKrylovSchur.getRitzValues().begin()->imagpart;
2563
2564 return -std::abs(std::sqrt(real * real + comp * comp));
2565
2566 // cout << "Anasazi::BlockKrylovSchur::solve() did not converge!" << std::endl;
2567 // return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2568 } else { // solverreturn==Anasazi::Converged
2569 double real = eigProb->getSolution().Evals.begin()->realpart;
2570 double comp = eigProb->getSolution().Evals.begin()->imagpart;
2571
2572 return std::abs(std::sqrt(real * real + comp * comp));
2573
2574 // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2575 // return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2576 }
2577}
2578
2602double computeSmallestMagEig(const RCP<const Thyra::LinearOpBase<double>> &A, double tol,
2603 bool isHermitian, int numBlocks, int restart, int verbosity) {
2604 typedef Thyra::LinearOpBase<double> OP;
2605 typedef Thyra::MultiVectorBase<double> MV;
2606
2607 int startVectors = 1;
2608
2609 // construct an initial guess
2610 const RCP<MV> ivec = Thyra::createMember(A->domain());
2611 Thyra::randomize(-1.0, 1.0, ivec.ptr());
2612
2613 RCP<Anasazi::BasicEigenproblem<double, MV, OP>> eigProb =
2614 rcp(new Anasazi::BasicEigenproblem<double, MV, OP>(A, ivec));
2615 eigProb->setNEV(1);
2616 eigProb->setHermitian(isHermitian);
2617
2618 // set the problem up
2619 if (not eigProb->setProblem()) {
2620 // big time failure!
2621 return Teuchos::ScalarTraits<double>::nan();
2622 }
2623
2624 // we want largert magnitude eigenvalue
2625 std::string which("SM"); // smallest magnitude
2626
2627 // Create the parameter list for the eigensolver
2628 Teuchos::ParameterList MyPL;
2629 MyPL.set("Verbosity", verbosity);
2630 MyPL.set("Which", which);
2631 MyPL.set("Block Size", startVectors);
2632 MyPL.set("Num Blocks", numBlocks);
2633 MyPL.set("Maximum Restarts", restart);
2634 MyPL.set("Convergence Tolerance", tol);
2635
2636 // build status test manager
2637 // RCP<Anasazi::StatusTestMaxIters<double,MV,OP> > statTest
2638 // = rcp(new Anasazi::StatusTestMaxIters<double,MV,OP>(10));
2639
2640 // Create the Block Krylov Schur solver
2641 // This takes as inputs the eigenvalue problem and the solver parameters
2642 Anasazi::BlockKrylovSchurSolMgr<double, MV, OP> MyBlockKrylovSchur(eigProb, MyPL);
2643
2644 // Solve the eigenvalue problem, and save the return code
2645 Anasazi::ReturnType solverreturn = MyBlockKrylovSchur.solve();
2646
2647 if (solverreturn == Anasazi::Unconverged) {
2648 // cout << "Anasazi::BlockKrylovSchur::solve() did not converge! " << std::endl;
2649 return -std::abs(MyBlockKrylovSchur.getRitzValues().begin()->realpart);
2650 } else { // solverreturn==Anasazi::Converged
2651 // cout << "Anasazi::BlockKrylovSchur::solve() converged!" << endl;
2652 return std::abs(eigProb->getSolution().Evals.begin()->realpart);
2653 }
2654}
2655
2664ModifiableLinearOp getDiagonalOp(const Teko::LinearOp &A, const DiagonalType &dt) {
2665 switch (dt) {
2666 case Diagonal: return getDiagonalOp(A);
2667 case Lumped: return getLumpedMatrix(A);
2668 case AbsRowSum: return getAbsRowSumMatrix(A);
2669 case NotDiag:
2670 default: TEUCHOS_TEST_FOR_EXCEPT(true);
2671 };
2672
2673 return Teuchos::null;
2674}
2675
2684ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp &A, const Teko::DiagonalType &dt) {
2685 switch (dt) {
2686 case Diagonal: return getInvDiagonalOp(A);
2687 case Lumped: return getInvLumpedMatrix(A);
2688 case AbsRowSum: return getAbsRowSumInvMatrix(A);
2689 case NotDiag:
2690 default: TEUCHOS_TEST_FOR_EXCEPT(true);
2691 };
2692
2693 return Teuchos::null;
2694}
2695
2702std::string getDiagonalName(const DiagonalType &dt) {
2703 switch (dt) {
2704 case Diagonal: return "Diagonal";
2705 case Lumped: return "Lumped";
2706 case AbsRowSum: return "AbsRowSum";
2707 case NotDiag: return "NotDiag";
2708 case BlkDiag: return "BlkDiag";
2709 };
2710
2711 return "<error>";
2712}
2713
2722DiagonalType getDiagonalType(std::string name) {
2723 if (name == "Diagonal") return Diagonal;
2724 if (name == "Lumped") return Lumped;
2725 if (name == "AbsRowSum") return AbsRowSum;
2726 if (name == "BlkDiag") return BlkDiag;
2727
2728 return NotDiag;
2729}
2730
2731#ifdef TEKO_HAVE_EPETRA
2732LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G, const LinearOp &Op) {
2733#ifdef Teko_ENABLE_Isorropia
2734 Teuchos::ParameterList probeList;
2735 Prober prober(G, probeList, true);
2736 Teuchos::RCP<Epetra_CrsMatrix> Mat = rcp(new Epetra_CrsMatrix(Copy, *G));
2737 Teko::Epetra::EpetraOperatorWrapper Mwrap(Op);
2738 prober.probe(Mwrap, *Mat);
2739 return Thyra::epetraLinearOp(Mat);
2740#else
2741 (void)G;
2742 (void)Op;
2743 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Probe requires Isorropia");
2744#endif
2745}
2746#endif
2747
2748double norm_1(const MultiVector &v, std::size_t col) {
2749 Teuchos::Array<double> n(v->domain()->dim());
2750 Thyra::norms_1<double>(*v, n);
2751
2752 return n[col];
2753}
2754
2755double norm_2(const MultiVector &v, std::size_t col) {
2756 Teuchos::Array<double> n(v->domain()->dim());
2757 Thyra::norms_2<double>(*v, n);
2758
2759 return n[col];
2760}
2761
2762#ifdef TEKO_HAVE_EPETRA
2763void putScalar(const ModifiableLinearOp &op, double scalar) {
2764 try {
2765 // get Epetra_Operator
2766 RCP<Epetra_Operator> eOp = Thyra::get_Epetra_Operator(*op);
2767
2768 // cast it to a CrsMatrix
2769 RCP<Epetra_CrsMatrix> eCrsOp = rcp_dynamic_cast<Epetra_CrsMatrix>(eOp, true);
2770
2771 eCrsOp->PutScalar(scalar);
2772 } catch (std::exception &e) {
2773 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
2774
2775 *out << "Teko: putScalar requires an Epetra_CrsMatrix\n";
2776 *out << " Could not extract an Epetra_Operator from a \"" << op->description() << std::endl;
2777 *out << " OR\n";
2778 *out << " Could not cast an Epetra_Operator to a Epetra_CrsMatrix\n";
2779 *out << std::endl;
2780 *out << "*** THROWN EXCEPTION ***\n";
2781 *out << e.what() << std::endl;
2782 *out << "************************\n";
2783
2784 throw e;
2785 }
2786}
2787#endif
2788
2789void clipLower(MultiVector &v, double lowerBound) {
2790 using Teuchos::RCP;
2791 using Teuchos::rcp_dynamic_cast;
2792
2793 // cast so entries are accessible
2794 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2795 // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2796
2797 for (Thyra::Ordinal i = 0; i < v->domain()->dim(); i++) {
2798 RCP<Thyra::SpmdVectorBase<double>> spmdVec =
2799 rcp_dynamic_cast<Thyra::SpmdVectorBase<double>>(v->col(i), true);
2800
2801 Teuchos::ArrayRCP<double> values;
2802 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2803 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2804 for (Teuchos::ArrayRCP<double>::size_type j = 0; j < values.size(); j++) {
2805 if (values[j] < lowerBound) values[j] = lowerBound;
2806 }
2807 }
2808}
2809
2810void clipUpper(MultiVector &v, double upperBound) {
2811 using Teuchos::RCP;
2812 using Teuchos::rcp_dynamic_cast;
2813
2814 // cast so entries are accessible
2815 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2816 // = rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(v);
2817 for (Thyra::Ordinal i = 0; i < v->domain()->dim(); i++) {
2818 RCP<Thyra::SpmdVectorBase<double>> spmdVec =
2819 rcp_dynamic_cast<Thyra::SpmdVectorBase<double>>(v->col(i), true);
2820
2821 Teuchos::ArrayRCP<double> values;
2822 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2823 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2824 for (Teuchos::ArrayRCP<double>::size_type j = 0; j < values.size(); j++) {
2825 if (values[j] > upperBound) values[j] = upperBound;
2826 }
2827 }
2828}
2829
2830void replaceValue(MultiVector &v, double currentValue, double newValue) {
2831 using Teuchos::RCP;
2832 using Teuchos::rcp_dynamic_cast;
2833
2834 // cast so entries are accessible
2835 // RCP<Thyra::SpmdMultiVectorBase<double> > spmdMVec
2836 // = rcp_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(v,true);
2837 for (Thyra::Ordinal i = 0; i < v->domain()->dim(); i++) {
2838 RCP<Thyra::SpmdVectorBase<double>> spmdVec =
2839 rcp_dynamic_cast<Thyra::SpmdVectorBase<double>>(v->col(i), true);
2840
2841 Teuchos::ArrayRCP<double> values;
2842 // spmdMVec->getNonconstLocalData(Teuchos::ptrFromRef(values),Teuchos::ptrFromRef(i));
2843 spmdVec->getNonconstLocalData(Teuchos::ptrFromRef(values));
2844 for (Teuchos::ArrayRCP<double>::size_type j = 0; j < values.size(); j++) {
2845 if (values[j] == currentValue) values[j] = newValue;
2846 }
2847 }
2848}
2849
2850void columnAverages(const MultiVector &v, std::vector<double> &averages) {
2851 averages.resize(v->domain()->dim());
2852
2853 // sum over each column
2854 Thyra::sums<double>(*v, averages);
2855
2856 // build averages
2857 Thyra::Ordinal rows = v->range()->dim();
2858 for (std::size_t i = 0; i < averages.size(); i++) averages[i] = averages[i] / rows;
2859}
2860
2861double average(const MultiVector &v) {
2862 Thyra::Ordinal rows = v->range()->dim();
2863 Thyra::Ordinal cols = v->domain()->dim();
2864
2865 std::vector<double> averages;
2866 columnAverages(v, averages);
2867
2868 double sum = 0.0;
2869 for (std::size_t i = 0; i < averages.size(); i++) sum += averages[i] * rows;
2870
2871 return sum / (rows * cols);
2872}
2873
2874bool isPhysicallyBlockedLinearOp(const LinearOp &op) {
2875 // See if the operator is a PBLO
2876 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> pblo =
2877 rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(op);
2878 if (!pblo.is_null()) return true;
2879
2880 // See if the operator is a wrapped PBLO
2881 ST scalar = 0.0;
2882 Thyra::EOpTransp transp = Thyra::NOTRANS;
2883 RCP<const Thyra::LinearOpBase<ST>> wrapped_op;
2884 Thyra::unwrap(op, &scalar, &transp, &wrapped_op);
2885 pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(wrapped_op);
2886 if (!pblo.is_null()) return true;
2887
2888 return false;
2889}
2890
2891RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> getPhysicallyBlockedLinearOp(
2892 const LinearOp &op, ST *scalar, bool *transp) {
2893 // If the operator is a TpetraLinearOp
2894 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> pblo =
2895 rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(op);
2896 if (!pblo.is_null()) {
2897 *scalar = 1.0;
2898 *transp = false;
2899 return pblo;
2900 }
2901
2902 // If the operator is a wrapped TpetraLinearOp
2903 RCP<const Thyra::LinearOpBase<ST>> wrapped_op;
2904 Thyra::EOpTransp eTransp = Thyra::NOTRANS;
2905 Thyra::unwrap(op, scalar, &eTransp, &wrapped_op);
2906 pblo = rcp_dynamic_cast<const Thyra::PhysicallyBlockedLinearOpBase<double>>(wrapped_op, true);
2907 if (!pblo.is_null()) {
2908 *transp = true;
2909 if (eTransp == Thyra::NOTRANS) *transp = false;
2910 return pblo;
2911 }
2912
2913 return Teuchos::null;
2914}
2915
2916std::string formatBlockName(const std::string &prefix, int i, int j, int nrow) {
2917 unsigned digits = 0;
2918 auto blockId = nrow - 1;
2919 do {
2920 blockId /= 10;
2921 digits++;
2922 } while (blockId);
2923
2924 std::ostringstream ss;
2925 ss << prefix << "_";
2926 ss << std::setfill('0') << std::setw(digits) << i;
2927 ss << "_";
2928 ss << std::setfill('0') << std::setw(digits) << j;
2929 ss << ".mm";
2930 return ss.str();
2931}
2932
2933void writeMatrix(const std::string &filename, const Teko::LinearOp &op) {
2934 using Teuchos::RCP;
2935 using Teuchos::rcp_dynamic_cast;
2936#ifdef TEKO_HAVE_EPETRA
2937 const RCP<const Thyra::EpetraLinearOp> eOp = rcp_dynamic_cast<const Thyra::EpetraLinearOp>(op);
2938#endif
2939
2940 if (Teko::TpetraHelpers::isTpetraLinearOp(op)) {
2941 const RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT>> tOp =
2942 rcp_dynamic_cast<const Thyra::TpetraLinearOp<ST, LO, GO, NT>>(op);
2943 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT>> crsOp =
2944 rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT>>(tOp->getConstTpetraOperator(),
2945 true);
2946 using Writer = Tpetra::MatrixMarket::Writer<Tpetra::CrsMatrix<ST, LO, GO, NT>>;
2947 Writer::writeMapFile(("rowmap_" + filename).c_str(), *(crsOp->getRowMap()));
2948 Writer::writeMapFile(("colmap_" + filename).c_str(), *(crsOp->getColMap()));
2949 Writer::writeMapFile(("domainmap_" + filename).c_str(), *(crsOp->getDomainMap()));
2950 Writer::writeMapFile(("rangemap_" + filename).c_str(), *(crsOp->getRangeMap()));
2951 Writer::writeSparseFile(filename.c_str(), crsOp);
2952 }
2953#ifdef TEKO_HAVE_EPETRA
2954 else if (eOp != Teuchos::null) {
2955 const RCP<const Epetra_CrsMatrix> crsOp =
2956 rcp_dynamic_cast<const Epetra_CrsMatrix>(eOp->epetra_op(), true);
2957 EpetraExt::BlockMapToMatrixMarketFile(("rowmap_" + filename).c_str(), crsOp->RowMap());
2958 EpetraExt::BlockMapToMatrixMarketFile(("colmap_" + filename).c_str(), crsOp->ColMap());
2959 EpetraExt::BlockMapToMatrixMarketFile(("domainmap_" + filename).c_str(), crsOp->DomainMap());
2960 EpetraExt::BlockMapToMatrixMarketFile(("rangemap_" + filename).c_str(), crsOp->RangeMap());
2961 EpetraExt::RowMatrixToMatrixMarketFile(filename.c_str(), *crsOp);
2962 }
2963#endif
2964 else if (isPhysicallyBlockedLinearOp(op)) {
2965 double scalar = 0.0;
2966 bool transp = false;
2967 RCP<const Thyra::PhysicallyBlockedLinearOpBase<double>> blocked_op =
2968 getPhysicallyBlockedLinearOp(op, &scalar, &transp);
2969
2970 int numRows = blocked_op->productRange()->numBlocks();
2971 int numCols = blocked_op->productDomain()->numBlocks();
2972
2973 for (int r = 0; r < numRows; ++r)
2974 for (int c = 0; c < numCols; ++c) {
2975 auto block = Teko::explicitScale(scalar, blocked_op->getBlock(r, c));
2976 if (transp) block = Teko::explicitTranspose(block);
2977 writeMatrix(formatBlockName(filename, r, c, numRows), block);
2978 }
2979 } else {
2980 TEUCHOS_ASSERT(false);
2981 }
2982}
2983
2984} // namespace Teko
@ BlkDiag
Specifies that a block diagonal approximation is to be used.
@ NotDiag
For user convenience, if Teko recieves this value, exceptions will be thrown.
@ AbsRowSum
Specifies that the diagonal entry is .
@ Diagonal
Specifies that just the diagonal is used.
@ Lumped
Specifies that row sum is used to form a diagonal.