314 void expandSearchSpace();
317 void applyOperators();
320 void updateProjections();
323 void solveProjectedEigenproblem();
326 void computeProjectedEigenvectors( Teuchos::SerialDenseMatrix<int,ScalarType> &X );
329 void scaleEigenvectors(
const Teuchos::SerialDenseMatrix<int,ScalarType> &X,
330 Teuchos::SerialDenseMatrix<int,ScalarType> &X_alpha,
331 Teuchos::SerialDenseMatrix<int,ScalarType> &X_beta );
334 void sortValues( std::vector<MagnitudeType> &realParts,
335 std::vector<MagnitudeType> &imagParts,
336 std::vector<int> &permVec,
340 void computeResidual();
343 void computeRitzIndex();
346 void computeRitzVectors();
349 RCP<Eigenproblem<ScalarType,MV,OP> > d_problem;
350 Teuchos::ParameterList d_pl;
359 int d_maxSubspaceDim;
362 std::string d_initType;
364 bool d_relativeConvergence;
367 RCP<OutputManager<ScalarType> > d_outputMan;
368 RCP<OrthoManager<ScalarType,MV> > d_orthoMan;
369 RCP<SortManager<MagnitudeType> > d_sortMan;
370 RCP<StatusTest<ScalarType,MV,OP> > d_tester;
373 std::vector< Value<ScalarType> > d_eigenvalues;
377 std::vector<MagnitudeType> d_resNorms;
383 RCP<MV> d_ritzVecSpace;
385 Teuchos::Array< RCP<const MV> > d_auxVecs;
388 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_VAV;
389 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_VBV;
390 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_S;
391 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_T;
392 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_Q;
393 RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_Z;
396 std::vector<MagnitudeType> d_alphar;
397 std::vector<MagnitudeType> d_alphai;
398 std::vector<MagnitudeType> d_betar;
399 std::vector<int> d_ritzIndex;
400 std::vector<int> d_convergedIndices;
401 std::vector<int> d_expansionIndices;
414 std::string d_testSpace;
422 bool d_ritzIndexValid;
423 bool d_ritzVectorsValid;
430template <
class MagnitudeType,
class MV,
class OP>
435 typedef std::complex<MagnitudeType> ScalarType;
442 Teuchos::ParameterList &pl)
445 MagnitudeType::this_class_is_missing_a_specialization();
456template <
class ScalarType,
class MV,
class OP>
463 Teuchos::ParameterList &pl )
465 TEUCHOS_TEST_FOR_EXCEPTION( problem == Teuchos::null, std::invalid_argument,
"No Eigenproblem given to solver." );
466 TEUCHOS_TEST_FOR_EXCEPTION( outputman == Teuchos::null, std::invalid_argument,
"No OutputManager given to solver." );
467 TEUCHOS_TEST_FOR_EXCEPTION( orthoman == Teuchos::null, std::invalid_argument,
"No OrthoManager given to solver." );
468 TEUCHOS_TEST_FOR_EXCEPTION( sortman == Teuchos::null, std::invalid_argument,
"No SortManager given to solver." );
469 TEUCHOS_TEST_FOR_EXCEPTION( tester == Teuchos::null, std::invalid_argument,
"No StatusTest given to solver." );
470 TEUCHOS_TEST_FOR_EXCEPTION( !problem->isProblemSet(), std::invalid_argument,
"Problem has not been set." );
474 TEUCHOS_TEST_FOR_EXCEPTION( problem->getA()==Teuchos::null && problem->getOperator()==Teuchos::null,
475 std::invalid_argument,
"Either A or Operator must be non-null on Eigenproblem");
476 d_A = problem->getA()!=Teuchos::null ? problem->getA() : problem->getOperator();
477 d_B = problem->getM();
478 d_P = problem->getPrec();
480 d_outputMan = outputman;
482 d_orthoMan = orthoman;
485 d_NEV = d_problem->getNEV();
486 d_initType = d_pl.get<std::string>(
"Initial Guess",
"Random");
487 d_numToPrint = d_pl.get<
int>(
"Print Number of Ritz Values",-1);
488 d_useLeading = d_pl.get<
bool>(
"Use Leading Vectors",
false);
490 if( d_B != Teuchos::null )
495 if( d_P != Teuchos::null )
500 d_testSpace = d_pl.get<std::string>(
"Test Space",
"V");
501 TEUCHOS_TEST_FOR_EXCEPTION( d_testSpace!=
"V" && d_testSpace!=
"AV" && d_testSpace!=
"BV", std::invalid_argument,
502 "Anasazi::GeneralizedDavidson: Test Space must be V, AV, or BV" );
503 TEUCHOS_TEST_FOR_EXCEPTION( d_testSpace==
"V" ?
false : !d_haveB, std::invalid_argument,
504 "Anasazi::GeneralizedDavidson: Test Space must be V for standard eigenvalue problem" );
507 int blockSize = d_pl.get<
int>(
"Block Size",1);
508 int maxSubDim = d_pl.get<
int>(
"Maximum Subspace Dimension",3*d_NEV*blockSize);
510 d_maxSubspaceDim = -1;
511 setSize( blockSize, maxSubDim );
512 d_relativeConvergence = d_pl.get<
bool>(
"Relative Convergence Tolerance",
false);
515 TEUCHOS_TEST_FOR_EXCEPTION( d_blockSize <= 0, std::invalid_argument,
"Block size must be positive");
516 TEUCHOS_TEST_FOR_EXCEPTION( d_maxSubspaceDim <= 0, std::invalid_argument,
"Maximum Subspace Dimension must be positive" );
517 TEUCHOS_TEST_FOR_EXCEPTION( d_problem->getNEV()+2 > pl.get<
int>(
"Maximum Subspace Dimension"),
518 std::invalid_argument,
"Maximum Subspace Dimension must be strictly greater than NEV");
520 std::invalid_argument,
"Maximum Subspace Dimension cannot exceed problem size");
526 d_ritzIndexValid =
false;
527 d_ritzVectorsValid =
false;
534template <
class ScalarType,
class MV,
class OP>
540 d_outputMan->stream(
Warnings) <<
"WARNING: GeneralizedDavidson::iterate called without first calling initialize" << std::endl;
541 d_outputMan->stream(
Warnings) <<
" Default initialization will be performed" << std::endl;
546 if( d_outputMan->isVerbosity(
Debug) )
555 while( d_tester->getStatus() !=
Passed && d_curDim+d_expansionSize <= d_maxSubspaceDim )
565 solveProjectedEigenproblem();
570 int numToSort = std::max(d_blockSize,d_NEV);
571 numToSort = std::min(numToSort,d_curDim);
577 if( d_outputMan->isVerbosity(
Debug) )
591template <
class ScalarType,
class MV,
class OP>
612template <
class ScalarType,
class MV,
class OP>
615 setSize(blockSize,d_maxSubspaceDim);
621template <
class ScalarType,
class MV,
class OP>
624 if( blockSize != d_blockSize || maxSubDim != d_maxSubspaceDim )
626 d_blockSize = blockSize;
627 d_maxSubspaceDim = maxSubDim;
628 d_initialized =
false;
630 d_outputMan->stream(
Debug) <<
" >> Anasazi::GeneralizedDavidson: Allocating eigenproblem"
631 <<
" state with block size of " << d_blockSize
632 <<
" and maximum subspace dimension of " << d_maxSubspaceDim << std::endl;
635 d_alphar.resize(d_maxSubspaceDim);
636 d_alphai.resize(d_maxSubspaceDim);
637 d_betar.resize(d_maxSubspaceDim);
640 int msd = d_maxSubspaceDim;
643 RCP<const MV> initVec = d_problem->getInitVec();
650 d_VAV = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
651 d_S = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
652 d_Q = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
653 d_Z = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
659 d_VBV = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
660 d_T = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
671 d_R =
MVT::Clone(*initVec,std::max(d_blockSize,d_NEV)+1);
672 d_ritzVecSpace =
MVT::Clone(*initVec,std::max(d_blockSize,d_NEV)+1);
687template <
class ScalarType,
class MV,
class OP>
692 d_curDim = (state.
curDim > 0 ? state.
curDim : d_blockSize );
694 d_outputMan->stream(
Debug) <<
" >> Anasazi::GeneralizedDavidson: Initializing state"
695 <<
" with subspace dimension " << d_curDim << std::endl;
698 std::vector<int> initInds(d_curDim);
699 for(
int i=0; i<d_curDim; ++i )
706 bool reset_V =
false;;
714 else if(
MVT::GetNumberVecs(*d_problem->getInitVec()) < d_blockSize || d_initType ==
"Random" )
722 RCP<const MV> initVec =
MVT::CloneView(*d_problem->getInitVec(),initInds);
730 int rank = d_orthoMan->projectAndNormalize( *V1, d_auxVecs );
731 TEUCHOS_TEST_FOR_EXCEPTION( rank < d_blockSize, std::logic_error,
732 "Anasazi::GeneralizedDavidson::initialize(): Error generating initial orthonormal basis" );
735 if( d_outputMan->isVerbosity(
Debug) )
737 d_outputMan->stream(
Debug) <<
" >> Anasazi::GeneralizedDavidson: Error in V^T V == I: "
738 << d_orthoMan->orthonormError( *V1 ) << std::endl;
748 if( state.
AV != d_AV )
759 Teuchos::SerialDenseMatrix<int,ScalarType> VAV1( Teuchos::View, *d_VAV, d_curDim, d_curDim );
762 if( !reset_V && state.
VAV != Teuchos::null && state.
VAV->numRows() >= d_curDim && state.
VAV->numCols() >= d_curDim )
764 if( state.
VAV != d_VAV )
766 Teuchos::SerialDenseMatrix<int,ScalarType> state_VAV( Teuchos::View, *state.
VAV, d_curDim, d_curDim );
767 VAV1.assign( state_VAV );
773 if( d_testSpace ==
"V" )
777 else if( d_testSpace ==
"AV" )
781 else if( d_testSpace ==
"BV" )
797 if( state.
BV != d_BV )
807 Teuchos::SerialDenseMatrix<int,ScalarType> VBV1( Teuchos::View, *d_VBV, d_curDim, d_curDim );
810 if( !reset_V && state.
VBV != Teuchos::null && state.
VBV->numRows() >= d_curDim && state.
VBV->numCols() >= d_curDim )
812 if( state.
VBV != d_VBV )
814 Teuchos::SerialDenseMatrix<int,ScalarType> state_VBV( Teuchos::View, *state.
VBV, d_curDim, d_curDim );
815 VBV1.assign( state_VBV );
821 if( d_testSpace ==
"V" )
825 else if( d_testSpace ==
"AV" )
829 else if( d_testSpace ==
"BV" )
837 solveProjectedEigenproblem();
840 int numToSort = std::max(d_blockSize,d_NEV);
841 numToSort = std::min(numToSort,d_curDim);
848 d_initialized =
true;
854template <
class ScalarType,
class MV,
class OP>
864template <
class ScalarType,
class MV,
class OP>
865std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
874template <
class ScalarType,
class MV,
class OP>
875std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
878 std::vector<int> resIndices(numWanted);
879 for(
int i=0; i<numWanted; ++i )
884 std::vector<MagnitudeType> resNorms;
885 d_orthoMan->norm( *R_checked, resNorms );
893template <
class ScalarType,
class MV,
class OP>
896 std::vector< Value<ScalarType> > ritzValues;
897 for(
int ival=0; ival<d_curDim; ++ival )
900 thisVal.
realpart = d_alphar[ival] / d_betar[ival];
901 if( d_betar[ival] != MT::zero() )
902 thisVal.
imagpart = d_alphai[ival] / d_betar[ival];
906 ritzValues.push_back( thisVal );
921template <
class ScalarType,
class MV,
class OP>
923 const Teuchos::Array< RCP<const MV> > &auxVecs )
928 typename Teuchos::Array< RCP<const MV> >::const_iterator arrItr;
930 for( arrItr=auxVecs.begin(); arrItr!=auxVecs.end(); ++arrItr )
935 d_initialized =
false;
941template <
class ScalarType,
class MV,
class OP>
945 std::vector<MagnitudeType> realRitz(d_curDim), imagRitz(d_curDim);
947 for(
int i=0; i<d_curDim; ++i )
949 realRitz[i] = ritzVals[i].realpart;
950 imagRitz[i] = ritzVals[i].imagpart;
953 std::vector<int> permVec;
954 sortValues( realRitz, imagRitz, permVec, d_curDim );
956 std::vector<int> sel(d_curDim,0);
957 for(
int ii=0; ii<numWanted; ++ii )
958 sel[ permVec[ii] ]=1;
965 int work_size=10*d_maxSubspaceDim+16;
966 std::vector<ScalarType> work(work_size);
972 Teuchos::LAPACK<int,ScalarType> lapack;
973 lapack.TGSEN( ijob, wantq, wantz, &sel[0], d_curDim, d_S->values(), d_S->stride(), d_T->values(), d_T->stride(),
974 &d_alphar[0], &d_alphai[0], &d_betar[0], d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(),
975 &sdim, 0, 0, 0, &work[0], work_size, &iwork, iwork_size, &info );
977 d_ritzIndexValid =
false;
978 d_ritzVectorsValid =
false;
980 std::stringstream ss;
981 ss <<
"Anasazi::GeneralizedDavidson: TGSEN returned error code " << info << std::endl;
982 TEUCHOS_TEST_FOR_EXCEPTION( info<0, std::runtime_error, ss.str() );
988 d_outputMan->stream(
Warnings) <<
"WARNING: " << ss.str() << std::endl;
989 d_outputMan->stream(
Warnings) <<
" Problem may not be correctly sorted" << std::endl;
993 char getCondNum =
'N';
996 int work_size = d_curDim;
997 std::vector<ScalarType> work(work_size);
1002 Teuchos::LAPACK<int,ScalarType> lapack;
1003 lapack.TRSEN (getCondNum, getQ, &sel[0], d_curDim, d_S->values (),
1004 d_S->stride (), d_Z->values (), d_Z->stride (),
1005 &d_alphar[0], &d_alphai[0], &subDim, 0, 0, &work[0],
1006 work_size, &iwork, iwork_size, &info );
1008 std::fill( d_betar.begin(), d_betar.end(), ST::one() );
1010 d_ritzIndexValid =
false;
1011 d_ritzVectorsValid =
false;
1013 TEUCHOS_TEST_FOR_EXCEPTION(
1014 info < 0, std::runtime_error,
"Anasazi::GeneralizedDavidson::"
1015 "sortProblem: LAPACK routine TRSEN returned error code INFO = "
1016 << info <<
" < 0. This indicates that argument " << -info
1017 <<
" (counting starts with one) of TRSEN had an illegal value.");
1023 TEUCHOS_TEST_FOR_EXCEPTION(
1024 info == 1, std::runtime_error,
"Anasazi::GeneralizedDavidson::"
1025 "sortProblem: LAPACK routine TRSEN returned error code INFO = 1. "
1026 "This indicates that the reordering failed because some eigenvalues "
1027 "are too close to separate (the problem is very ill-conditioned).");
1029 TEUCHOS_TEST_FOR_EXCEPTION(
1030 info > 1, std::logic_error,
"Anasazi::GeneralizedDavidson::"
1031 "sortProblem: LAPACK routine TRSEN returned error code INFO = "
1032 << info <<
" > 1. This should not be possible. It may indicate an "
1033 "error either in LAPACK itself, or in Teuchos' LAPACK wrapper.");
1045template <
class ScalarType,
class MV,
class OP>
1046void GeneralizedDavidson<ScalarType,MV,OP>::expandSearchSpace()
1050 std::vector<int> newIndices(d_expansionSize);
1051 for(
int i=0; i<d_expansionSize; ++i )
1053 newIndices[i] = d_curDim+i;
1057 std::vector<int> curIndices(d_curDim);
1058 for(
int i=0; i<d_curDim; ++i )
1062 RCP<MV> V_new = MVT::CloneViewNonConst( *d_V, newIndices);
1063 RCP<const MV> V_cur = MVT::CloneView( *d_V, curIndices);
1064 RCP<const MV> R_active = MVT::CloneView( *d_R, d_expansionIndices);
1069 OPT::Apply( *d_P, *R_active, *V_new );
1074 MVT::SetBlock( *R_active, newIndices, *d_V );
1078 Teuchos::Array< RCP<const MV> > against = d_auxVecs;
1079 against.push_back( V_cur );
1080 int rank = d_orthoMan->projectAndNormalize(*V_new,against);
1082 if( d_outputMan->isVerbosity(
Debug) )
1084 std::vector<int> allIndices(d_curDim+d_expansionSize);
1085 for(
int i=0; i<d_curDim+d_expansionSize; ++i )
1088 RCP<const MV> V_all = MVT::CloneView( *d_V, allIndices );
1090 d_outputMan->stream(
Debug) <<
" >> Anasazi::GeneralizedDavidson: Error in V^T V == I: "
1091 << d_orthoMan->orthonormError( *V_all ) << std::endl;
1094 TEUCHOS_TEST_FOR_EXCEPTION( rank != d_expansionSize, std::runtime_error,
1095 "Anasazi::GeneralizedDavidson::ExpandSearchSpace(): Orthonormalization of new vectors failed" );
1102template <
class ScalarType,
class MV,
class OP>
1103void GeneralizedDavidson<ScalarType,MV,OP>::applyOperators()
1106 std::vector<int> newIndices(d_expansionSize);
1107 for(
int i=0; i<d_expansionSize; ++i )
1108 newIndices[i] = d_curDim+i;
1111 RCP<const MV> V_new = MVT::CloneView( *d_V, newIndices );
1112 RCP<MV> AV_new = MVT::CloneViewNonConst( *d_AV, newIndices );
1115 OPT::Apply( *d_A, *V_new, *AV_new );
1116 d_opApplies += MVT::GetNumberVecs( *V_new );
1121 RCP<MV> BV_new = MVT::CloneViewNonConst( *d_BV, newIndices );
1122 OPT::Apply( *d_B, *V_new, *BV_new );
1129template <
class ScalarType,
class MV,
class OP>
1130void GeneralizedDavidson<ScalarType,MV,OP>::updateProjections()
1133 std::vector<int> newIndices(d_expansionSize);
1134 for(
int i=0; i<d_expansionSize; ++i )
1135 newIndices[i] = d_curDim+i;
1137 std::vector<int> curIndices(d_curDim);
1138 for(
int i=0; i<d_curDim; ++i )
1141 std::vector<int> allIndices(d_curDim+d_expansionSize);
1142 for(
int i=0; i<d_curDim+d_expansionSize; ++i )
1146 RCP<const MV> W_new, W_all;
1147 if( d_testSpace ==
"V" )
1149 W_new = MVT::CloneView(*d_V, newIndices );
1150 W_all = MVT::CloneView(*d_V, allIndices );
1152 else if( d_testSpace ==
"AV" )
1154 W_new = MVT::CloneView(*d_AV, newIndices );
1155 W_all = MVT::CloneView(*d_AV, allIndices );
1157 else if( d_testSpace ==
"BV" )
1159 W_new = MVT::CloneView(*d_BV, newIndices );
1160 W_all = MVT::CloneView(*d_BV, allIndices );
1164 RCP<const MV> AV_new = MVT::CloneView(*d_AV, newIndices);
1165 RCP<const MV> AV_current = MVT::CloneView(*d_AV, curIndices);
1168 Teuchos::SerialDenseMatrix<int,ScalarType> VAV_lastrow( Teuchos::View, *d_VAV, d_expansionSize, d_curDim, d_curDim, 0 );
1169 MVT::MvTransMv( ST::one(), *W_new, *AV_current, VAV_lastrow );
1172 Teuchos::SerialDenseMatrix<int,ScalarType> VAV_lastcol( Teuchos::View, *d_VAV, d_curDim+d_expansionSize, d_expansionSize, 0, d_curDim );
1173 MVT::MvTransMv( ST::one(), *W_all, *AV_new, VAV_lastcol );
1178 RCP<const MV> BV_new = MVT::CloneView(*d_BV, newIndices);
1179 RCP<const MV> BV_current = MVT::CloneView(*d_BV, curIndices);
1182 Teuchos::SerialDenseMatrix<int,ScalarType> VBV_lastrow( Teuchos::View, *d_VBV, d_expansionSize, d_curDim, d_curDim, 0 );
1183 MVT::MvTransMv( ST::one(), *W_new, *BV_current, VBV_lastrow );
1186 Teuchos::SerialDenseMatrix<int,ScalarType> VBV_lastcol( Teuchos::View, *d_VBV, d_curDim+d_expansionSize, d_expansionSize, 0, d_curDim );
1187 MVT::MvTransMv( ST::one(), *W_all, *BV_new, VBV_lastcol );
1191 d_curDim += d_expansionSize;
1193 d_ritzIndexValid =
false;
1194 d_ritzVectorsValid =
false;
1200template <
class ScalarType,
class MV,
class OP>
1201void GeneralizedDavidson<ScalarType,MV,OP>::solveProjectedEigenproblem()
1208 d_S->assign(*d_VAV);
1209 d_T->assign(*d_VBV);
1212 char leftVecs =
'V';
1213 char rightVecs =
'V';
1214 char sortVals =
'N';
1217 Teuchos::LAPACK<int,ScalarType> lapack;
1221 std::vector<ScalarType> work(1);
1222 lapack.GGES( leftVecs, rightVecs, sortVals, NULL, d_curDim, d_S->values(), d_S->stride(),
1223 d_T->values(), d_T->stride(), &sdim, &d_alphar[0], &d_alphai[0], &d_betar[0],
1224 d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(), &work[0], work_size, 0, &info );
1226 work_size = work[0];
1227 work.resize(work_size);
1228 lapack.GGES( leftVecs, rightVecs, sortVals, NULL, d_curDim, d_S->values(), d_S->stride(),
1229 d_T->values(), d_T->stride(), &sdim, &d_alphar[0], &d_alphai[0], &d_betar[0],
1230 d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(), &work[0], work_size, 0, &info );
1232 d_ritzIndexValid =
false;
1233 d_ritzVectorsValid =
false;
1235 std::stringstream ss;
1236 ss <<
"Anasazi::GeneralizedDavidson: GGES returned error code " << info << std::endl;
1237 TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1243 d_S->assign(*d_VAV);
1248 int work_size = 3*d_curDim;
1249 std::vector<ScalarType> work(work_size);
1252 Teuchos::LAPACK<int,ScalarType> lapack;
1253 lapack.GEES( vecs, d_curDim, d_S->values(), d_S->stride(), &sdim, &d_alphar[0], &d_alphai[0],
1254 d_Z->values(), d_Z->stride(), &work[0], work_size, 0, 0, &info);
1256 std::fill( d_betar.begin(), d_betar.end(), ST::one() );
1258 d_ritzIndexValid =
false;
1259 d_ritzVectorsValid =
false;
1261 std::stringstream ss;
1262 ss <<
"Anasazi::GeneralizedDavidson: GEES returned error code " << info << std::endl;
1263 TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1275template <
class ScalarType,
class MV,
class OP>
1276void GeneralizedDavidson<ScalarType,MV,OP>::computeRitzIndex()
1278 if( d_ritzIndexValid )
1281 d_ritzIndex.resize( d_curDim );
1283 while( i < d_curDim )
1285 if( d_alphai[i] == ST::zero() )
1293 d_ritzIndex[i+1] = -1;
1297 d_ritzIndexValid =
true;
1308template <
class ScalarType,
class MV,
class OP>
1309void GeneralizedDavidson<ScalarType,MV,OP>::computeRitzVectors()
1311 if( d_ritzVectorsValid )
1318 std::vector<int> checkedIndices(d_residualSize);
1319 for(
int ii=0; ii<d_residualSize; ++ii )
1320 checkedIndices[ii] = ii;
1323 Teuchos::SerialDenseMatrix<int,ScalarType> X(Teuchos::Copy,*d_Z,d_curDim,d_curDim);
1324 computeProjectedEigenvectors( X );
1327 Teuchos::SerialDenseMatrix<int,ScalarType> X_wanted(Teuchos::View,X,d_curDim,d_residualSize);
1330 d_ritzVecs = MVT::CloneViewNonConst( *d_ritzVecSpace, checkedIndices );
1332 std::vector<int> curIndices(d_curDim);
1333 for(
int i=0; i<d_curDim; ++i )
1336 RCP<const MV> V_current = MVT::CloneView( *d_V, curIndices );
1339 MVT::MvTimesMatAddMv(ST::one(),*V_current,X_wanted,ST::zero(),*d_ritzVecs);
1342 std::vector<MagnitudeType> scale(d_residualSize);
1343 MVT::MvNorm( *d_ritzVecs, scale );
1344 Teuchos::LAPACK<int,ScalarType> lapack;
1345 for(
int i=0; i<d_residualSize; ++i )
1347 if( d_ritzIndex[i] == 0 )
1349 scale[i] = 1.0/scale[i];
1351 else if( d_ritzIndex[i] == 1 )
1353 MagnitudeType nrm = lapack.LAPY2(scale[i],scale[i+1]);
1355 scale[i+1] = 1.0/nrm;
1358 MVT::MvScale( *d_ritzVecs, scale );
1360 d_ritzVectorsValid =
true;
1367template <
class ScalarType,
class MV,
class OP>
1368void GeneralizedDavidson<ScalarType,MV,OP>::sortValues( std::vector<MagnitudeType> &realParts,
1369 std::vector<MagnitudeType> &imagParts,
1370 std::vector<int> &permVec,
1375 TEUCHOS_TEST_FOR_EXCEPTION( (
int) realParts.size()<N, std::runtime_error,
1376 "Anasazi::GeneralizedDavidson::SortValues: Number of requested sorted values greater than vector length." );
1377 TEUCHOS_TEST_FOR_EXCEPTION( (
int) imagParts.size()<N, std::runtime_error,
1378 "Anasazi::GeneralizedDavidson::SortValues: Number of requested sorted values greater than vector length." );
1380 RCP< std::vector<int> > rcpPermVec = Teuchos::rcpFromRef(permVec);
1382 d_sortMan->sort( realParts, imagParts, rcpPermVec, N );
1384 d_ritzIndexValid =
false;
1385 d_ritzVectorsValid =
false;
1399template <
class ScalarType,
class MV,
class OP>
1400void GeneralizedDavidson<ScalarType,MV,OP>::computeProjectedEigenvectors(
1401 Teuchos::SerialDenseMatrix<int,ScalarType> &X )
1403 int N = X.numRows();
1406 Teuchos::SerialDenseMatrix<int,ScalarType> S(Teuchos::Copy, *d_S, N, N);
1407 Teuchos::SerialDenseMatrix<int,ScalarType> T(Teuchos::Copy, *d_T, N, N);
1408 Teuchos::SerialDenseMatrix<int,ScalarType> VL(Teuchos::Copy, *d_Q, N, N);
1410 char whichVecs =
'R';
1412 int work_size = 6*d_maxSubspaceDim;
1413 std::vector<ScalarType> work(work_size,ST::zero());
1417 Teuchos::LAPACK<int,ScalarType> lapack;
1418 lapack.TGEVC( whichVecs, howMany, 0, N, S.values(), S.stride(), T.values(), T.stride(),
1419 VL.values(), VL.stride(), X.values(), X.stride(), N, &M, &work[0], &info );
1421 std::stringstream ss;
1422 ss <<
"Anasazi::GeneralizedDavidson: TGEVC returned error code " << info << std::endl;
1423 TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1427 Teuchos::SerialDenseMatrix<int,ScalarType> S(Teuchos::Copy, *d_S, N, N);
1428 Teuchos::SerialDenseMatrix<int,ScalarType> VL(Teuchos::Copy, *d_Z, N, N);
1430 char whichVecs =
'R';
1433 std::vector<ScalarType> work(3*N);
1437 Teuchos::LAPACK<int,ScalarType> lapack;
1439 lapack.TREVC( whichVecs, howMany, &sel, N, S.values(), S.stride(), VL.values(), VL.stride(),
1440 X.values(), X.stride(), N, &m, &work[0], &info );
1442 std::stringstream ss;
1443 ss <<
"Anasazi::GeneralizedDavidson: TREVC returned error code " << info << std::endl;
1444 TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1451template <
class ScalarType,
class MV,
class OP>
1452void GeneralizedDavidson<ScalarType,MV,OP>::scaleEigenvectors(
1453 const Teuchos::SerialDenseMatrix<int,ScalarType> &X,
1454 Teuchos::SerialDenseMatrix<int,ScalarType> &X_alpha,
1455 Teuchos::SerialDenseMatrix<int,ScalarType> &X_beta )
1457 int Nr = X.numRows();
1458 int Nc = X.numCols();
1460 TEUCHOS_TEST_FOR_EXCEPTION( Nr>d_curDim, std::logic_error,
1461 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: Matrix size exceeds current dimension");
1462 TEUCHOS_TEST_FOR_EXCEPTION( Nc>d_curDim, std::logic_error,
1463 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: Matrix size exceeds current dimension");
1464 TEUCHOS_TEST_FOR_EXCEPTION( X_alpha.numRows()!=Nr, std::logic_error,
1465 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of rows in Xalpha does not match X");
1466 TEUCHOS_TEST_FOR_EXCEPTION( X_alpha.numCols()!=Nc, std::logic_error,
1467 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of cols in Xalpha does not match X");
1468 TEUCHOS_TEST_FOR_EXCEPTION( X_beta.numRows()!=Nr, std::logic_error,
1469 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of rows in Xbeta does not match X");
1470 TEUCHOS_TEST_FOR_EXCEPTION( X_beta.numCols()!=Nc, std::logic_error,
1471 "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of cols in Xbeta does not match X");
1475 Teuchos::SerialDenseMatrix<int,ScalarType> Alpha(Nc,Nc,
true);
1476 Teuchos::SerialDenseMatrix<int,ScalarType> Beta(Nc,Nc,
true);
1480 for(
int i=0; i<Nc; ++i )
1482 if( d_ritzIndex[i] == 0 )
1484 Alpha(i,i) = d_alphar[i];
1485 Beta(i,i) = d_betar[i];
1487 else if( d_ritzIndex[i] == 1 )
1489 Alpha(i,i) = d_alphar[i];
1490 Alpha(i,i+1) = d_alphai[i];
1491 Beta(i,i) = d_betar[i];
1495 Alpha(i,i-1) = d_alphai[i];
1496 Alpha(i,i) = d_alphar[i];
1497 Beta(i,i) = d_betar[i];
1504 err = X_alpha.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, ST::one(), X, Alpha, ST::zero() );
1505 std::stringstream astream;
1506 astream <<
"GeneralizedDavidson::ScaleEigenvectors: multiply returned error code " << err;
1507 TEUCHOS_TEST_FOR_EXCEPTION( err!=0, std::runtime_error, astream.str() );
1510 err = X_beta.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, ST::one(), X, Beta, ST::zero() );
1511 std::stringstream bstream;
1512 bstream <<
"GeneralizedDavidson::ScaleEigenvectors: multiply returned error code " << err;
1513 TEUCHOS_TEST_FOR_EXCEPTION( err!=0, std::runtime_error, bstream.str() );
1519template <
class ScalarType,
class MV,
class OP>
1520void GeneralizedDavidson<ScalarType,MV,OP>::computeResidual()
1525 d_residualSize = std::max( d_blockSize, d_NEV );
1526 if( d_curDim < d_residualSize )
1528 d_residualSize = d_curDim;
1530 else if( d_ritzIndex[d_residualSize-1] == 1 )
1536 std::vector<int> residualIndices(d_residualSize);
1537 for(
int i=0; i<d_residualSize; ++i )
1538 residualIndices[i] = i;
1541 Teuchos::SerialDenseMatrix<int,ScalarType> X(Teuchos::Copy,*d_Z,d_curDim,d_curDim);
1544 computeProjectedEigenvectors( X );
1547 Teuchos::SerialDenseMatrix<int,ScalarType> X_alpha(d_curDim,d_residualSize);
1548 Teuchos::SerialDenseMatrix<int,ScalarType> X_beta(d_curDim,d_residualSize);
1551 Teuchos::SerialDenseMatrix<int,ScalarType> X_wanted(Teuchos::View, X, d_curDim, d_residualSize);
1554 scaleEigenvectors( X_wanted, X_alpha, X_beta );
1557 RCP<MV> R_active = MVT::CloneViewNonConst( *d_R, residualIndices );
1560 std::vector<int> activeIndices(d_curDim);
1561 for(
int i=0; i<d_curDim; ++i )
1565 RCP<const MV> AV_active = MVT::CloneView( *d_AV, activeIndices );
1566 MVT::MvTimesMatAddMv(ST::one(),*AV_active, X_beta, ST::zero(),*R_active);
1570 RCP<const MV> BV_active = MVT::CloneView( *d_BV, activeIndices );
1571 MVT::MvTimesMatAddMv(ST::one(),*BV_active, X_alpha,-ST::one(), *R_active);
1575 RCP<const MV> V_active = MVT::CloneView( *d_V, activeIndices );
1576 MVT::MvTimesMatAddMv(ST::one(),*V_active, X_alpha,-ST::one(), *R_active);
1591 Teuchos::LAPACK<int,ScalarType> lapack;
1592 Teuchos::BLAS<int,ScalarType> blas;
1593 std::vector<MagnitudeType> resScaling(d_residualSize);
1594 for(
int icol=0; icol<d_residualSize; ++icol )
1596 if( d_ritzIndex[icol] == 0 )
1598 MagnitudeType Xnrm = blas.NRM2( d_curDim, X_wanted[icol], 1);
1599 MagnitudeType ABscaling = d_relativeConvergence ? d_alphar[icol] : d_betar[icol];
1600 resScaling[icol] = MT::one() / (Xnrm * ABscaling);
1602 else if( d_ritzIndex[icol] == 1 )
1604 MagnitudeType Xnrm1 = blas.NRM2( d_curDim, X_wanted[icol], 1 );
1605 MagnitudeType Xnrm2 = blas.NRM2( d_curDim, X_wanted[icol+1], 1 );
1606 MagnitudeType Xnrm = lapack.LAPY2(Xnrm1,Xnrm2);
1607 MagnitudeType ABscaling = d_relativeConvergence ? lapack.LAPY2(d_alphar[icol],d_alphai[icol])
1609 resScaling[icol] = MT::one() / (Xnrm * ABscaling);
1610 resScaling[icol+1] = MT::one() / (Xnrm * ABscaling);
1613 MVT::MvScale( *R_active, resScaling );
1616 d_resNorms.resize(d_residualSize);
1617 MVT::MvNorm(*R_active,d_resNorms);
1624 for(
int i=0; i<d_residualSize; ++i )
1626 if( d_ritzIndex[i] == 1 )
1628 MagnitudeType nrm = lapack.LAPY2(d_resNorms[i],d_resNorms[i+1]);
1629 d_resNorms[i] = nrm;
1630 d_resNorms[i+1] = nrm;
1635 d_tester->checkStatus(
this);
1638 if( d_useLeading || d_blockSize >= d_NEV )
1640 d_expansionSize=d_blockSize;
1641 if( d_ritzIndex[d_blockSize-1]==1 )
1644 d_expansionIndices.resize(d_expansionSize);
1645 for(
int i=0; i<d_expansionSize; ++i )
1646 d_expansionIndices[i] = i;
1650 std::vector<int> convergedVectors = d_tester->whichVecs();
1654 for( startVec=0; startVec<d_residualSize; ++startVec )
1656 if( std::find(convergedVectors.begin(),convergedVectors.end(),startVec)==convergedVectors.end() )
1662 int endVec = startVec + d_blockSize - 1;
1663 if( endVec > (d_residualSize-1) )
1665 endVec = d_residualSize-1;
1666 startVec = d_residualSize-d_blockSize;
1670 if( d_ritzIndex[startVec]==-1 )
1676 if( d_ritzIndex[endVec]==1 )
1679 d_expansionSize = 1+endVec-startVec;
1680 d_expansionIndices.resize(d_expansionSize);
1681 for(
int i=0; i<d_expansionSize; ++i )
1682 d_expansionIndices[i] = startVec+i;
1689template <
class ScalarType,
class MV,
class OP>
1694 myout.setf(std::ios::scientific, std::ios::floatfield);
1697 myout <<
"================================================================================" << endl;
1699 myout <<
" GeneralizedDavidson Solver Solver Status" << endl;
1701 myout <<
"The solver is "<<(d_initialized ?
"initialized." :
"not initialized.") << endl;
1702 myout <<
"The number of iterations performed is " << d_iteration << endl;
1703 myout <<
"The number of operator applies performed is " << d_opApplies << endl;
1704 myout <<
"The block size is " << d_expansionSize << endl;
1705 myout <<
"The current basis size is " << d_curDim << endl;
1706 myout <<
"The number of requested eigenvalues is " << d_NEV << endl;
1707 myout <<
"The number of converged values is " << d_tester->howMany() << endl;
1710 myout.setf(std::ios_base::right, std::ios_base::adjustfield);
1714 myout <<
"CURRENT RITZ VALUES" << endl;
1716 myout << std::setw(24) <<
"Ritz Value"
1717 << std::setw(30) <<
"Residual Norm" << endl;
1718 myout <<
"--------------------------------------------------------------------------------" << endl;
1719 if( d_residualSize > 0 )
1721 std::vector<MagnitudeType> realRitz(d_curDim), imagRitz(d_curDim);
1722 std::vector< Value<ScalarType> > ritzVals =
getRitzValues();
1723 for(
int i=0; i<d_curDim; ++i )
1725 realRitz[i] = ritzVals[i].realpart;
1726 imagRitz[i] = ritzVals[i].imagpart;
1728 std::vector<int> permvec;
1729 sortValues( realRitz, imagRitz, permvec, d_curDim );
1731 int numToPrint = std::max( d_numToPrint, d_NEV );
1732 numToPrint = std::min( d_curDim, numToPrint );
1739 if( d_ritzIndex[permvec[numToPrint-1]] != 0 )
1743 while( i<numToPrint )
1745 if( imagRitz[i] == ST::zero() )
1747 myout << std::setw(15) << realRitz[i];
1748 myout <<
" + i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1749 if( i < d_residualSize )
1750 myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1752 myout <<
" Not Computed" << endl;
1759 myout << std::setw(15) << realRitz[i];
1760 myout <<
" + i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1761 if( i < d_residualSize )
1762 myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1764 myout <<
" Not Computed" << endl;
1767 myout << std::setw(15) << realRitz[i];
1768 myout <<
" - i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1769 if( i < d_residualSize )
1770 myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1772 myout <<
" Not Computed" << endl;
1780 myout << std::setw(20) <<
"[ NONE COMPUTED ]" << endl;
1784 myout <<
"================================================================================" << endl;