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