15#ifndef ANASAZI_BASIC_ORTHOMANAGER_HPP
16#define ANASAZI_BASIC_ORTHOMANAGER_HPP
29#include "Teuchos_TimeMonitor.hpp"
30#include "Teuchos_LAPACK.hpp"
31#include "Teuchos_BLAS.hpp"
32#ifdef ANASAZI_BASIC_ORTHO_DEBUG
33# include <Teuchos_FancyOStream.hpp>
38 template<
class ScalarType,
class MV,
class OP>
42 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
43 typedef Teuchos::ScalarTraits<ScalarType> SCT;
53 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa = 1.41421356 ,
54 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType eps = 0.0,
55 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType tol = 0.20 );
107 Teuchos::Array<Teuchos::RCP<const MV> > Q,
108 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C
109 = Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix<int,ScalarType> >(Teuchos::null)),
110 Teuchos::RCP<MV> MX = Teuchos::null,
111 Teuchos::Array<Teuchos::RCP<const MV> > MQ = Teuchos::tuple(Teuchos::RCP<const MV>(Teuchos::null))
155 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B = Teuchos::null,
156 Teuchos::RCP<MV> MX = Teuchos::null
221 Teuchos::Array<Teuchos::RCP<const MV> > Q,
222 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C
223 = Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix<int,ScalarType> >(Teuchos::null)),
224 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B = Teuchos::null,
225 Teuchos::RCP<MV> MX = Teuchos::null,
226 Teuchos::Array<Teuchos::RCP<const MV> > MQ = Teuchos::tuple(Teuchos::RCP<const MV>(Teuchos::null))
238 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
245 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
246 orthogErrorMat(
const MV &X1,
const MV &X2, Teuchos::RCP<const MV> MX1, Teuchos::RCP<const MV> MX2)
const;
254 void setKappa(
typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa ) { kappa_ = kappa; }
257 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
getKappa()
const {
return kappa_; }
263 MagnitudeType kappa_;
268 int findBasis(MV &X, Teuchos::RCP<MV> MX,
269 Teuchos::SerialDenseMatrix<int,ScalarType> &B,
270 bool completeBasis,
int howMany = -1 )
const;
275 Teuchos::RCP<Teuchos::Time> timerReortho_;
282 template<
class ScalarType,
class MV,
class OP>
284 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa,
285 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType eps,
286 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType tol ) :
287 MatOrthoManager<ScalarType,MV,OP>(Op), kappa_(kappa), eps_(eps), tol_(tol)
288#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
289 , timerReortho_(Teuchos::TimeMonitor::getNewTimer(
"Anasazi::BasicOrthoManager::Re-orthogonalization"))
292 TEUCHOS_TEST_FOR_EXCEPTION(eps_ < SCT::magnitude(SCT::zero()),std::invalid_argument,
293 "Anasazi::BasicOrthoManager::BasicOrthoManager(): argument \"eps\" must be non-negative.");
295 Teuchos::LAPACK<int,MagnitudeType> lapack;
296 eps_ = lapack.LAMCH(
'E');
297 eps_ = Teuchos::ScalarTraits<MagnitudeType>::pow(eps_,.75);
299 TEUCHOS_TEST_FOR_EXCEPTION(
300 tol_ < SCT::magnitude(SCT::zero()) || tol_ > SCT::magnitude(SCT::one()),
301 std::invalid_argument,
302 "Anasazi::BasicOrthoManager::BasicOrthoManager(): argument \"tol\" must be in [0,1].");
309 template<
class ScalarType,
class MV,
class OP>
310 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
312 const ScalarType ONE = SCT::one();
314 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank);
316 for (
int i=0; i<rank; i++) {
319 return xTx.normFrobenius();
326 template<
class ScalarType,
class MV,
class OP>
327 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
331 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r1,r2);
333 return xTx.normFrobenius();
339 template<
class ScalarType,
class MV,
class OP>
342 Teuchos::Array<Teuchos::RCP<const MV> > Q,
343 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
345 Teuchos::Array<Teuchos::RCP<const MV> > MQ
362#ifdef ANASAZI_BASIC_ORTHO_DEBUG
364 Teuchos::RCP<Teuchos::FancyOStream>
365 out = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
366 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
367 *out <<
"Entering Anasazi::BasicOrthoManager::projectMat(...)\n";
370 ScalarType ONE = SCT::one();
375 std::vector<int> qcs(nq);
377 if (nq == 0 || xc == 0 || xr == 0) {
378#ifdef ANASAZI_BASIC_ORTHO_DEBUG
379 *out <<
"Leaving Anasazi::BasicOrthoManager::projectMat(...)\n";
392#ifdef ANASAZI_BASIC_ORTHO_DEBUG
393 *out <<
"Allocating MX...\n";
395 if (MX == Teuchos::null) {
404 MX = Teuchos::rcpFromRef(X);
410 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
411 "Anasazi::BasicOrthoManager::projectMat(): MVT returned negative dimensions for X,MX" );
413 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr || xr!=qr, std::invalid_argument,
414 "Anasazi::BasicOrthoManager::projectMat(): Size of X not consistent with MX,Q" );
417 for (
int i=0; i<nq; i++) {
419 "Anasazi::BasicOrthoManager::projectMat(): Q lengths not mutually consistent" );
421 TEUCHOS_TEST_FOR_EXCEPTION( qr <
static_cast<ptrdiff_t
>(qcs[i]), std::invalid_argument,
422 "Anasazi::BasicOrthoManager::projectMat(): Q has less rows than columns" );
424 if ( C[i] == Teuchos::null ) {
425 C[i] = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) );
428 TEUCHOS_TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc , std::invalid_argument,
429 "Anasazi::BasicOrthoManager::projectMat(): Size of Q not consistent with size of C" );
436 std::vector<ScalarType> oldDot( xc );
438#ifdef ANASAZI_BASIC_ORTHO_DEBUG
439 *out <<
"oldDot = { ";
440 std::copy(oldDot.begin(), oldDot.end(), std::ostream_iterator<ScalarType>(*out,
" "));
446 for (
int i=0; i<nq; i++) {
450#ifdef ANASAZI_BASIC_ORTHO_DEBUG
451 *out <<
"Applying projector P_Q[" << i <<
"]...\n";
458 if (MQ[i] == Teuchos::null) {
459#ifdef ANASAZI_BASIC_ORTHO_DEBUG
460 *out <<
"Updating MX via M*X...\n";
466#ifdef ANASAZI_BASIC_ORTHO_DEBUG
467 *out <<
"Updating MX via M*Q...\n";
475 std::vector<ScalarType> newDot(xc);
477#ifdef ANASAZI_BASIC_ORTHO_DEBUG
478 *out <<
"newDot = { ";
479 std::copy(newDot.begin(), newDot.end(), std::ostream_iterator<ScalarType>(*out,
" "));
484 for (
int j = 0; j < xc; ++j) {
486 if ( SCT::magnitude(kappa_*newDot[j]) < SCT::magnitude(oldDot[j]) ) {
487#ifdef ANASAZI_BASIC_ORTHO_DEBUG
488 *out <<
"kappa_*newDot[" <<j<<
"] == " << kappa_*newDot[j] <<
"... another step of Gram-Schmidt.\n";
490#ifdef ANASAZI_TEUCHOS_TIME_MONITOR
491 Teuchos::TimeMonitor lcltimer( *timerReortho_ );
493 for (
int i=0; i<nq; i++) {
494 Teuchos::SerialDenseMatrix<int,ScalarType> C2(C[i]->numRows(), C[i]->numCols());
499#ifdef ANASAZI_BASIC_ORTHO_DEBUG
500 *out <<
"Applying projector P_Q[" << i <<
"]...\n";
506 if (MQ[i] == Teuchos::null) {
507#ifdef ANASAZI_BASIC_ORTHO_DEBUG
508 *out <<
"Updating MX via M*X...\n";
514#ifdef ANASAZI_BASIC_ORTHO_DEBUG
515 *out <<
"Updating MX via M*Q...\n";
525#ifdef ANASAZI_BASIC_ORTHO_DEBUG
526 *out <<
"Leaving Anasazi::BasicOrthoManager::projectMat(...)\n";
534 template<
class ScalarType,
class MV,
class OP>
537 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
538 Teuchos::RCP<MV> MX)
const {
548 if (MX == Teuchos::null) {
558 if ( B == Teuchos::null ) {
559 B = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
566 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
567 "Anasazi::BasicOrthoManager::normalizeMat(): X must be non-empty" );
568 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
569 "Anasazi::BasicOrthoManager::normalizeMat(): Size of X not consistent with size of B" );
570 TEUCHOS_TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument,
571 "Anasazi::BasicOrthoManager::normalizeMat(): Size of X not consistent with size of MX" );
572 TEUCHOS_TEST_FOR_EXCEPTION(
static_cast<ptrdiff_t
>(xc) > xr, std::invalid_argument,
573 "Anasazi::BasicOrthoManager::normalizeMat(): Size of X not feasible for normalization" );
575 return findBasis(X, MX, *B,
true );
582 template<
class ScalarType,
class MV,
class OP>
585 Teuchos::Array<Teuchos::RCP<const MV> > Q,
586 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
587 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
589 Teuchos::Array<Teuchos::RCP<const MV> > MQ
592#ifdef ANASAZI_BASIC_ORTHO_DEBUG
594 Teuchos::RCP<Teuchos::FancyOStream>
595 out = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
596 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
597 *out <<
"Entering Anasazi::BasicOrthoManager::projectAndNormalizeMat(...)\n";
608 if ( B == Teuchos::null ) {
609 B = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
614 if (MX == Teuchos::null) {
616#ifdef ANASAZI_BASIC_ORTHO_DEBUG
617 *out <<
"Allocating MX...\n";
626 MX = Teuchos::rcpFromRef(X);
632 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
"Anasazi::BasicOrthoManager::projectAndNormalizeMat(): X must be non-empty" );
634 ptrdiff_t numbas = 0;
635 for (
int i=0; i<nq; i++) {
640 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
641 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Size of X must be consistent with size of B" );
643 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
644 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): MVT returned negative dimensions for X,MX" );
646 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument,
647 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Size of X must be consistent with size of MX" );
649 TEUCHOS_TEST_FOR_EXCEPTION( numbas+xc > xr, std::invalid_argument,
650 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Orthogonality constraints not feasible" );
653#ifdef ANASAZI_BASIC_ORTHO_DEBUG
654 *out <<
"Orthogonalizing X against Q...\n";
658 Teuchos::SerialDenseMatrix<int,ScalarType> oldCoeff(xc,1);
664 int curxsize = xc - rank;
669#ifdef ANASAZI_BASIC_ORTHO_DEBUG
670 *out <<
"Attempting to find orthonormal basis for X...\n";
672 rank = findBasis(X,MX,*B,
false,curxsize);
674 if (oldrank != -1 && rank != oldrank) {
680 for (
int i=0; i<xc; i++) {
681 (*B)(i,oldrank) = oldCoeff(i,0);
686 if (rank != oldrank) {
694 for (
int i=0; i<xc; i++) {
695 oldCoeff(i,0) = (*B)(i,rank);
702#ifdef ANASAZI_BASIC_ORTHO_DEBUG
703 *out <<
"Finished computing basis.\n";
708 TEUCHOS_TEST_FOR_EXCEPTION( rank < oldrank,
OrthoError,
709 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): basis lost rank; this shouldn't happen");
711 if (rank != oldrank) {
727#ifdef ANASAZI_BASIC_ORTHO_DEBUG
728 *out <<
"Randomizing X[" << rank <<
"]...\n";
730 Teuchos::RCP<MV> curX, curMX;
731 std::vector<int> ind(1);
736#ifdef ANASAZI_BASIC_ORTHO_DEBUG
737 *out <<
"Applying operator to random vector.\n";
749 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > dummyC(0);
756 TEUCHOS_TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error,
757 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Debug error in rank variable." );
759#ifdef ANASAZI_BASIC_ORTHO_DEBUG
760 *out <<
"Leaving Anasazi::BasicOrthoManager::projectAndNormalizeMat(...)\n";
771 template<
class ScalarType,
class MV,
class OP>
772 int BasicOrthoManager<ScalarType, MV, OP>::findBasis(
773 MV &X, Teuchos::RCP<MV> MX,
774 Teuchos::SerialDenseMatrix<int,ScalarType> &B,
775 bool completeBasis,
int howMany )
const {
792#ifdef ANASAZI_BASIC_ORTHO_DEBUG
794 Teuchos::RCP<Teuchos::FancyOStream>
795 out = Teuchos::getFancyOStream(Teuchos::rcpFromRef(std::cout));
796 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
797 *out <<
"Entering Anasazi::BasicOrthoManager::findBasis(...)\n";
800 const ScalarType ONE = SCT::one();
801 const MagnitudeType ZERO = SCT::magnitude(SCT::zero());
803 int xc = MVT::GetNumberVecs( X );
812 TEUCHOS_TEST_FOR_EXCEPTION(this->_hasOp ==
true && MX == Teuchos::null, std::logic_error,
813 "Anasazi::BasicOrthoManager::findBasis(): calling routine did not specify MS.");
814 TEUCHOS_TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::logic_error,
815 "Anasazi::BasicOrthoManager::findBasis(): Invalid howMany parameter" );
820 int xstart = xc - howMany;
822 for (
int j = xstart; j < xc; j++) {
831 for (
int i=j+1; i<xc; ++i) {
836 std::vector<int> index(1);
838 Teuchos::RCP<MV> Xj = MVT::CloneViewNonConst( X, index );
839 Teuchos::RCP<MV> MXj;
840 if ((this->_hasOp)) {
842 MXj = MVT::CloneViewNonConst( *MX, index );
850 std::vector<int> prev_idx( numX );
851 Teuchos::RCP<const MV> prevX, prevMX;
854 for (
int i=0; i<numX; ++i) prev_idx[i] = i;
855 prevX = MVT::CloneViewNonConst( X, prev_idx );
857 prevMX = MVT::CloneViewNonConst( *MX, prev_idx );
866 for (
int numTrials = 0; numTrials < 10; numTrials++) {
867#ifdef ANASAZI_BASIC_ORTHO_DEBUG
868 *out <<
"Trial " << numTrials <<
" for vector " << j <<
"\n";
872 Teuchos::SerialDenseMatrix<int,ScalarType> product(numX, 1);
873 std::vector<MagnitudeType> origNorm(1), newNorm(1), newNorm2(1);
878 Teuchos::RCP<MV> oldMXj = MVT::CloneCopy( *MXj );
880#ifdef ANASAZI_BASIC_ORTHO_DEBUG
881 *out <<
"origNorm = " << origNorm[0] <<
"\n";
892#ifdef ANASAZI_BASIC_ORTHO_DEBUG
893 *out <<
"Orthogonalizing X[" << j <<
"]...\n";
895 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *Xj );
902#ifdef ANASAZI_BASIC_ORTHO_DEBUG
903 *out <<
"Updating MX[" << j <<
"]...\n";
905 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *MXj );
910 MagnitudeType product_norm = product.normOne();
912#ifdef ANASAZI_BASIC_ORTHO_DEBUG
913 *out <<
"newNorm = " << newNorm[0] <<
"\n";
914 *out <<
"prodoct_norm = " << product_norm <<
"\n";
918 if ( product_norm/newNorm[0] >= tol_ || newNorm[0] < eps_*origNorm[0]) {
919#ifdef ANASAZI_BASIC_ORTHO_DEBUG
920 if (product_norm/newNorm[0] >= tol_) {
921 *out <<
"product_norm/newNorm == " << product_norm/newNorm[0] <<
"... another step of Gram-Schmidt.\n";
924 *out <<
"eps*origNorm == " << eps_*origNorm[0] <<
"... another step of Gram-Schmidt.\n";
929 Teuchos::SerialDenseMatrix<int,ScalarType> P2(numX,1);
932#ifdef ANASAZI_BASIC_ORTHO_DEBUG
933 *out <<
"Orthogonalizing X[" << j <<
"]...\n";
935 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj );
936 if ((this->_hasOp)) {
937#ifdef ANASAZI_BASIC_ORTHO_DEBUG
938 *out <<
"Updating MX[" << j <<
"]...\n";
940 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj );
944 product_norm = P2.normOne();
945#ifdef ANASAZI_BASIC_ORTHO_DEBUG
946 *out <<
"newNorm2 = " << newNorm2[0] <<
"\n";
947 *out <<
"product_norm = " << product_norm <<
"\n";
949 if ( product_norm/newNorm2[0] >= tol_ || newNorm2[0] < eps_*origNorm[0] ) {
951#ifdef ANASAZI_BASIC_ORTHO_DEBUG
952 if (product_norm/newNorm2[0] >= tol_) {
953 *out <<
"product_norm/newNorm2 == " << product_norm/newNorm2[0] <<
"... setting vector to zero.\n";
955 else if (newNorm[0] < newNorm2[0]) {
956 *out <<
"newNorm2 > newNorm... setting vector to zero.\n";
959 *out <<
"eps*origNorm == " << eps_*origNorm[0] <<
"... setting vector to zero.\n";
962 MVT::MvInit(*Xj,ZERO);
963 if ((this->_hasOp)) {
964#ifdef ANASAZI_BASIC_ORTHO_DEBUG
965 *out <<
"Setting MX[" << j <<
"] to zero as well.\n";
967 MVT::MvInit(*MXj,ZERO);
974 if (numTrials == 0) {
975 for (
int i=0; i<numX; i++) {
976 B(i,j) = product(i,0);
982 if ( newNorm[0] != ZERO && newNorm[0] > SCT::sfmin() ) {
983#ifdef ANASAZI_BASIC_ORTHO_DEBUG
984 *out <<
"Normalizing X[" << j <<
"], norm(X[" << j <<
"]) = " << newNorm[0] <<
"\n";
988 MVT::MvScale( *Xj, ONE/newNorm[0]);
990#ifdef ANASAZI_BASIC_ORTHO_DEBUG
991 *out <<
"Normalizing M*X[" << j <<
"]...\n";
994 MVT::MvScale( *MXj, ONE/newNorm[0]);
998 if (numTrials == 0) {
1007#ifdef ANASAZI_BASIC_ORTHO_DEBUG
1008 *out <<
"Not normalizing M*X[" << j <<
"]...\n";
1015 if (completeBasis) {
1017#ifdef ANASAZI_BASIC_ORTHO_DEBUG
1018 *out <<
"Inserting random vector in X[" << j <<
"]...\n";
1020 MVT::MvRandom( *Xj );
1022#ifdef ANASAZI_BASIC_ORTHO_DEBUG
1023 *out <<
"Updating M*X[" << j <<
"]...\n";
1025 OPT::Apply( *(this->_Op), *Xj, *MXj );
1026 this->_OpCounter += MVT::GetNumberVecs(*Xj);
1037 if (rankDef ==
true) {
1038 TEUCHOS_TEST_FOR_EXCEPTION( rankDef && completeBasis,
OrthoError,
1039 "Anasazi::BasicOrthoManager::findBasis(): Unable to complete basis" );
1040#ifdef ANASAZI_BASIC_ORTHO_DEBUG
1041 *out <<
"Returning early, rank " << j <<
" from Anasazi::BasicOrthoManager::findBasis(...)\n";
1048#ifdef ANASAZI_BASIC_ORTHO_DEBUG
1049 *out <<
"Returning " << xc <<
" from Anasazi::BasicOrthoManager::findBasis(...)\n";
Anasazi header file which uses auto-configuration information to include necessary C++ headers.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Declaration of basic traits for the multivector type.
Virtual base class which defines basic traits for the operator type.
int normalizeMat(MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null) const
This method takes a multivector X and attempts to compute an orthonormal basis for ,...
Teuchos::ScalarTraits< ScalarType >::magnitudeType getKappa() const
Return parameter for re-orthogonalization threshold.
~BasicOrthoManager()
Destructor.
int projectAndNormalizeMat(MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const
Given a set of bases Q[i] and a multivector X, this method computes an orthonormal basis for .
void setKappa(typename Teuchos::ScalarTraits< ScalarType >::magnitudeType kappa)
Set parameter for re-orthogonalization threshold.
void projectMat(MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const
Given a list of mutually orthogonal and internally orthonormal bases Q, this method projects a multiv...
Teuchos::ScalarTraits< ScalarType >::magnitudeType orthonormErrorMat(const MV &X, Teuchos::RCP< const MV > MX=Teuchos::null) const
This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of ...
Teuchos::ScalarTraits< ScalarType >::magnitudeType orthogErrorMat(const MV &X1, const MV &X2, Teuchos::RCP< const MV > MX1, Teuchos::RCP< const MV > MX2) const
This method computes the error in orthogonality of two multivectors, measured as the Frobenius norm o...
BasicOrthoManager(Teuchos::RCP< const OP > Op=Teuchos::null, typename Teuchos::ScalarTraits< ScalarType >::magnitudeType kappa=1.41421356, typename Teuchos::ScalarTraits< ScalarType >::magnitudeType eps=0.0, typename Teuchos::ScalarTraits< ScalarType >::magnitudeType tol=0.20)
Constructor specifying re-orthogonalization tolerance.
MatOrthoManager(Teuchos::RCP< const OP > Op=Teuchos::null)
Default constructor.
void normMat(const MV &X, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, Teuchos::RCP< const MV > MX=Teuchos::null) const
Provides the norm induced by the matrix-based inner product.
void innerProdMat(const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z, Teuchos::RCP< const MV > MX=Teuchos::null, Teuchos::RCP< const MV > MY=Teuchos::null) const
Provides a matrix-based inner product.
Traits class which defines basic operations on multivectors.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void MvDot(const MV &mv, const MV &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
Virtual base class which defines basic traits for the operator type.
static void Apply(const OP &Op, const MV &x, MV &y)
Application method which performs operation y = Op*x. An OperatorError exception is thrown if there i...
Exception thrown to signal error in an orthogonalization manager method.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package.