55 #ifndef __Teko_Utilities_hpp__ 56 #define __Teko_Utilities_hpp__ 58 #include "Epetra_CrsMatrix.h" 61 #include "Teuchos_VerboseObject.hpp" 64 #include "Thyra_LinearOpBase.hpp" 65 #include "Thyra_PhysicallyBlockedLinearOpBase.hpp" 66 #include "Thyra_ProductVectorSpaceBase.hpp" 67 #include "Thyra_VectorSpaceBase.hpp" 68 #include "Thyra_ProductMultiVectorBase.hpp" 69 #include "Thyra_MultiVectorStdOps.hpp" 70 #include "Thyra_MultiVectorBase.hpp" 71 #include "Thyra_VectorBase.hpp" 72 #include "Thyra_VectorStdOps.hpp" 73 #include "Thyra_DefaultBlockedLinearOp.hpp" 74 #include "Thyra_DefaultMultipliedLinearOp.hpp" 75 #include "Thyra_DefaultScaledAdjointLinearOp.hpp" 76 #include "Thyra_DefaultAddedLinearOp.hpp" 77 #include "Thyra_DefaultIdentityLinearOp.hpp" 78 #include "Thyra_DefaultZeroLinearOp.hpp" 81 #ifndef _MSC_EXTENSIONS 82 #define _MSC_EXTENSIONS 83 #define TEKO_DEFINED_MSC_EXTENSIONS 89 #define Teko_DEBUG_INT 5 93 using Thyra::multiply;
96 using Thyra::identity;
98 using Thyra::block2x2;
99 using Thyra::block2x1;
100 using Thyra::block1x2;
120 Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(
int dim,
double * coords,
const Epetra_CrsMatrix & stencil);
144 Teuchos::RCP<Epetra_CrsMatrix> buildGraphLaplacian(
double * x,
double * y,
double * z,
int stride,
const Epetra_CrsMatrix & stencil);
152 const Teuchos::RCP<Teuchos::FancyOStream> getOutputStream();
156 #ifndef Teko_DEBUG_OFF 157 #define Teko_DEBUG_EXPR(str) str 158 #define Teko_DEBUG_MSG(str,level) if(level<=Teko_DEBUG_INT) { \ 159 Teuchos::RCP<Teuchos::FancyOStream> out = Teko::getOutputStream(); \ 160 *out << "Teko: " << str << std::endl; } 161 #define Teko_DEBUG_MSG_BEGIN(level) if(level<=Teko_DEBUG_INT) { \ 162 Teko::getOutputStream()->pushTab(3); \ 163 *Teko::getOutputStream() << "Teko: Begin debug MSG\n"; \ 164 std::ostream & DEBUG_STREAM = *Teko::getOutputStream(); \ 165 Teko::getOutputStream()->pushTab(3); 166 #define Teko_DEBUG_MSG_END() Teko::getOutputStream()->popTab(); \ 167 *Teko::getOutputStream() << "Teko: End debug MSG\n"; \ 168 Teko::getOutputStream()->popTab(); } 169 #define Teko_DEBUG_PUSHTAB() Teko::getOutputStream()->pushTab(3) 170 #define Teko_DEBUG_POPTAB() Teko::getOutputStream()->popTab() 172 struct __DebugScope__ {
173 __DebugScope__(
const std::string & str,
int level)
174 : str_(str), level_(level)
175 { Teko_DEBUG_MSG(
"BEGIN "+str_,level_); Teko_DEBUG_PUSHTAB(); }
177 { Teko_DEBUG_POPTAB(); Teko_DEBUG_MSG(
"END "+str_,level_); }
178 std::string str_;
int level_; };
179 #define Teko_DEBUG_SCOPE(str,level) __DebugScope__ __dbgScope__(str,level); 181 #define Teko_DEBUG_EXPR(str) 182 #define Teko_DEBUG_MSG(str,level) 183 #define Teko_DEBUG_MSG_BEGIN(level) if(false) { \ 184 std::ostream & DEBUG_STREAM = *Teko::getOutputStream(); 185 #define Teko_DEBUG_MSG_END() } 186 #define Teko_DEBUG_PUSHTAB() 187 #define Teko_DEBUG_POPTAB() 188 #define Teko_DEBUG_SCOPE(str,level) 192 typedef Teuchos::RCP<const Thyra::VectorSpaceBase<double> > VectorSpace;
199 typedef Teuchos::RCP<Thyra::ProductMultiVectorBase<double> > BlockedMultiVector;
200 typedef Teuchos::RCP<Thyra::MultiVectorBase<double> > MultiVector;
206 inline const MultiVector
toMultiVector(
const BlockedMultiVector & bmv) {
return bmv; }
210 {
return Teuchos::rcp_dynamic_cast<Thyra::ProductMultiVectorBase<double> >(bmv); }
214 {
return bmv->productSpace()->numBlocks(); }
217 inline MultiVector
getBlock(
int i,
const BlockedMultiVector & bmv)
218 {
return Teuchos::rcp_const_cast<Thyra::MultiVectorBase<double> >(bmv->getMultiVectorBlock(i)); }
222 {
return v->clone_mv(); }
225 inline MultiVector
copyAndInit(
const MultiVector & v,
double scalar)
226 { MultiVector mv = v->clone_mv(); Thyra::assign(mv.ptr(),scalar);
return mv; }
229 inline BlockedMultiVector
deepcopy(
const BlockedMultiVector & v)
245 inline MultiVector
datacopy(
const MultiVector & src,MultiVector & dst)
247 if(dst==Teuchos::null)
250 bool rangeCompat = src->range()->isCompatible(*dst->range());
251 bool domainCompat = src->domain()->isCompatible(*dst->domain());
253 if(not (rangeCompat && domainCompat))
257 Thyra::assign<double>(dst.ptr(),*src);
274 inline BlockedMultiVector
datacopy(
const BlockedMultiVector & src,BlockedMultiVector & dst)
276 if(dst==Teuchos::null)
279 bool rangeCompat = src->range()->isCompatible(*dst->range());
280 bool domainCompat = src->domain()->isCompatible(*dst->domain());
282 if(not (rangeCompat && domainCompat))
286 Thyra::assign<double>(dst.ptr(),*src);
291 BlockedMultiVector buildBlockedMultiVector(
const std::vector<MultiVector> & mvs);
303 Teuchos::RCP<Thyra::VectorBase<double> > indicatorVector(
304 const std::vector<int> & indices,
305 const VectorSpace & vs,
306 double onValue=1.0,
double offValue=0.0);
314 typedef Teuchos::RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > BlockedLinearOp;
315 typedef Teuchos::RCP<const Thyra::LinearOpBase<double> > LinearOp;
316 typedef Teuchos::RCP<Thyra::LinearOpBase<double> > InverseLinearOp;
317 typedef Teuchos::RCP<Thyra::LinearOpBase<double> > ModifiableLinearOp;
320 inline LinearOp
zero(
const VectorSpace & vs)
321 {
return Thyra::zero<double>(vs,vs); }
324 void putScalar(
const ModifiableLinearOp & op,
double scalar);
328 {
return lo->range(); }
332 {
return lo->domain(); }
337 Teuchos::RCP<Thyra::LinearOpBase<double> > lo = Teuchos::rcp_const_cast<Thyra::LinearOpBase<double> >(clo);
338 return Teuchos::rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> > (lo);
344 Teuchos::RCP<Thyra::LinearOpBase<double> > lo = Teuchos::rcp_const_cast<Thyra::LinearOpBase<double> >(clo);
345 return Teuchos::rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> > (lo);
349 inline LinearOp
toLinearOp(BlockedLinearOp & blo) {
return blo; }
352 inline const LinearOp
toLinearOp(
const BlockedLinearOp & blo) {
return blo; }
355 inline LinearOp
toLinearOp(ModifiableLinearOp & blo) {
return blo; }
358 inline const LinearOp
toLinearOp(
const ModifiableLinearOp & blo) {
return blo; }
362 {
return blo->productRange()->numBlocks(); }
366 {
return blo->productDomain()->numBlocks(); }
369 inline LinearOp
getBlock(
int i,
int j,
const BlockedLinearOp & blo)
370 {
return blo->getBlock(i,j); }
373 inline void setBlock(
int i,
int j,BlockedLinearOp & blo,
const LinearOp & lo)
374 {
return blo->setBlock(i,j,lo); }
378 {
return rcp(
new Thyra::DefaultBlockedLinearOp<double>()); }
392 { blo->beginBlockFill(rowCnt,colCnt); }
404 { blo->beginBlockFill(); }
408 { blo->endBlockFill(); }
411 BlockedLinearOp getUpperTriBlocks(
const BlockedLinearOp & blo,
bool callEndBlockFill=
true);
414 BlockedLinearOp getLowerTriBlocks(
const BlockedLinearOp & blo,
bool callEndBlockFill=
true);
435 BlockedLinearOp zeroBlockedOp(
const BlockedLinearOp & blo);
438 bool isZeroOp(
const LinearOp op);
448 ModifiableLinearOp getAbsRowSumMatrix(
const LinearOp & op);
458 ModifiableLinearOp getAbsRowSumInvMatrix(
const LinearOp & op);
467 ModifiableLinearOp getLumpedMatrix(
const LinearOp & op);
477 ModifiableLinearOp getInvLumpedMatrix(
const LinearOp & op);
502 void applyOp(
const LinearOp & A,
const MultiVector & x,MultiVector & y,
double alpha=1.0,
double beta=0.0);
522 inline void applyOp(
const LinearOp & A,
const BlockedMultiVector & x,BlockedMultiVector & y,
double alpha=1.0,
double beta=0.0)
524 applyOp(A,x_mv,y_mv,alpha,beta); }
535 void update(
double alpha,
const MultiVector & x,
double beta,MultiVector & y);
538 inline void update(
double alpha,
const BlockedMultiVector & x,
double beta,BlockedMultiVector & y)
540 update(alpha,x_mv,beta,y_mv); }
543 inline void scale(
const double alpha,MultiVector & x) { Thyra::scale<double>(alpha,x.ptr()); }
546 inline void scale(
const double alpha,BlockedMultiVector & x)
550 inline LinearOp
scale(
const double alpha,ModifiableLinearOp & a)
551 {
return Thyra::nonconstScale(alpha,a); }
554 inline LinearOp
adjoint(ModifiableLinearOp & a)
555 {
return Thyra::nonconstAdjoint(a); }
573 const ModifiableLinearOp getDiagonalOp(
const LinearOp & op);
580 const MultiVector getDiagonal(
const LinearOp & op);
593 const ModifiableLinearOp getInvDiagonalOp(
const LinearOp & op);
607 const LinearOp explicitMultiply(
const LinearOp & opl,
const LinearOp & opm,
const LinearOp & opr);
623 const ModifiableLinearOp explicitMultiply(
const LinearOp & opl,
const LinearOp & opm,
const LinearOp & opr,
624 const ModifiableLinearOp & destOp);
636 const LinearOp explicitMultiply(
const LinearOp & opl,
const LinearOp & opr);
651 const ModifiableLinearOp explicitMultiply(
const LinearOp & opl,
const LinearOp & opr,
652 const ModifiableLinearOp & destOp);
664 const LinearOp explicitAdd(
const LinearOp & opl,
const LinearOp & opr);
678 const ModifiableLinearOp explicitAdd(
const LinearOp & opl,
const LinearOp & opr,
679 const ModifiableLinearOp & destOp);
683 const ModifiableLinearOp explicitSum(
const LinearOp & opl,
684 const ModifiableLinearOp & destOp);
689 const LinearOp explicitTranspose(
const LinearOp & op);
694 const LinearOp buildDiagonal(
const MultiVector & v,
const std::string & lbl=
"ANYM");
699 const LinearOp buildInvDiagonal(
const MultiVector & v,
const std::string & lbl=
"ANYM");
726 double computeSpectralRad(
const Teuchos::RCP<
const Thyra::LinearOpBase<double> > & A,
double tol,
727 bool isHermitian=
false,
int numBlocks=5,
int restart=0,
int verbosity=0);
752 double computeSmallestMagEig(
const Teuchos::RCP<
const Thyra::LinearOpBase<double> > & A,
double tol,
753 bool isHermitian=
false,
int numBlocks=5,
int restart=0,
int verbosity=0);
771 ModifiableLinearOp getDiagonalOp(
const Teko::LinearOp & A,
const DiagonalType & dt);
781 ModifiableLinearOp getInvDiagonalOp(
const Teko::LinearOp & A,
const DiagonalType & dt);
788 const MultiVector getDiagonal(
const LinearOp & op,
const DiagonalType & dt);
796 std::string getDiagonalName(
const DiagonalType & dt);
806 DiagonalType getDiagonalType(std::string name);
808 LinearOp probe(Teuchos::RCP<const Epetra_CrsGraph> &G,
const LinearOp & Op);
812 double norm_1(
const MultiVector & v,std::size_t col);
816 double norm_2(
const MultiVector & v,std::size_t col);
821 void clipLower(MultiVector & v,
double lowerBound);
826 void clipUpper(MultiVector & v,
double upperBound);
831 void replaceValue(MultiVector & v,
double currentValue,
double newValue);
835 void columnAverages(
const MultiVector & v,std::vector<double> & averages);
839 double average(
const MultiVector & v);
844 #ifdef TEKO_DEFINED_MSC_EXTENSIONS 845 #undef _MSC_EXTENSIONS LinearOp adjoint(ModifiableLinearOp &a)
Construct an implicit adjoint of the linear operators.
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.
Specifies that a block diagonal approximation is to be used.
Specifies that just the diagonal is used.
MultiVector deepcopy(const MultiVector &v)
Perform a deep copy of the vector.
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
Specifies that row sum is used to form a diagonal.
void scale(const double alpha, MultiVector &x)
Scale a multivector by a constant.
BlockedLinearOp toBlockedLinearOp(LinearOp &clo)
Converts a LinearOp to a BlockedLinearOp.
void endBlockFill(BlockedLinearOp &blo)
Notify the blocked operator that the fill stage is completed.
MultiVector copyAndInit(const MultiVector &v, double scalar)
Perform a deep copy of the vector.
VectorSpace domainSpace(const LinearOp &lo)
Get the domain space of a linear operator.
void beginBlockFill(BlockedLinearOp &blo, int rowCnt, int colCnt)
Let the blocked operator know that you are going to set the sub blocks.
const BlockedMultiVector toBlockedMultiVector(const MultiVector &bmv)
Convert to a BlockedMultiVector from a MultiVector.
VectorSpace rangeSpace(const LinearOp &lo)
Get the range space of a linear operator.
LinearOp toLinearOp(BlockedLinearOp &blo)
Convert to a LinearOp from a BlockedLinearOp.
LinearOp zero(const VectorSpace &vs)
Build a square zero operator from a single vector space.
MultiVector datacopy(const MultiVector &src, MultiVector &dst)
Copy the contents of a multivector to a destination vector.
MultiVector getBlock(int i, const BlockedMultiVector &bmv)
Get the ith block from a BlockedMultiVector object.
For user convenience, if Teko recieves this value, exceptions will be thrown.
DiagonalType
Type describing the type of diagonal to construct.
int blockColCount(const BlockedLinearOp &blo)
Get the column count in a block linear operator.
BlockedLinearOp createBlockedOp()
Build a new blocked linear operator.
void setBlock(int i, int j, BlockedLinearOp &blo, const LinearOp &lo)
Set the i,j block in a BlockedLinearOp object.
Specifies that the diagonal entry is .