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
00034 #ifndef ANASAZI_BASIC_ORTHOMANAGER_HPP
00035 #define ANASAZI_BASIC_ORTHOMANAGER_HPP
00036
00044
00045
00046 #include "AnasaziConfigDefs.hpp"
00047 #include "AnasaziMultiVecTraits.hpp"
00048 #include "AnasaziOperatorTraits.hpp"
00049 #include "AnasaziMatOrthoManager.hpp"
00050 #include "Teuchos_TimeMonitor.hpp"
00051
00052 namespace Anasazi {
00053
00054 template<class ScalarType, class MV, class OP>
00055 class BasicOrthoManager : public MatOrthoManager<ScalarType,MV,OP> {
00056
00057 private:
00058 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
00059 typedef Teuchos::ScalarTraits<ScalarType> SCT;
00060 typedef MultiVecTraits<ScalarType,MV> MVT;
00061 typedef OperatorTraits<ScalarType,MV,OP> OPT;
00062
00063 public:
00064
00066
00067
00068 BasicOrthoManager( Teuchos::RCP<const OP> Op = Teuchos::null, typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa = 1.5625 );
00069
00070
00072 ~BasicOrthoManager() {}
00074
00075
00077
00078
00080 void setKappa( typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa ) { kappa_ = kappa; }
00081
00083 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType getKappa() const { return kappa_; }
00084
00086
00087
00089
00090
00091
00130 void projectMat (
00131 MV &X,
00132 Teuchos::RCP<MV> MX = Teuchos::null,
00133 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C = Teuchos::tuple(Teuchos::null),
00134 Teuchos::Array<Teuchos::RCP<const MV> > Q = Teuchos::tuple(Teuchos::null) ) const;
00135
00136
00175 int normalizeMat (
00176 MV &X,
00177 Teuchos::RCP<MV> MX = Teuchos::null,
00178 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B = Teuchos::tuple(Teuchos::null) ) const;
00179
00180
00240 int projectAndNormalizeMat (
00241 MV &X,
00242 Teuchos::RCP<MV> MX = Teuchos::null,
00243 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C = Teuchos::tuple(Teuchos::null),
00244 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B = Teuchos::null,
00245 Teuchos::Array<Teuchos::RCP<const MV> > Q = Teuchos::tuple(Teuchos::null) ) const;
00246
00248
00250
00251
00256 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00257 orthonormErrorMat(const MV &X, Teuchos::RCP<const MV> MX = Teuchos::null) const;
00258
00263 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00264 orthogErrorMat(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const;
00265
00267
00268 private:
00269
00271 MagnitudeType kappa_;
00272
00273
00274 int findBasis(MV &X, Teuchos::RCP<MV> MX,
00275 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > C,
00276 bool completeBasis, int howMany = -1 ) const;
00277
00278
00279
00280
00281 Teuchos::RCP<Teuchos::Time> timerReortho_;
00282
00283 };
00284
00285
00287
00288 template<class ScalarType, class MV, class OP>
00289 BasicOrthoManager<ScalarType,MV,OP>::BasicOrthoManager( Teuchos::RCP<const OP> Op,
00290 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType kappa ) :
00291 MatOrthoManager<ScalarType,MV,OP>(Op),
00292 kappa_(kappa),
00293 timerReortho_(Teuchos::TimeMonitor::getNewTimer("BasicOrthoManager::Re-orthogonalization"))
00294 {}
00295
00296
00298
00299 template<class ScalarType, class MV, class OP>
00300 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00301 BasicOrthoManager<ScalarType,MV,OP>::orthonormErrorMat(const MV &X, Teuchos::RCP<const MV> MX) const {
00302 const ScalarType ONE = SCT::one();
00303 int rank = MVT::GetNumberVecs(X);
00304 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank);
00305 innerProdMat(X,X,MX,xTx);
00306 for (int i=0; i<rank; i++) {
00307 xTx(i,i) -= ONE;
00308 }
00309 return xTx.normFrobenius();
00310 }
00311
00313
00314 template<class ScalarType, class MV, class OP>
00315 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
00316 BasicOrthoManager<ScalarType,MV,OP>::orthogErrorMat(const MV &X1, Teuchos::RCP<const MV> MX1, const MV &X2) const {
00317 int r1 = MVT::GetNumberVecs(X1);
00318 int r2 = MVT::GetNumberVecs(X2);
00319 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r2,r1);
00320 innerProdMat(X2,X1,MX1,xTx);
00321 return xTx.normFrobenius();
00322 }
00323
00325
00326 template<class ScalarType, class MV, class OP>
00327 int BasicOrthoManager<ScalarType, MV, OP>::projectAndNormalizeMat(
00328 MV &X, Teuchos::RCP<MV> MX,
00329 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00330 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00331 Teuchos::Array<Teuchos::RCP<const MV> > Q ) const {
00332
00333 int nq = Q.length();
00334 int xc = MVT::GetNumberVecs( X );
00335 int xr = MVT::GetVecLength( X );
00336 int rank;
00337
00338
00339
00340
00341 if ( B == Teuchos::null ) {
00342 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
00343 }
00344
00345
00346 if (this->_hasOp) {
00347 if (MX == Teuchos::null) {
00348
00349 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
00350 OPT::Apply(*(this->_Op),X,*MX);
00351 this->_OpCounter += MVT::GetNumberVecs(X);
00352 }
00353 }
00354 else {
00355
00356 MX = Teuchos::rcp( &X, false );
00357 }
00358
00359 int mxc = MVT::GetNumberVecs( *MX );
00360 int mxr = MVT::GetVecLength( *MX );
00361
00362
00363 TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument, "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): X must be non-empty" );
00364
00365 int numbas = 0;
00366 for (int i=0; i<nq; i++) {
00367 numbas += MVT::GetNumberVecs( *Q[i] );
00368 }
00369
00370
00371 TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
00372 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Size of X must be consistant with size of B" );
00373
00374 TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
00375 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): MVT returned negative dimensions for X,MX" );
00376
00377 TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument,
00378 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Size of X must be consistant with size of MX" );
00379
00380 TEST_FOR_EXCEPTION( numbas+xc > xr, std::invalid_argument,
00381 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Orthogonality constraints not feasible" );
00382
00383
00384 projectMat(X,MX,C,Q);
00385
00386
00387 Teuchos::SerialDenseMatrix<int,ScalarType> oldCoeff(xc,1);
00388
00389
00390 rank = 0;
00391 int numTries = 10;
00392 int oldrank = -1;
00393 do {
00394 int curxsize = xc - rank;
00395
00396
00397
00398
00399 rank = findBasis(X,MX,B,false,curxsize);
00400
00401 if (rank < xc && numTries == 10) {
00402
00403
00404 for (int i=0; i<xc; i++) {
00405 oldCoeff(i,0) = (*B)(i,rank);
00406 }
00407 }
00408
00409 if (oldrank != -1 && rank != oldrank) {
00410
00411 for (int i=0; i<xc; i++) {
00412 (*B)(i,oldrank) = oldCoeff(i,0);
00413 }
00414 }
00415
00416 if (rank == xc) {
00417
00418 break;
00419 }
00420 else {
00421 TEST_FOR_EXCEPTION( rank < oldrank, OrthoError,
00422 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): basis lost rank; this shouldn't happen");
00423
00424 if (rank != oldrank) {
00425
00426 numTries = 10;
00427 }
00428
00429
00430 oldrank = rank;
00431
00432
00433 if (numTries <= 0) {
00434 break;
00435 }
00436
00437 numTries--;
00438
00439
00440 #ifdef ANASAZI_BASICORTHO_DEBUG
00441 cout << "Random for column " << rank << endl;
00442 #endif
00443 Teuchos::RCP<MV> curX, curMX;
00444 std::vector<int> ind(1);
00445 ind[0] = rank;
00446 curX = MVT::CloneView(X,ind);
00447 MVT::MvRandom(*curX);
00448 if (this->_hasOp) {
00449 curMX = MVT::CloneView(*MX,ind);
00450 OPT::Apply( *(this->_Op), *curX, *curMX );
00451 this->_OpCounter += MVT::GetNumberVecs(*curX);
00452 }
00453
00454
00455
00456
00457 projectMat(*curX,curMX,Teuchos::null,Q);
00458 }
00459 } while (1);
00460
00461
00462 TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error,
00463 "Anasazi::BasicOrthoManager::projectAndNormalizeMat(): Debug error in rank variable." );
00464 return rank;
00465 }
00466
00467
00468
00470
00471 template<class ScalarType, class MV, class OP>
00472 int BasicOrthoManager<ScalarType, MV, OP>::normalizeMat(
00473 MV &X, Teuchos::RCP<MV> MX,
00474 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B ) const {
00475
00476 return findBasis(X, MX, B, true );
00477 }
00478
00479
00480
00482 template<class ScalarType, class MV, class OP>
00483 void BasicOrthoManager<ScalarType, MV, OP>::projectMat(
00484 MV &X, Teuchos::RCP<MV> MX,
00485 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
00486 Teuchos::Array<Teuchos::RCP<const MV> > Q) const {
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496
00497
00498
00499
00500
00501
00502 ScalarType ONE = SCT::one();
00503
00504 int xc = MVT::GetNumberVecs( X );
00505 int xr = MVT::GetVecLength( X );
00506 int nq = Q.length();
00507 std::vector<int> qcs(nq);
00508
00509 if (nq == 0 || xc == 0 || xr == 0) {
00510 return;
00511 }
00512 int qr = MVT::GetVecLength ( *Q[0] );
00513
00514
00515
00516 C.resize(nq);
00517
00518
00519
00520 if (this->_hasOp) {
00521 if (MX == Teuchos::null) {
00522
00523 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
00524 OPT::Apply(*(this->_Op),X,*MX);
00525 this->_OpCounter += MVT::GetNumberVecs(X);
00526 }
00527 }
00528 else {
00529
00530 MX = Teuchos::rcp( &X, false );
00531 }
00532 int mxc = MVT::GetNumberVecs( *MX );
00533 int mxr = MVT::GetVecLength( *MX );
00534
00535
00536 TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
00537 "Anasazi::BasicOrthoManager::projectMat(): MVT returned negative dimensions for X,MX" );
00538
00539 TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr || xr!=qr, std::invalid_argument,
00540 "Anasazi::BasicOrthoManager::projectMat(): Size of X not consistant with MX,Q" );
00541
00542
00543 int baslen = 0;
00544 for (int i=0; i<nq; i++) {
00545 TEST_FOR_EXCEPTION( MVT::GetVecLength( *Q[i] ) != qr, std::invalid_argument,
00546 "Anasazi::BasicOrthoManager::projectMat(): Q lengths not mutually consistant" );
00547 qcs[i] = MVT::GetNumberVecs( *Q[i] );
00548 TEST_FOR_EXCEPTION( qr < qcs[i], std::invalid_argument,
00549 "Anasazi::BasicOrthoManager::projectMat(): Q has less rows than columns" );
00550 baslen += qcs[i];
00551
00552
00553 if ( C[i] == Teuchos::null ) {
00554 C[i] = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) );
00555 }
00556 else {
00557 TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc , std::invalid_argument,
00558 "Anasazi::BasicOrthoManager::projectMat(): Size of Q not consistant with size of C" );
00559 }
00560 }
00561
00562
00563
00564
00565 std::vector<ScalarType> oldDot( xc );
00566 MVT::MvDot( X, *MX, &oldDot );
00567
00568 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq);
00569
00570 for (int i=0; i<nq; i++) {
00571
00572 innerProdMat(*Q[i],X,MX,*C[i]);
00573
00574 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X );
00575
00576
00577 if (this->_hasOp) {
00578 if (xc <= qcs[i]) {
00579 OPT::Apply( *(this->_Op), X, *MX);
00580 this->_OpCounter += MVT::GetNumberVecs(X);
00581 }
00582 else {
00583
00584 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
00585 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
00586 this->_OpCounter += MVT::GetNumberVecs(*Q[i]);
00587 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX );
00588 }
00589 }
00590 }
00591
00592
00593 std::vector<ScalarType> newDot(xc);
00594 MVT::MvDot( X, *MX, &newDot );
00595
00596
00597 for (int j = 0; j < xc; ++j) {
00598
00599 if ( SCT::magnitude(kappa_*newDot[j]) < SCT::magnitude(oldDot[j]) ) {
00600 Teuchos::TimeMonitor lcltimer( *timerReortho_ );
00601 for (int i=0; i<nq; i++) {
00602 Teuchos::SerialDenseMatrix<int,ScalarType> C2(*C[i]);
00603
00604
00605 innerProdMat(*Q[i],X,MX,C2);
00606 *C[i] += C2;
00607 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X );
00608
00609
00610 if (this->_hasOp) {
00611 if (MQ[i].get()) {
00612
00613 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX );
00614 }
00615 else if (xc <= qcs[i]) {
00616
00617 OPT::Apply( *(this->_Op), X, *MX);
00618 this->_OpCounter += MVT::GetNumberVecs(X);
00619 }
00620 }
00621 }
00622 break;
00623 }
00624 }
00625 }
00626
00627
00629
00630
00631 template<class ScalarType, class MV, class OP>
00632 int BasicOrthoManager<ScalarType, MV, OP>::findBasis(
00633 MV &X, Teuchos::RCP<MV> MX,
00634 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
00635 bool completeBasis, int howMany ) const {
00636
00637 using std::cout;
00638 using std::endl;
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659 const ScalarType ONE = SCT::one();
00660 const MagnitudeType ZERO = SCT::magnitude(SCT::zero());
00661 const ScalarType EPS = SCT::eps();
00662
00663 int xc = MVT::GetNumberVecs( X );
00664 int xr = MVT::GetVecLength( X );
00665
00666 if (howMany == -1) {
00667 howMany = xc;
00668 }
00669
00670
00671
00672
00673
00674
00675
00676 if (this->_hasOp) {
00677 if (MX == Teuchos::null) {
00678
00679 MX = MVT::Clone(X,xc);
00680 OPT::Apply(*(this->_Op),X,*MX);
00681 this->_OpCounter += MVT::GetNumberVecs(X);
00682 }
00683 }
00684
00685
00686
00687
00688 if ( B == Teuchos::null ) {
00689 B = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
00690 }
00691
00692 int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc;
00693 int mxr = (this->_hasOp) ? MVT::GetVecLength( *MX ) : xr;
00694
00695
00696 TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
00697 "Anasazi::BasicOrthoManager::findBasis(): X must be non-empty" );
00698 TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
00699 "Anasazi::BasicOrthoManager::findBasis(): Size of X not consistant with size of B" );
00700 TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument,
00701 "Anasazi::BasicOrthoManager::findBasis(): Size of X not consistant with size of MX" );
00702 TEST_FOR_EXCEPTION( xc > xr, std::invalid_argument,
00703 "Anasazi::BasicOrthoManager::findBasis(): Size of X not feasible for normalization" );
00704 TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::invalid_argument,
00705 "Anasazi::BasicOrthoManager::findBasis(): Invalid howMany parameter" );
00706
00707
00708
00709
00710 int xstart = xc - howMany;
00711
00712 for (int j = xstart; j < xc; j++) {
00713
00714
00715 int numX = j;
00716
00717
00718
00719
00720
00721 for (int i=j+1; i<xc; ++i) {
00722 (*B)(i,j) = ZERO;
00723 }
00724
00725
00726 std::vector<int> index(1);
00727 index[0] = j;
00728 Teuchos::RCP<MV> Xj = MVT::CloneView( X, index );
00729 Teuchos::RCP<MV> MXj;
00730 if ((this->_hasOp)) {
00731
00732 MXj = MVT::CloneView( *MX, index );
00733 }
00734 else {
00735
00736 MXj = Xj;
00737 }
00738
00739
00740 std::vector<int> prev_idx( numX );
00741 Teuchos::RCP<const MV> prevX, prevMX;
00742
00743 if (numX > 0) {
00744 for (int i=0; i<numX; ++i) prev_idx[i] = i;
00745 prevX = MVT::CloneView( X, prev_idx );
00746 if (this->_hasOp) {
00747 prevMX = MVT::CloneView( *MX, prev_idx );
00748 }
00749 }
00750
00751 bool rankDef = true;
00752
00753
00754
00755
00756 for (int numTrials = 0; numTrials < 10; numTrials++) {
00757
00758
00759 Teuchos::SerialDenseMatrix<int,ScalarType> product(numX, 1);
00760 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
00761
00762
00763
00764
00765 Teuchos::RCP<MV> oldMXj = MVT::CloneCopy( *MXj );
00766 MVT::MvDot( *Xj, *MXj, &oldDot );
00767
00768 TEST_FOR_EXCEPTION( SCT::real(oldDot[0]) < ZERO, OrthoError,
00769 "Anasazi::BasicOrthoManager::findBasis(): Negative definiteness discovered in inner product" );
00770
00771 if (numX > 0) {
00772
00773
00774
00775 innerProdMat(*prevX,*Xj,MXj,product);
00776
00777
00778
00779 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *Xj );
00780
00781
00782 if (this->_hasOp) {
00783
00784
00785
00786 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *MXj );
00787 }
00788
00789
00790 MVT::MvDot( *Xj, *MXj, &newDot );
00791
00792
00793 if ( SCT::magnitude(kappa_*newDot[0]) < SCT::magnitude(oldDot[0]) ) {
00794
00795
00796 Teuchos::SerialDenseMatrix<int,ScalarType> P2(numX,1);
00797
00798 innerProdMat(*prevX,*Xj,MXj,P2);
00799 product += P2;
00800 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj );
00801 if ((this->_hasOp)) {
00802 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj );
00803 }
00804 }
00805
00806 }
00807
00808
00809 MVT::MvDot( *Xj, *oldMXj, &newDot );
00810
00811
00812 if (numTrials == 0) {
00813 for (int i=0; i<numX; i++) {
00814 (*B)(i,j) = product(i,0);
00815 }
00816 }
00817
00818
00819 #ifdef ANASAZI_BASICORTHO_DEBUG
00820 cout << "olddot: " << SCT::magnitude(oldDot[0]) << " newdot: " << SCT::magnitude(newDot[0]);
00821 #endif
00822 if ( SCT::magnitude(newDot[0]) > SCT::magnitude(oldDot[0]*EPS*EPS) && SCT::real(newDot[0]) > ZERO ) {
00823 #ifdef ANASAZI_BASICORTHO_DEBUG
00824 cout << " ACCEPTED" << endl;
00825 #endif
00826
00827
00828 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
00829
00830 MVT::MvAddMv( ONE/diag, *Xj, ZERO, *Xj, *Xj );
00831 if (this->_hasOp) {
00832
00833 MVT::MvAddMv( ONE/diag, *MXj, ZERO, *MXj, *MXj );
00834 }
00835
00836
00837 if (numTrials == 0) {
00838 (*B)(j,j) = diag;
00839 }
00840
00841
00842 rankDef = false;
00843 break;
00844 }
00845 else {
00846 #ifdef ANASAZI_BASICORTHO_DEBUG
00847 cout << " REJECTED" << endl;
00848 #endif
00849
00850
00851
00852 (*B)(j,j) = ZERO;
00853
00854 if (completeBasis) {
00855
00856 #ifdef ANASAZI_BASICORTHO_DEBUG
00857 cout << "Random for column " << j << endl;
00858 #endif
00859 MVT::MvRandom( *Xj );
00860 if (this->_hasOp) {
00861 OPT::Apply( *(this->_Op), *Xj, *MXj );
00862 this->_OpCounter += MVT::GetNumberVecs(*Xj);
00863 }
00864 }
00865 else {
00866 rankDef = true;
00867 break;
00868 }
00869
00870 }
00871
00872 }
00873
00874
00875 if (rankDef == true) {
00876 MVT::MvInit( *Xj, ZERO );
00877 if (this->_hasOp) {
00878 MVT::MvInit( *MXj, ZERO );
00879 }
00880 TEST_FOR_EXCEPTION( completeBasis, OrthoError,
00881 "Anasazi::BasicOrthoManager::findBasis(): Unable to complete basis" );
00882 return j;
00883 }
00884
00885 }
00886
00887 return xc;
00888 }
00889
00890 }
00891
00892 #endif // ANASAZI_BASIC_ORTHOMANAGER_HPP
00893