00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00035 #ifndef ANASAZI_SVQB_ORTHOMANAGER_HPP
00036 #define ANASAZI_SVQB_ORTHOMANAGER_HPP
00037
00045 #include "AnasaziConfigDefs.hpp"
00046 #include "AnasaziMultiVecTraits.hpp"
00047 #include "AnasaziOperatorTraits.hpp"
00048 #include "AnasaziMatOrthoManager.hpp"
00049 #include "Teuchos_LAPACK.hpp"
00050
00051 namespace Anasazi {
00052
00053 template<class ScalarType, class MV, class OP>
00054 class SVQBOrthoManager : public MatOrthoManager<ScalarType,MV,OP> {
00055
00056 private:
00057 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00058 typedef Teuchos::ScalarTraits<ScalarType> SCT;
00059 typedef Teuchos::ScalarTraits<MagnitudeType> SCTM;
00060 typedef MultiVecTraits<ScalarType,MV> MVT;
00061 typedef OperatorTraits<ScalarType,MV,OP> OPT;
00062 std::string dbgstr;
00063
00064
00065 public:
00066
00067 using MatOrthoManager<ScalarType,MV,OP>::norm;
00068
00070
00071
00072 SVQBOrthoManager( Teuchos::RefCountPtr<const OP> Op = Teuchos::null, bool debug = false );
00073
00075 ~SVQBOrthoManager() {};
00077
00078
00080
00081
00082
00106 void project ( MV &X, Teuchos::RefCountPtr<MV> MX,
00107 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00108 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q) const;
00109
00110
00113 void project ( MV &X,
00114 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00115 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q) const {
00116 project(X,Teuchos::null,C,Q);
00117 }
00118
00119
00147 int normalize ( MV &X, Teuchos::RefCountPtr<MV> MX,
00148 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const;
00149
00150
00153 int normalize ( MV &X, Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const {
00154 return normalize(X,Teuchos::null,B);
00155 }
00156
00157
00193 int projectAndNormalize ( MV &X, Teuchos::RefCountPtr<MV> MX,
00194 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00195 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00196 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q ) const;
00197
00200 int projectAndNormalize ( MV &X,
00201 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00202 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00203 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q ) const {
00204 return projectAndNormalize(X,Teuchos::null,C,B,Q);
00205 }
00206
00208
00210
00211
00215 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00216 orthonormError(const MV &X) const {
00217 return orthonormError(X,Teuchos::null);
00218 }
00219
00224 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00225 orthonormError(const MV &X, Teuchos::RefCountPtr<const MV> MX) const;
00226
00230 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00231 orthogError(const MV &X1, const MV &X2) const {
00232 return orthogError(X1,Teuchos::null,X2);
00233 }
00234
00239 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00240 orthogError(const MV &X1, Teuchos::RefCountPtr<const MV> MX1, const MV &X2) const;
00241
00243
00244 private:
00245
00246 MagnitudeType eps_;
00247 bool debug_;
00248
00249
00250 int findBasis( MV &X, Teuchos::RefCountPtr<MV> MX,
00251 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00252 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00253 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q,
00254 bool normalize ) const;
00255 };
00256
00257
00259
00260 template<class ScalarType, class MV, class OP>
00261 SVQBOrthoManager<ScalarType,MV,OP>::SVQBOrthoManager( Teuchos::RefCountPtr<const OP> Op, bool debug)
00262 : MatOrthoManager<ScalarType,MV,OP>(Op), dbgstr(" *** "), debug_(debug) {
00263
00264 Teuchos::LAPACK<int,MagnitudeType> lapack;
00265 eps_ = lapack.LAMCH('E');
00266 if (debug_) {
00267 std::cout << "eps_ == " << eps_ << endl;
00268 }
00269 }
00270
00271
00273
00274 template<class ScalarType, class MV, class OP>
00275 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00276 SVQBOrthoManager<ScalarType,MV,OP>::orthonormError(const MV &X, Teuchos::RefCountPtr<const MV> MX) const {
00277 const ScalarType ONE = SCT::one();
00278 int rank = MVT::GetNumberVecs(X);
00279 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank);
00280 innerProd(X,X,MX,xTx);
00281 for (int i=0; i<rank; i++) {
00282 xTx(i,i) -= ONE;
00283 }
00284 return xTx.normFrobenius();
00285 }
00286
00288
00289 template<class ScalarType, class MV, class OP>
00290 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00291 SVQBOrthoManager<ScalarType,MV,OP>::orthogError(const MV &X1, Teuchos::RefCountPtr<const MV> MX1, const MV &X2) const {
00292 int r1 = MVT::GetNumberVecs(X1);
00293 int r2 = MVT::GetNumberVecs(X2);
00294 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r2,r1);
00295 innerProd(X2,X1,MX1,xTx);
00296 return xTx.normFrobenius();
00297 }
00298
00299
00301
00302 template<class ScalarType, class MV, class OP>
00303 int SVQBOrthoManager<ScalarType, MV, OP>::projectAndNormalize(
00304 MV &X, Teuchos::RefCountPtr<MV> MX,
00305 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00306 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00307 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q ) const {
00308
00309 return findBasis(X,MX,C,B,Q,true);
00310 }
00311
00312
00313
00315
00316 template<class ScalarType, class MV, class OP>
00317 int SVQBOrthoManager<ScalarType, MV, OP>::normalize(
00318 MV &X, Teuchos::RefCountPtr<MV> MX,
00319 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const {
00320 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C;
00321 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q;
00322 return findBasis(X,MX,C,B,Q,true);
00323 }
00324
00325
00326
00328 template<class ScalarType, class MV, class OP>
00329 void SVQBOrthoManager<ScalarType, MV, OP>::project(
00330 MV &X, Teuchos::RefCountPtr<MV> MX,
00331 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00332 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q) const {
00333 findBasis(X,MX,C,Teuchos::null,Q,false);
00334 }
00335
00336
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363 template<class ScalarType, class MV, class OP>
00364 int SVQBOrthoManager<ScalarType, MV, OP>::findBasis(
00365 MV &X, Teuchos::RefCountPtr<MV> MX,
00366 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00367 Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00368 Teuchos::Array<Teuchos::RefCountPtr<const MV> > Q,
00369 bool normalize) const {
00370
00371 const ScalarType ONE = SCT::one();
00372 const MagnitudeType MONE = SCTM::one();
00373 const MagnitudeType ZERO = SCTM::zero();
00374
00375 int numGS = 0,
00376 numSVQB = 0,
00377 numRand = 0;
00378
00379
00380 int xc = MVT::GetNumberVecs(X);
00381 int xr = MVT::GetVecLength( X );
00382
00383
00384 int nq = Q.length();
00385 int qr = (nq == 0) ? 0 : MVT::GetVecLength(*Q[0]);
00386 int qsize = 0;
00387 std::vector<int> qcs(nq);
00388 for (int i=0; i<nq; i++) {
00389 qcs[i] = MVT::GetNumberVecs(*Q[i]);
00390 qsize += qcs[i];
00391 }
00392
00393 if (normalize == true && qsize + xc > xr) {
00394
00395 TEST_FOR_EXCEPTION( true, std::invalid_argument,
00396 "Anasazi::SVQBOrthoManager::findBasis(): Orthogonalization constraints not feasible" );
00397 }
00398
00399
00400 if (normalize == false && (qsize == 0 || xc == 0)) {
00401
00402 return 0;
00403 }
00404 else if (normalize == true && (xc == 0 || xr == 0)) {
00405
00406 TEST_FOR_EXCEPTION( true, std::invalid_argument,
00407 "Anasazi::SVQBOrthoManager::findBasis(): X must be non-empty" );
00408 }
00409
00410
00411 TEST_FOR_EXCEPTION( qsize != 0 && qr != xr , std::invalid_argument,
00412 "Anasazi::SVQBOrthoManager::findBasis(): Size of X not consistant with size of Q" );
00413
00414
00415
00416
00417
00418 C.resize(nq);
00419 Teuchos::Array<Teuchos::RefCountPtr<Teuchos::SerialDenseMatrix<int,ScalarType> > > newC(nq);
00420
00421 for (int i=0; i<nq; i++) {
00422
00423 TEST_FOR_EXCEPTION( MVT::GetVecLength( *Q[i] ) != qr, std::invalid_argument,
00424 "Anasazi::SVQBOrthoManager::findBasis(): Size of Q not mutually consistant" );
00425 TEST_FOR_EXCEPTION( qr < qcs[i], std::invalid_argument,
00426 "Anasazi::SVQBOrthoManager::findBasis(): Q has less rows than columns" );
00427
00428 if ( C[i] == Teuchos::null ) {
00429 C[i] = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) );
00430 }
00431 else {
00432 TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc, std::invalid_argument,
00433 "Anasazi::SVQBOrthoManager::findBasis(): Size of Q not consistant with C" );
00434 }
00435
00436 C[i]->putScalar(ZERO);
00437 newC[i] = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(*C[i]) );
00438 }
00439
00440
00442
00443
00444
00445
00446 if (normalize == true) {
00447 if ( B == Teuchos::null ) {
00448 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
00449 }
00450 TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
00451 "Anasazi::SVQBOrthoManager::findBasis(): Size of B not consistant with X" );
00452
00453 B->putScalar(ZERO);
00454 for (int i=0; i<xc; i++) {
00455 (*B)(i,i) = MONE;
00456 }
00457 }
00458
00459
00460
00461
00462
00463
00464
00465
00466 Teuchos::RefCountPtr<MV> workX;
00467 if (normalize) {
00468 workX = MVT::Clone(X,xc);
00469 }
00470 if (this->_hasOp) {
00471 if (MX == Teuchos::null) {
00472
00473 MX = MVT::Clone(X,xc);
00474 OPT::Apply(*(this->_Op),X,*MX);
00475 }
00476 }
00477 else {
00478 MX = Teuchos::rcp(&X,false);
00479 }
00480 std::vector<MagnitudeType> normX(xc), invnormX(xc);
00481 Teuchos::SerialDenseMatrix<int,ScalarType> XtMX(xc,xc), workU(1,1);
00482 Teuchos::LAPACK<int,ScalarType> lapack;
00483
00484
00485
00486
00487 std::vector<ScalarType> work;
00488 std::vector<MagnitudeType> lambda, lambdahi, rwork;
00489 if (normalize) {
00490
00491 int lwork = lapack.ILAENV(1,"hetrd","VU",xc,-1,-1,-1);
00492
00493
00494 TEST_FOR_EXCEPTION( lwork < 0, OrthoError,
00495 "Anasazi::SVQBOrthoManager::findBasis(): Error code from ILAENV" );
00496
00497 lwork = (lwork+2)*xc;
00498 work.resize(lwork);
00499
00500 lwork = (3*xc-2 > 1) ? 3*xc - 2 : 1;
00501 rwork.resize(lwork);
00502
00503 lambda.resize(xc);
00504 lambdahi.resize(xc);
00505 workU.reshape(xc,xc);
00506 }
00507
00508
00509 int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc;
00510 int mxr = (this->_hasOp) ? MVT::GetVecLength( *MX ) : xr;
00511 TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument,
00512 "Anasazi::SVQBOrthoManager::findBasis(): Size of X not consistant with MX" );
00513
00514
00515 bool doGramSchmidt = true;
00516
00517 MagnitudeType tolerance = MONE/SCTM::squareroot(eps_);
00518
00519
00520 while (doGramSchmidt) {
00521
00523
00524 if (qsize > 0) {
00525
00526 numGS++;
00527
00528
00529 norm(X,MX,&normX);
00530
00531 Teuchos::RefCountPtr<MV> Xi,MXi;
00532 std::vector<int> ind(1);
00533 for (int i=0; i<xc; i++) {
00534 invnormX[i] = (normX[i] == ZERO) ? ZERO : MONE/normX[i];
00535 ind[0] = i;
00536 Xi = MVT::CloneView(X,ind);
00537 MVT::MvAddMv(ZERO,*Xi,invnormX[i],*Xi,*Xi);
00538 Xi = Teuchos::null;
00539 if (this->_hasOp) {
00540 MXi = MVT::CloneView(*MX,ind);
00541 MVT::MvAddMv(ZERO,*MXi,invnormX[i],*MXi,*MXi);
00542 MXi = Teuchos::null;
00543 }
00544 }
00545
00546 if (debug_) {
00547 std::vector<MagnitudeType> nrm2(xc);
00548 cout << dbgstr << "max post-scale norm: (with/without MX) : ";
00549 MagnitudeType maxpsnw = ZERO, maxpsnwo = ZERO;
00550 norm(X,MX,&nrm2);
00551 for (int i=0; i<xc; i++) {
00552 maxpsnw = (nrm2[i] > maxpsnw ? nrm2[i] : maxpsnw);
00553 }
00554 norm(X,&nrm2);
00555 for (int i=0; i<xc; i++) {
00556 maxpsnwo = (nrm2[i] > maxpsnwo ? nrm2[i] : maxpsnwo);
00557 }
00558 cout << "(" << maxpsnw << "," << maxpsnwo << ")" << endl;
00559 }
00560
00561 for (int i=0; i<nq; i++) {
00562 innerProd(*Q[i],X,MX,*newC[i]);
00563 }
00564
00565 for (int i=0; i<nq; i++) {
00566 MVT::MvTimesMatAddMv(-ONE,*Q[i],*newC[i],ONE,X);
00567 }
00568
00569 for (int i=0; i<xc; i++) {
00570 ind[0] = i;
00571 Xi = MVT::CloneView(X,ind);
00572 MVT::MvAddMv(ZERO,*Xi,normX[i],*Xi,*Xi);
00573 Xi = Teuchos::null;
00574 }
00575
00576 if (this->_hasOp) {
00577 OPT::Apply(*(this->_Op),X,*MX);
00578 }
00579
00580
00581
00582
00583
00584
00585 MagnitudeType maxNorm = ZERO;
00586 for (int j=0; j<xc; j++) {
00587 MagnitudeType sum = ZERO;
00588 for (int k=0; k<nq; k++) {
00589 for (int i=0; i<qcs[k]; i++) {
00590 sum += SCT::magnitude((*newC[k])(i,j))*SCT::magnitude((*newC[k])(i,j));
00591 }
00592 }
00593 maxNorm = (sum > maxNorm) ? sum : maxNorm;
00594 }
00595
00596
00597 if (maxNorm < 0.36) {
00598 doGramSchmidt = false;
00599 }
00600
00601
00602 for (int k=0; k<nq; k++) {
00603 for (int j=0; j<xc; j++) {
00604 for (int i=0; i<qcs[k]; i++) {
00605 (*newC[k])(i,j) *= normX[j];
00606 }
00607 }
00608 }
00609
00610 if (normalize) {
00611
00612 int info;
00613 for (int i=0; i<nq; i++) {
00614 info = C[i]->multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,ONE,*newC[i],*B,ONE);
00615 TEST_FOR_EXCEPTION(info != 0, std::logic_error, "Anasazi::SVQBOrthoManager::findBasis(): Input error to SerialDenseMatrix::multiply.");
00616 }
00617 }
00618 else {
00619
00620 for (int i=0; i<nq; i++) {
00621 (*C[i]) += *newC[i];
00622 }
00623 }
00624 }
00625 else {
00626
00627 doGramSchmidt = false;
00628 }
00629
00630
00632
00633 if (normalize) {
00634
00635 MagnitudeType condT = tolerance;
00636
00637 while (condT >= tolerance) {
00638
00639 numSVQB++;
00640
00641
00642 innerProd(X,X,MX,XtMX);
00643
00644
00645 std::vector<MagnitudeType> Dh(xc), Dhi(xc);
00646 for (int i=0; i<xc; i++) {
00647 Dh[i] = SCT::magnitude(SCT::squareroot(XtMX(i,i)));
00648 Dhi[i] = (Dh[i] == ZERO ? ZERO : MONE/Dh[i]);
00649 }
00650
00651 for (int i=0; i<xc; i++) {
00652 for (int j=0; j<xc; j++) {
00653 XtMX(i,j) *= Dhi[i]*Dhi[j];
00654 }
00655 }
00656
00657
00658 int info;
00659 lapack.HEEV('V', 'U', xc, XtMX.values(), XtMX.stride(), &lambda[0], &work[0], work.size(), &rwork[0], &info);
00660 TEST_FOR_EXCEPTION( info != 0, OrthoError,
00661 "Anasazi::SVQBOrthoManager::findBasis(): Error code from HEEV" );
00662 if (debug_) {
00663 cout << dbgstr << "eigenvalues of XtMX: (";
00664 for (int i=0; i<xc-1; i++) {
00665 cout << lambda[i] << ",";
00666 }
00667 cout << lambda[xc-1] << ")" << endl;
00668 }
00669
00670
00671
00672 MagnitudeType maxLambda = lambda[xc-1],
00673 minLambda = lambda[0];
00674 int iZeroMax = -1;
00675 for (int i=0; i<xc; i++) {
00676 if (lambda[i] < 10*eps_*maxLambda) {
00677 iZeroMax = i;
00678 lambda[i] = ZERO;
00679 lambdahi[i] = ZERO;
00680 }
00681
00682
00683
00684
00685
00686
00687 else {
00688 lambda[i] = SCTM::squareroot(lambda[i]);
00689 lambdahi[i] = MONE/lambda[i];
00690 }
00691 }
00692
00693
00694
00695
00696 std::vector<int> ind(xc);
00697 for (int i=0; i<xc; i++) {ind[i] = i;}
00698 MVT::SetBlock(X,ind,*workX);
00699
00700
00701 workU.assign(XtMX);
00702 for (int j=0; j<xc; j++) {
00703 for (int i=0; i<xc; i++) {
00704 workU(i,j) *= Dhi[i]*lambdahi[j];
00705 }
00706 }
00707
00708 MVT::MvTimesMatAddMv(ONE,*workX,workU,ZERO,X);
00709
00710
00711
00712
00713 if (this->_hasOp) {
00714 if (maxLambda >= tolerance * minLambda) {
00715
00716 OPT::Apply(*(this->_Op),X,*MX);
00717 }
00718 else {
00719
00720
00721 MVT::SetBlock(*MX,ind,*workX);
00722
00723
00724 MVT::MvTimesMatAddMv(ONE,*workX,workU,ZERO,*MX);
00725 }
00726 }
00727
00728
00729
00730 for (int j=0; j<xc; j++) {
00731 for (int i=0; i<xc; i++) {
00732 workU(i,j) = Dh[i] * (*B)(i,j);
00733 }
00734 }
00735 info = B->multiply(Teuchos::CONJ_TRANS,Teuchos::NO_TRANS,ONE,XtMX,workU,ZERO);
00736 TEST_FOR_EXCEPTION(info != 0, std::logic_error, "Anasazi::SVQBOrthoManager::findBasis(): Input error to SerialDenseMatrix::multiply.");
00737 for (int j=0; j<xc ;j++) {
00738 for (int i=0; i<xc; i++) {
00739 (*B)(i,j) *= lambda[i];
00740 }
00741 }
00742
00743
00744 if (iZeroMax >= 0) {
00745 if (debug_) {
00746 cout << dbgstr << "augmenting multivec with " << iZeroMax+1 << " random directions" << endl;
00747 }
00748
00749 numRand++;
00750
00751 std::vector<int> ind(iZeroMax+1);
00752 for (int i=0; i<iZeroMax+1; i++) {
00753 ind[i] = i;
00754 }
00755 Teuchos::RefCountPtr<MV> Xnull,MXnull;
00756 Xnull = MVT::CloneView(X,ind);
00757 MVT::MvRandom(*Xnull);
00758 if (this->_hasOp) {
00759 MXnull = MVT::CloneView(*MX,ind);
00760 OPT::Apply(*(this->_Op),*Xnull,*MXnull);
00761 MXnull = Teuchos::null;
00762 }
00763 Xnull = Teuchos::null;
00764 condT = tolerance;
00765 doGramSchmidt = true;
00766 break;
00767 }
00768
00769 condT = SCTM::magnitude(maxLambda / minLambda);
00770 if (debug_) {
00771 cout << dbgstr << "condT: " << condT << endl;
00772 }
00773
00774 }
00775
00776 if ((doGramSchmidt == false) && (condT > SCTM::squareroot(tolerance))) {
00777 doGramSchmidt = true;
00778 }
00779 }
00780
00781
00782 }
00783
00784 if (debug_) {
00785 cout << dbgstr << "(numGS,numSVQB,numRand) : "
00786 << "(" << numGS
00787 << "," << numSVQB
00788 << "," << numRand
00789 << ")" << endl;
00790 }
00791
00792 return xc;
00793 }
00794
00795 }
00796
00797 #endif // ANASAZI_SVQB_ORTHOMANAGER_HPP
00798