Teko Version of the Day
Loading...
Searching...
No Matches
Teko_Utilities.hpp
Go to the documentation of this file.
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
17
18#ifndef __Teko_Utilities_hpp__
19#define __Teko_Utilities_hpp__
20
21#include "Teko_ConfigDefs.hpp"
22
23#ifdef TEKO_HAVE_EPETRA
24#include "Epetra_CrsMatrix.h"
25#endif
26
27#include "Tpetra_CrsMatrix.hpp"
28
29// Teuchos includes
30#include "Teuchos_VerboseObject.hpp"
31
32// Thyra includes
33#include "Thyra_LinearOpBase.hpp"
34#include "Thyra_PhysicallyBlockedLinearOpBase.hpp"
35#include "Thyra_ProductVectorSpaceBase.hpp"
36#include "Thyra_VectorSpaceBase.hpp"
37#include "Thyra_ProductMultiVectorBase.hpp"
38#include "Thyra_MultiVectorStdOps.hpp"
39#include "Thyra_MultiVectorBase.hpp"
40#include "Thyra_VectorBase.hpp"
41#include "Thyra_VectorStdOps.hpp"
42#include "Thyra_DefaultBlockedLinearOp.hpp"
43#include "Thyra_DefaultMultipliedLinearOp.hpp"
44#include "Thyra_DefaultScaledAdjointLinearOp.hpp"
45#include "Thyra_DefaultAddedLinearOp.hpp"
46#include "Thyra_DefaultIdentityLinearOp.hpp"
47#include "Thyra_DefaultZeroLinearOp.hpp"
48
49#ifdef _MSC_VER
50#ifndef _MSC_EXTENSIONS
51#define _MSC_EXTENSIONS
52#define TEKO_DEFINED_MSC_EXTENSIONS
53#endif
54#include <iso646.h> // For C alternative tokens
55#endif
56
57// #define Teko_DEBUG_OFF
58#define Teko_DEBUG_INT 5
59
60namespace Teko {
61
62using Thyra::add;
63using Thyra::block1x2;
64using Thyra::block2x1;
65using Thyra::block2x2;
66using Thyra::identity;
67using Thyra::multiply;
68using Thyra::scale;
69using Thyra::zero; // make it to take one argument (square matrix)
70
89#ifdef TEKO_HAVE_EPETRA
90Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(int dim, double *coords,
91 const Epetra_CrsMatrix &stencil);
92#endif
93
94Teuchos::RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > buildGraphLaplacian(
95 int dim, ST *coords, const Tpetra::CrsMatrix<ST, LO, GO, NT> &stencil);
96
119#ifdef TEKO_HAVE_EPETRA
120Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(double *x, double *y, double *z, int stride,
121 const Epetra_CrsMatrix &stencil);
122#endif
123
124Teuchos::RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > buildGraphLaplacian(
125 ST *x, ST *y, ST *z, GO stride, const Tpetra::CrsMatrix<ST, LO, GO, NT> &stencil);
126
133const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream();
134// inline const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream();
135// { return Teuchos::VerboseObjectBase::getDefaultOStream(); }
136
137#ifndef Teko_DEBUG_OFF
138//#if 0
139#define Teko_DEBUG_EXPR(str) str
140#define Teko_DEBUG_MSG(str, level) \
141 if (level <= Teko_DEBUG_INT) { \
142 Teuchos::RCP<Teuchos::FancyOStream> out = Teko::getOutputStream(); \
143 *out << "Teko: " << str << std::endl; \
144 }
145#define Teko_DEBUG_MSG_BEGIN(level) \
146 if (level <= Teko_DEBUG_INT) { \
147 Teko::getOutputStream()->pushTab(3); \
148 *Teko::getOutputStream() << "Teko: Begin debug MSG\n"; \
149 std::ostream &DEBUG_STREAM = *Teko::getOutputStream(); \
150 Teko::getOutputStream()->pushTab(3);
151#define Teko_DEBUG_MSG_END() \
152 Teko::getOutputStream()->popTab(); \
153 *Teko::getOutputStream() << "Teko: End debug MSG\n"; \
154 Teko::getOutputStream()->popTab(); \
155 }
156#define Teko_DEBUG_PUSHTAB() Teko::getOutputStream()->pushTab(3)
157#define Teko_DEBUG_POPTAB() Teko::getOutputStream()->popTab()
158#define Teko_DEBUG_SCOPE(str, level)
159
160// struct __DebugScope__ {
161// __DebugScope__(const std::string & str,int level)
162// : str_(str), level_(level)
163// { Teko_DEBUG_MSG("BEGIN "+str_,level_); Teko_DEBUG_PUSHTAB(); }
164// ~__DebugScope__()
165// { Teko_DEBUG_POPTAB(); Teko_DEBUG_MSG("END "+str_,level_); }
166// std::string str_; int level_; };
167// #define Teko_DEBUG_SCOPE(str,level) __DebugScope__ __dbgScope__(str,level);
168#else
169#define Teko_DEBUG_EXPR(str)
170#define Teko_DEBUG_MSG(str, level)
171#define Teko_DEBUG_MSG_BEGIN(level) \
172 if (false) { \
173 std::ostream &DEBUG_STREAM = *Teko::getOutputStream();
174#define Teko_DEBUG_MSG_END() }
175#define Teko_DEBUG_PUSHTAB()
176#define Teko_DEBUG_POPTAB()
177#define Teko_DEBUG_SCOPE(str, level)
178#endif
179
180// typedefs for increased simplicity
181typedef Teuchos::RCP<const Thyra::VectorSpaceBase<double> > VectorSpace;
182
183// ----------------------------------------------------------------------------
184
186
187
188typedef Teuchos::RCP<Thyra::ProductMultiVectorBase<double> > BlockedMultiVector;
189typedef Teuchos::RCP<Thyra::MultiVectorBase<double> > MultiVector;
190
192inline MultiVector toMultiVector(BlockedMultiVector &bmv) { return bmv; }
193
195inline const MultiVector toMultiVector(const BlockedMultiVector &bmv) { return bmv; }
196
198inline const BlockedMultiVector toBlockedMultiVector(const MultiVector &bmv) {
199 return Teuchos::rcp_dynamic_cast<Thyra::ProductMultiVectorBase<double> >(bmv);
200}
201
203inline int blockCount(const BlockedMultiVector &bmv) { return bmv->productSpace()->numBlocks(); }
204
206inline MultiVector getBlock(int i, const BlockedMultiVector &bmv) {
207 return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(bmv->getMultiVectorBlock(i));
208}
209
211inline MultiVector deepcopy(const MultiVector &v) { return v->clone_mv(); }
212
214inline MultiVector copyAndInit(const MultiVector &v, double scalar) {
215 MultiVector mv = v->clone_mv();
216 Thyra::assign(mv.ptr(), scalar);
217 return mv;
218}
219
221inline BlockedMultiVector deepcopy(const BlockedMultiVector &v) {
222 return toBlockedMultiVector(v->clone_mv());
223}
224
238inline MultiVector datacopy(const MultiVector &src, MultiVector &dst) {
239 if (dst == Teuchos::null) return deepcopy(src);
240
241 bool rangeCompat = src->range()->isCompatible(*dst->range());
242 bool domainCompat = src->domain()->isCompatible(*dst->domain());
243
244 if (not(rangeCompat && domainCompat)) return deepcopy(src);
245
246 // perform data copy
247 Thyra::assign<double>(dst.ptr(), *src);
248 return dst;
249}
250
264inline BlockedMultiVector datacopy(const BlockedMultiVector &src, BlockedMultiVector &dst) {
265 if (dst == Teuchos::null) return deepcopy(src);
266
267 bool rangeCompat = src->range()->isCompatible(*dst->range());
268 bool domainCompat = src->domain()->isCompatible(*dst->domain());
269
270 if (not(rangeCompat && domainCompat)) return deepcopy(src);
271
272 // perform data copy
273 Thyra::assign<double>(dst.ptr(), *src);
274 return dst;
275}
276
278BlockedMultiVector buildBlockedMultiVector(const std::vector<MultiVector> &mvs);
279
290Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(const std::vector<int> &indices,
291 const VectorSpace &vs,
292 double onValue = 1.0,
293 double offValue = 0.0);
294
296
297// ----------------------------------------------------------------------------
298
300
301typedef Teuchos::RCP<Thyra::PhysicallyBlockedLinearOpBase<ST> > BlockedLinearOp;
302typedef Teuchos::RCP<const Thyra::LinearOpBase<ST> > LinearOp;
303typedef Teuchos::RCP<Thyra::LinearOpBase<ST> > InverseLinearOp;
304typedef Teuchos::RCP<Thyra::LinearOpBase<ST> > ModifiableLinearOp;
305
307inline LinearOp zero(const VectorSpace &vs) { return Thyra::zero<ST>(vs, vs); }
308
309inline LinearOp zero(const VectorSpace &range, const VectorSpace &domain) {
310 return Thyra::zero<ST>(range, domain);
311}
312
314#ifdef TEKO_HAVE_EPETRA
315void putScalar(const ModifiableLinearOp &op, double scalar);
316#endif
317
319inline VectorSpace rangeSpace(const LinearOp &lo) { return lo->range(); }
320
322inline VectorSpace domainSpace(const LinearOp &lo) { return lo->domain(); }
323
325inline BlockedLinearOp toBlockedLinearOp(LinearOp &clo) {
326 Teuchos::RCP<Thyra::LinearOpBase<double> > lo =
327 Teuchos::rcp_const_cast<Thyra::LinearOpBase<double> >(clo);
328 return Teuchos::rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(lo);
329}
330
332inline const BlockedLinearOp toBlockedLinearOp(const LinearOp &clo) {
333 Teuchos::RCP<Thyra::LinearOpBase<double> > lo =
334 Teuchos::rcp_const_cast<Thyra::LinearOpBase<double> >(clo);
335 return Teuchos::rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(lo);
336}
337
339inline LinearOp toLinearOp(BlockedLinearOp &blo) { return blo; }
340
342inline const LinearOp toLinearOp(const BlockedLinearOp &blo) { return blo; }
343
345inline LinearOp toLinearOp(ModifiableLinearOp &blo) { return blo; }
346
348inline const LinearOp toLinearOp(const ModifiableLinearOp &blo) { return blo; }
349
351inline int blockRowCount(const BlockedLinearOp &blo) { return blo->productRange()->numBlocks(); }
352
354inline int blockColCount(const BlockedLinearOp &blo) { return blo->productDomain()->numBlocks(); }
355
357inline LinearOp getBlock(int i, int j, const BlockedLinearOp &blo) { return blo->getBlock(i, j); }
358
360inline void setBlock(int i, int j, BlockedLinearOp &blo, const LinearOp &lo) {
361 return blo->setBlock(i, j, lo);
362}
363
365inline BlockedLinearOp createBlockedOp() {
366 return rcp(new Thyra::DefaultBlockedLinearOp<double>());
367}
368
380inline void beginBlockFill(BlockedLinearOp &blo, int rowCnt, int colCnt) {
381 blo->beginBlockFill(rowCnt, colCnt);
382}
383
393inline void beginBlockFill(BlockedLinearOp &blo) { blo->beginBlockFill(); }
394
396inline void endBlockFill(BlockedLinearOp &blo) { blo->endBlockFill(); }
397
399BlockedLinearOp getUpperTriBlocks(const BlockedLinearOp &blo, bool callEndBlockFill = true);
400
402BlockedLinearOp getLowerTriBlocks(const BlockedLinearOp &blo, bool callEndBlockFill = true);
403
423BlockedLinearOp zeroBlockedOp(const BlockedLinearOp &blo);
424
426bool isZeroOp(const LinearOp op);
427
436ModifiableLinearOp getAbsRowSumMatrix(const LinearOp &op);
437
446ModifiableLinearOp getAbsRowSumInvMatrix(const LinearOp &op);
447
455ModifiableLinearOp getLumpedMatrix(const LinearOp &op);
456
465ModifiableLinearOp getInvLumpedMatrix(const LinearOp &op);
466
468
470
471
490void applyOp(const LinearOp &A, const MultiVector &x, MultiVector &y, double alpha = 1.0,
491 double beta = 0.0);
492
511void applyTransposeOp(const LinearOp &A, const MultiVector &x, MultiVector &y, double alpha = 1.0,
512 double beta = 0.0);
513
532inline void applyOp(const LinearOp &A, const BlockedMultiVector &x, BlockedMultiVector &y,
533 double alpha = 1.0, double beta = 0.0) {
534 const MultiVector x_mv = toMultiVector(x);
535 MultiVector y_mv = toMultiVector(y);
536 applyOp(A, x_mv, y_mv, alpha, beta);
537}
538
557inline void applyTransposeOp(const LinearOp &A, const BlockedMultiVector &x, BlockedMultiVector &y,
558 double alpha = 1.0, double beta = 0.0) {
559 const MultiVector x_mv = toMultiVector(x);
560 MultiVector y_mv = toMultiVector(y);
561 applyTransposeOp(A, x_mv, y_mv, alpha, beta);
562}
563
573void update(double alpha, const MultiVector &x, double beta, MultiVector &y);
574
576inline void update(double alpha, const BlockedMultiVector &x, double beta, BlockedMultiVector &y) {
577 MultiVector x_mv = toMultiVector(x);
578 MultiVector y_mv = toMultiVector(y);
579 update(alpha, x_mv, beta, y_mv);
580}
581
583inline void scale(const double alpha, MultiVector &x) { Thyra::scale<double>(alpha, x.ptr()); }
584
586inline void scale(const double alpha, BlockedMultiVector &x) {
587 MultiVector x_mv = toMultiVector(x);
588 scale(alpha, x_mv);
589}
590
592inline LinearOp scale(const double alpha, ModifiableLinearOp &a) {
593 return Thyra::nonconstScale(alpha, a);
594}
595
597inline LinearOp adjoint(ModifiableLinearOp &a) { return Thyra::nonconstAdjoint(a); }
598
600
602
603
615const ModifiableLinearOp getDiagonalOp(const LinearOp &op);
616
622const MultiVector getDiagonal(const LinearOp &op);
623
635const ModifiableLinearOp getInvDiagonalOp(const LinearOp &op);
636
649const LinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opm, const LinearOp &opr);
650
665const ModifiableLinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opm,
666 const LinearOp &opr, const ModifiableLinearOp &destOp);
667
678const LinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opr);
679
693const ModifiableLinearOp explicitMultiply(const LinearOp &opl, const LinearOp &opr,
694 const ModifiableLinearOp &destOp);
695
706const LinearOp explicitAdd(const LinearOp &opl, const LinearOp &opr);
707
720const ModifiableLinearOp explicitAdd(const LinearOp &opl, const LinearOp &opr,
721 const ModifiableLinearOp &destOp);
722
725const ModifiableLinearOp explicitSum(const LinearOp &opl, const ModifiableLinearOp &destOp);
726
730const LinearOp explicitTranspose(const LinearOp &op);
731
734const LinearOp explicitScale(double scalar, const LinearOp &op);
735
738double frobeniusNorm(const LinearOp &op);
739double oneNorm(const LinearOp &op);
740double infNorm(const LinearOp &op);
741
745const LinearOp buildDiagonal(const MultiVector &v, const std::string &lbl = "ANYM");
746
750const LinearOp buildInvDiagonal(const MultiVector &v, const std::string &lbl = "ANYM");
751
753
777double computeSpectralRad(const Teuchos::RCP<const Thyra::LinearOpBase<double> > &A, double tol,
778 bool isHermitian = false, int numBlocks = 5, int restart = 0,
779 int verbosity = 0);
780
804double computeSmallestMagEig(const Teuchos::RCP<const Thyra::LinearOpBase<double> > &A, double tol,
805 bool isHermitian = false, int numBlocks = 5, int restart = 0,
806 int verbosity = 0);
807
809typedef enum {
811 ,
813 ,
815 ,
817 ,
819} DiagonalType;
820
829ModifiableLinearOp getDiagonalOp(const Teko::LinearOp &A, const DiagonalType &dt);
830
839ModifiableLinearOp getInvDiagonalOp(const Teko::LinearOp &A, const DiagonalType &dt);
840
846const MultiVector getDiagonal(const LinearOp &op, const DiagonalType &dt);
847
854std::string getDiagonalName(const DiagonalType &dt);
855
864DiagonalType getDiagonalType(std::string name);
865
866#ifdef TEKO_HAVE_EPETRA
867LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G, const LinearOp &Op);
868#endif
869
872double norm_1(const MultiVector &v, std::size_t col);
873
876double norm_2(const MultiVector &v, std::size_t col);
877
881void clipLower(MultiVector &v, double lowerBound);
882
886void clipUpper(MultiVector &v, double upperBound);
887
891void replaceValue(MultiVector &v, double currentValue, double newValue);
892
895void columnAverages(const MultiVector &v, std::vector<double> &averages);
896
899double average(const MultiVector &v);
900
903bool isPhysicallyBlockedLinearOp(const LinearOp &op);
904
907Teuchos::RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(
908 const LinearOp &op, ST *scalar, bool *transp);
909
911std::string formatBlockName(const std::string &prefix, int i, int j, int nrow);
912
914void writeMatrix(const std::string &filename, const Teko::LinearOp &op);
915
916} // end namespace Teko
917
918#ifdef _MSC_VER
919#ifdef TEKO_DEFINED_MSC_EXTENSIONS
920#undef _MSC_EXTENSIONS
921#endif
922#endif
923
924#endif
MultiVector datacopy(const MultiVector &src, MultiVector &dst)
Copy the contents of a multivector to a destination vector.
@ 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.
const BlockedMultiVector toBlockedMultiVector(const MultiVector &bmv)
Convert to a BlockedMultiVector from a MultiVector.
MultiVector toMultiVector(BlockedMultiVector &bmv)
Convert to a MultiVector from a BlockedMultiVector.
int blockCount(const BlockedMultiVector &bmv)
Get the column count in a block linear operator.
MultiVector getBlock(int i, const BlockedMultiVector &bmv)
Get the ith block from a BlockedMultiVector object.
MultiVector copyAndInit(const MultiVector &v, double scalar)
Perform a deep copy of the vector.
MultiVector deepcopy(const MultiVector &v)
Perform a deep copy of the vector.