18#ifndef __Teko_Utilities_hpp__
19#define __Teko_Utilities_hpp__
21#include "Teko_ConfigDefs.hpp"
23#ifdef TEKO_HAVE_EPETRA
24#include "Epetra_CrsMatrix.h"
27#include "Tpetra_CrsMatrix.hpp"
30#include "Teuchos_VerboseObject.hpp"
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"
50#ifndef _MSC_EXTENSIONS
51#define _MSC_EXTENSIONS
52#define TEKO_DEFINED_MSC_EXTENSIONS
58#define Teko_DEBUG_INT 5
89#ifdef TEKO_HAVE_EPETRA
90Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(
int dim,
double *coords,
91 const Epetra_CrsMatrix &stencil);
94Teuchos::RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > buildGraphLaplacian(
95 int dim, ST *coords,
const Tpetra::CrsMatrix<ST, LO, GO, NT> &stencil);
119#ifdef TEKO_HAVE_EPETRA
120Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(
double *x,
double *y,
double *z,
int stride,
121 const Epetra_CrsMatrix &stencil);
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);
133const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream();
137#ifndef Teko_DEBUG_OFF
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; \
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(); \
156#define Teko_DEBUG_PUSHTAB() Teko::getOutputStream()->pushTab(3)
157#define Teko_DEBUG_POPTAB() Teko::getOutputStream()->popTab()
158#define Teko_DEBUG_SCOPE(str, level)
169#define Teko_DEBUG_EXPR(str)
170#define Teko_DEBUG_MSG(str, level)
171#define Teko_DEBUG_MSG_BEGIN(level) \
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)
181typedef Teuchos::RCP<const Thyra::VectorSpaceBase<double> > VectorSpace;
188typedef Teuchos::RCP<Thyra::ProductMultiVectorBase<double> > BlockedMultiVector;
189typedef Teuchos::RCP<Thyra::MultiVectorBase<double> > MultiVector;
195inline const MultiVector
toMultiVector(
const BlockedMultiVector &bmv) {
return bmv; }
199 return Teuchos::rcp_dynamic_cast<Thyra::ProductMultiVectorBase<double> >(bmv);
203inline int blockCount(
const BlockedMultiVector &bmv) {
return bmv->productSpace()->numBlocks(); }
206inline MultiVector
getBlock(
int i,
const BlockedMultiVector &bmv) {
207 return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(bmv->getMultiVectorBlock(i));
211inline MultiVector
deepcopy(
const MultiVector &v) {
return v->clone_mv(); }
214inline MultiVector
copyAndInit(
const MultiVector &v,
double scalar) {
215 MultiVector mv = v->clone_mv();
216 Thyra::assign(mv.ptr(), scalar);
221inline BlockedMultiVector
deepcopy(
const BlockedMultiVector &v) {
238inline MultiVector
datacopy(
const MultiVector &src, MultiVector &dst) {
239 if (dst == Teuchos::null)
return deepcopy(src);
241 bool rangeCompat = src->range()->isCompatible(*dst->range());
242 bool domainCompat = src->domain()->isCompatible(*dst->domain());
244 if (not(rangeCompat && domainCompat))
return deepcopy(src);
247 Thyra::assign<double>(dst.ptr(), *src);
264inline BlockedMultiVector
datacopy(
const BlockedMultiVector &src, BlockedMultiVector &dst) {
265 if (dst == Teuchos::null)
return deepcopy(src);
267 bool rangeCompat = src->range()->isCompatible(*dst->range());
268 bool domainCompat = src->domain()->isCompatible(*dst->domain());
270 if (not(rangeCompat && domainCompat))
return deepcopy(src);
273 Thyra::assign<double>(dst.ptr(), *src);
278BlockedMultiVector buildBlockedMultiVector(
const std::vector<MultiVector> &mvs);
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);
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;
307inline LinearOp zero(
const VectorSpace &vs) {
return Thyra::zero<ST>(vs, vs); }
309inline LinearOp zero(
const VectorSpace &range,
const VectorSpace &domain) {
310 return Thyra::zero<ST>(range, domain);
314#ifdef TEKO_HAVE_EPETRA
315void putScalar(
const ModifiableLinearOp &op,
double scalar);
319inline VectorSpace rangeSpace(
const LinearOp &lo) {
return lo->range(); }
322inline VectorSpace domainSpace(
const LinearOp &lo) {
return lo->domain(); }
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);
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);
339inline LinearOp toLinearOp(BlockedLinearOp &blo) {
return blo; }
342inline const LinearOp toLinearOp(
const BlockedLinearOp &blo) {
return blo; }
345inline LinearOp toLinearOp(ModifiableLinearOp &blo) {
return blo; }
348inline const LinearOp toLinearOp(
const ModifiableLinearOp &blo) {
return blo; }
351inline int blockRowCount(
const BlockedLinearOp &blo) {
return blo->productRange()->numBlocks(); }
354inline int blockColCount(
const BlockedLinearOp &blo) {
return blo->productDomain()->numBlocks(); }
357inline LinearOp
getBlock(
int i,
int j,
const BlockedLinearOp &blo) {
return blo->getBlock(i, j); }
360inline void setBlock(
int i,
int j, BlockedLinearOp &blo,
const LinearOp &lo) {
361 return blo->setBlock(i, j, lo);
365inline BlockedLinearOp createBlockedOp() {
366 return rcp(
new Thyra::DefaultBlockedLinearOp<double>());
380inline void beginBlockFill(BlockedLinearOp &blo,
int rowCnt,
int colCnt) {
381 blo->beginBlockFill(rowCnt, colCnt);
393inline void beginBlockFill(BlockedLinearOp &blo) { blo->beginBlockFill(); }
396inline void endBlockFill(BlockedLinearOp &blo) { blo->endBlockFill(); }
399BlockedLinearOp getUpperTriBlocks(
const BlockedLinearOp &blo,
bool callEndBlockFill =
true);
402BlockedLinearOp getLowerTriBlocks(
const BlockedLinearOp &blo,
bool callEndBlockFill =
true);
423BlockedLinearOp zeroBlockedOp(
const BlockedLinearOp &blo);
426bool isZeroOp(
const LinearOp op);
436ModifiableLinearOp getAbsRowSumMatrix(
const LinearOp &op);
446ModifiableLinearOp getAbsRowSumInvMatrix(
const LinearOp &op);
455ModifiableLinearOp getLumpedMatrix(
const LinearOp &op);
465ModifiableLinearOp getInvLumpedMatrix(
const LinearOp &op);
490void applyOp(
const LinearOp &A,
const MultiVector &x, MultiVector &y,
double alpha = 1.0,
511void applyTransposeOp(
const LinearOp &A,
const MultiVector &x, MultiVector &y,
double alpha = 1.0,
532inline void applyOp(
const LinearOp &A,
const BlockedMultiVector &x, BlockedMultiVector &y,
533 double alpha = 1.0,
double beta = 0.0) {
536 applyOp(A, x_mv, y_mv, alpha, beta);
557inline void applyTransposeOp(
const LinearOp &A,
const BlockedMultiVector &x, BlockedMultiVector &y,
558 double alpha = 1.0,
double beta = 0.0) {
561 applyTransposeOp(A, x_mv, y_mv, alpha, beta);
573void update(
double alpha,
const MultiVector &x,
double beta, MultiVector &y);
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);
583inline void scale(
const double alpha, MultiVector &x) { Thyra::scale<double>(alpha, x.ptr()); }
586inline void scale(
const double alpha, BlockedMultiVector &x) {
592inline LinearOp scale(
const double alpha, ModifiableLinearOp &a) {
593 return Thyra::nonconstScale(alpha, a);
597inline LinearOp adjoint(ModifiableLinearOp &a) {
return Thyra::nonconstAdjoint(a); }
615const ModifiableLinearOp getDiagonalOp(
const LinearOp &op);
622const MultiVector getDiagonal(
const LinearOp &op);
635const ModifiableLinearOp getInvDiagonalOp(
const LinearOp &op);
649const LinearOp explicitMultiply(
const LinearOp &opl,
const LinearOp &opm,
const LinearOp &opr);
665const ModifiableLinearOp explicitMultiply(
const LinearOp &opl,
const LinearOp &opm,
666 const LinearOp &opr,
const ModifiableLinearOp &destOp);
678const LinearOp explicitMultiply(
const LinearOp &opl,
const LinearOp &opr);
693const ModifiableLinearOp explicitMultiply(
const LinearOp &opl,
const LinearOp &opr,
694 const ModifiableLinearOp &destOp);
706const LinearOp explicitAdd(
const LinearOp &opl,
const LinearOp &opr);
720const ModifiableLinearOp explicitAdd(
const LinearOp &opl,
const LinearOp &opr,
721 const ModifiableLinearOp &destOp);
725const ModifiableLinearOp explicitSum(
const LinearOp &opl,
const ModifiableLinearOp &destOp);
730const LinearOp explicitTranspose(
const LinearOp &op);
734const LinearOp explicitScale(
double scalar,
const LinearOp &op);
738double frobeniusNorm(
const LinearOp &op);
739double oneNorm(
const LinearOp &op);
740double infNorm(
const LinearOp &op);
745const LinearOp buildDiagonal(
const MultiVector &v,
const std::string &lbl =
"ANYM");
750const LinearOp buildInvDiagonal(
const MultiVector &v,
const std::string &lbl =
"ANYM");
777double computeSpectralRad(
const Teuchos::RCP<
const Thyra::LinearOpBase<double> > &A,
double tol,
778 bool isHermitian =
false,
int numBlocks = 5,
int restart = 0,
804double computeSmallestMagEig(
const Teuchos::RCP<
const Thyra::LinearOpBase<double> > &A,
double tol,
805 bool isHermitian =
false,
int numBlocks = 5,
int restart = 0,
829ModifiableLinearOp getDiagonalOp(
const Teko::LinearOp &A,
const DiagonalType &dt);
839ModifiableLinearOp getInvDiagonalOp(
const Teko::LinearOp &A,
const DiagonalType &dt);
846const MultiVector getDiagonal(
const LinearOp &op,
const DiagonalType &dt);
854std::string getDiagonalName(
const DiagonalType &dt);
864DiagonalType getDiagonalType(std::string name);
866#ifdef TEKO_HAVE_EPETRA
867LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,
const LinearOp &Op);
872double norm_1(
const MultiVector &v, std::size_t col);
876double norm_2(
const MultiVector &v, std::size_t col);
881void clipLower(MultiVector &v,
double lowerBound);
886void clipUpper(MultiVector &v,
double upperBound);
891void replaceValue(MultiVector &v,
double currentValue,
double newValue);
895void columnAverages(
const MultiVector &v, std::vector<double> &averages);
899double average(
const MultiVector &v);
903bool isPhysicallyBlockedLinearOp(
const LinearOp &op);
907Teuchos::RCP<const Thyra::PhysicallyBlockedLinearOpBase<double> > getPhysicallyBlockedLinearOp(
908 const LinearOp &op, ST *scalar,
bool *transp);
911std::string formatBlockName(
const std::string &prefix,
int i,
int j,
int nrow);
914void writeMatrix(
const std::string &filename,
const Teko::LinearOp &op);
919#ifdef TEKO_DEFINED_MSC_EXTENSIONS
920#undef _MSC_EXTENSIONS
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.