47 #ifndef BELOS_ICGS_ORTHOMANAGER_HPP
48 #define BELOS_ICGS_ORTHOMANAGER_HPP
65 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
66 #ifdef BELOS_TEUCHOS_TIME_MONITOR
68 #endif // BELOS_TEUCHOS_TIME_MONITOR
73 template<
class ScalarType,
class MV,
class OP>
77 template<
class ScalarType,
class MV,
class OP>
80 template<
class ScalarType,
class MV,
class OP>
103 max_ortho_steps_( max_ortho_steps ),
105 sing_tol_( sing_tol ),
108 #ifdef BELOS_TEUCHOS_TIME_MONITOR
109 std::stringstream ss;
110 ss << label_ +
": ICGS[" << max_ortho_steps_ <<
"]";
112 std::string orthoLabel = ss.str() +
": Orthogonalization";
115 std::string updateLabel = ss.str() +
": Ortho (Update)";
118 std::string normLabel = ss.str() +
": Ortho (Norm)";
121 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
128 const std::string& label =
"Belos",
138 #ifdef BELOS_TEUCHOS_TIME_MONITOR
139 std::stringstream ss;
140 ss << label_ +
": ICGS[" << max_ortho_steps_ <<
"]";
142 std::string orthoLabel = ss.str() +
": Orthogonalization";
145 std::string updateLabel = ss.str() +
": Ortho (Update)";
148 std::string normLabel = ss.str() +
": Ortho (Norm)";
151 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
168 using Teuchos::parameterList;
172 RCP<ParameterList> params;
174 params = parameterList (*defaultParams);
189 int maxNumOrthogPasses;
190 MagnitudeType blkTol;
191 MagnitudeType singTol;
194 maxNumOrthogPasses = params->
get<
int> (
"maxNumOrthogPasses");
195 }
catch (InvalidParameterName&) {
196 maxNumOrthogPasses = defaultParams->get<
int> (
"maxNumOrthogPasses");
197 params->set (
"maxNumOrthogPasses", maxNumOrthogPasses);
208 blkTol = params->get<MagnitudeType> (
"blkTol");
209 }
catch (InvalidParameterName&) {
211 blkTol = params->get<MagnitudeType> (
"depTol");
214 params->remove (
"depTol");
215 }
catch (InvalidParameterName&) {
216 blkTol = defaultParams->get<MagnitudeType> (
"blkTol");
218 params->set (
"blkTol", blkTol);
222 singTol = params->get<MagnitudeType> (
"singTol");
223 }
catch (InvalidParameterName&) {
224 singTol = defaultParams->get<MagnitudeType> (
"singTol");
225 params->set (
"singTol", singTol);
228 max_ortho_steps_ = maxNumOrthogPasses;
238 if (defaultParams_.
is_null()) {
239 defaultParams_ = Belos::getICGSDefaultParameters<ScalarType, MV, OP>();
242 return defaultParams_;
256 using Teuchos::parameterList;
261 RCP<ParameterList> params = parameterList (*defaultParams);
282 params->
set (
"blkTol", blk_tol);
296 params->
set (
"singTol", sing_tol);
298 sing_tol_ = sing_tol;
481 void setLabel(
const std::string& label);
485 const std::string&
getLabel()
const {
return label_; }
511 int max_ortho_steps_;
513 MagnitudeType blk_tol_;
515 MagnitudeType sing_tol_;
519 #ifdef BELOS_TEUCHOS_TIME_MONITOR
521 #endif // BELOS_TEUCHOS_TIME_MONITOR
529 bool completeBasis,
int howMany = -1 )
const;
561 template<
class ScalarType,
class MV,
class OP>
564 template<
class ScalarType,
class MV,
class OP>
565 const typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType
570 template<
class ScalarType,
class MV,
class OP>
571 const typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType
575 template<
class ScalarType,
class MV,
class OP>
578 template<
class ScalarType,
class MV,
class OP>
579 const typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType
583 template<
class ScalarType,
class MV,
class OP>
584 const typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType
590 template<
class ScalarType,
class MV,
class OP>
593 if (label != label_) {
595 #ifdef BELOS_TEUCHOS_TIME_MONITOR
596 std::stringstream ss;
597 ss << label_ +
": ICGS[" << max_ortho_steps_ <<
"]";
599 std::string orthoLabel = ss.str() +
": Orthogonalization";
602 std::string updateLabel = ss.str() +
": Ortho (Update)";
605 std::string normLabel = ss.str() +
": Ortho (Norm)";
608 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
616 template<
class ScalarType,
class MV,
class OP>
619 const ScalarType ONE = SCT::one();
620 int rank = MVT::GetNumberVecs(X);
623 for (
int i=0; i<rank; i++) {
631 template<
class ScalarType,
class MV,
class OP>
634 int r1 = MVT::GetNumberVecs(X1);
635 int r2 = MVT::GetNumberVecs(X2);
643 template<
class ScalarType,
class MV,
class OP>
659 typedef typename Array< RCP< const MV > >::size_type size_type;
661 #ifdef BELOS_TEUCHOS_TIME_MONITOR
665 ScalarType ONE = SCT::one();
666 const MagnitudeType ZERO = MGT::zero();
669 int xc = MVT::GetNumberVecs( X );
670 ptrdiff_t xr = MVT::GetGlobalLength( X );
677 B =
rcp (
new serial_dense_matrix_type (xc, xc));
687 for (size_type k = 0; k < nq; ++k)
689 const int numRows = MVT::GetNumberVecs (*Q[k]);
690 const int numCols = xc;
693 C[k] =
rcp (
new serial_dense_matrix_type (numRows, numCols));
694 else if (C[k]->numRows() != numRows || C[k]->numCols() != numCols)
696 int err = C[k]->reshape (numRows, numCols);
698 "IMGS orthogonalization: failed to reshape "
699 "C[" << k <<
"] (the array of block "
700 "coefficients resulting from projecting X "
701 "against Q[1:" << nq <<
"]).");
707 if (MX == Teuchos::null) {
709 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
710 OPT::Apply(*(this->_Op),X,*MX);
718 int mxc = MVT::GetNumberVecs( *MX );
719 ptrdiff_t mxr = MVT::GetGlobalLength( *MX );
722 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
"Belos::ICGSOrthoManager::projectAndNormalize(): X must be non-empty" );
725 for (
int i=0; i<nq; i++) {
726 numbas += MVT::GetNumberVecs( *Q[i] );
731 "Belos::ICGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of B" );
734 "Belos::ICGSOrthoManager::projectAndNormalize(): MVT returned negative dimensions for X,MX" );
737 "Belos::ICGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of MX" );
743 bool dep_flg =
false;
749 dep_flg = blkOrtho1( X, MX, C, Q );
752 if ( B == Teuchos::null ) {
755 std::vector<ScalarType> diag(xc);
757 #ifdef BELOS_TEUCHOS_TIME_MONITOR
760 MVT::MvDot( X, *MX, diag );
762 (*B)(0,0) = SCT::squareroot(SCT::magnitude(diag[0]));
764 if (SCT::magnitude((*B)(0,0)) > ZERO) {
766 MVT::MvScale( X, ONE/(*B)(0,0) );
769 MVT::MvScale( *MX, ONE/(*B)(0,0) );
777 tmpX = MVT::CloneCopy(X);
779 tmpMX = MVT::CloneCopy(*MX);
783 dep_flg = blkOrtho( X, MX, C, Q );
789 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q );
792 MVT::Assign( *tmpX, X );
794 MVT::Assign( *tmpMX, *MX );
799 rank = findBasis( X, MX, B,
false );
804 rank = blkOrthoSing( *tmpX, tmpMX, C, B, Q );
807 MVT::Assign( *tmpX, X );
809 MVT::Assign( *tmpMX, *MX );
817 "Belos::ICGSOrthoManager::projectAndNormalize(): Debug error in rank variable." );
827 template<
class ScalarType,
class MV,
class OP>
832 #ifdef BELOS_TEUCHOS_TIME_MONITOR
837 return findBasis(X, MX, B,
true);
844 template<
class ScalarType,
class MV,
class OP>
864 #ifdef BELOS_TEUCHOS_TIME_MONITOR
868 int xc = MVT::GetNumberVecs( X );
869 ptrdiff_t xr = MVT::GetGlobalLength( X );
871 std::vector<int> qcs(nq);
873 if (nq == 0 || xc == 0 || xr == 0) {
876 ptrdiff_t qr = MVT::GetGlobalLength ( *Q[0] );
885 if (MX == Teuchos::null) {
887 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
888 OPT::Apply(*(this->_Op),X,*MX);
895 int mxc = MVT::GetNumberVecs( *MX );
896 ptrdiff_t mxr = MVT::GetGlobalLength( *MX );
900 "Belos::ICGSOrthoManager::project(): MVT returned negative dimensions for X,MX" );
903 "Belos::ICGSOrthoManager::project(): Size of X not consistant with MX,Q" );
907 for (
int i=0; i<nq; i++) {
909 "Belos::ICGSOrthoManager::project(): Q lengths not mutually consistant" );
910 qcs[i] = MVT::GetNumberVecs( *Q[i] );
912 "Belos::ICGSOrthoManager::project(): Q has less rows than columns" );
916 if ( C[i] == Teuchos::null ) {
921 "Belos::ICGSOrthoManager::project(): Size of Q not consistant with size of C" );
926 blkOrtho( X, MX, C, Q );
933 template<
class ScalarType,
class MV,
class OP>
956 const ScalarType ONE = SCT::one ();
957 const MagnitudeType ZERO = SCT::magnitude (SCT::zero ());
959 const int xc = MVT::GetNumberVecs (X);
960 const ptrdiff_t xr = MVT::GetGlobalLength (X);
973 if (MX == Teuchos::null) {
975 MX = MVT::Clone(X,xc);
976 OPT::Apply(*(this->_Op),X,*MX);
983 if ( B == Teuchos::null ) {
987 const int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc;
988 const ptrdiff_t mxr = (this->_hasOp) ? MVT::GetGlobalLength( *MX ) : xr;
992 "Belos::ICGSOrthoManager::findBasis(): X must be non-empty" );
994 "Belos::ICGSOrthoManager::findBasis(): Size of X not consistant with size of B" );
996 "Belos::ICGSOrthoManager::findBasis(): Size of X not consistant with size of MX" );
998 "Belos::ICGSOrthoManager::findBasis(): Size of X not feasible for normalization" );
1000 "Belos::ICGSOrthoManager::findBasis(): Invalid howMany parameter" );
1005 int xstart = xc - howMany;
1007 for (
int j = xstart; j < xc; j++) {
1013 bool addVec =
false;
1016 std::vector<int> index(1);
1022 MXj = MVT::CloneViewNonConst( *MX, index );
1030 std::vector<int> prev_idx( numX );
1035 for (
int i=0; i<numX; i++) {
1038 prevX = MVT::CloneView( X, prev_idx );
1040 prevMX = MVT::CloneView( *MX, prev_idx );
1043 oldMXj = MVT::CloneCopy( *MXj );
1048 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
1053 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1056 MVT::MvDot( *Xj, *MXj, oldDot );
1060 "Belos::ICGSOrthoManager::findBasis(): Negative definiteness discovered in inner product" );
1066 for (
int i=0; i<max_ortho_steps_; ++i) {
1070 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1079 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1082 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj );
1090 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1093 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj );
1107 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1110 MVT::MvDot( *Xj, *oldMXj, newDot );
1113 newDot[0] = oldDot[0];
1117 if (completeBasis) {
1121 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(sing_tol_*oldDot[0]) ) {
1126 std::cout <<
"Belos::ICGSOrthoManager::findBasis() --> Random for column " << numX << std::endl;
1131 MVT::MvRandom( *tempXj );
1133 tempMXj = MVT::Clone( X, 1 );
1134 OPT::Apply( *(this->_Op), *tempXj, *tempMXj );
1140 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1143 MVT::MvDot( *tempXj, *tempMXj, oldDot );
1146 for (
int num_orth=0; num_orth<max_ortho_steps_; num_orth++){
1148 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1154 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1157 MVT::MvTimesMatAddMv( -ONE, *prevX, product, ONE, *tempXj );
1160 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1163 MVT::MvTimesMatAddMv( -ONE, *prevMX, product, ONE, *tempMXj );
1168 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1171 MVT::MvDot( *tempXj, *tempMXj, newDot );
1174 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) {
1176 MVT::Assign( *tempXj, *Xj );
1178 MVT::Assign( *tempMXj, *MXj );
1190 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*blk_tol_) ) {
1198 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1200 MVT::MvScale( *Xj, ONE/diag );
1203 MVT::MvScale( *MXj, ONE/diag );
1217 for (
int i=0; i<numX; i++) {
1218 (*B)(i,j) = product(i,0);
1229 template<
class ScalarType,
class MV,
class OP>
1231 ICGSOrthoManager<ScalarType, MV, OP>::blkOrtho1 ( MV &X,
Teuchos::RCP<MV> MX,
1236 int xc = MVT::GetNumberVecs( X );
1237 const ScalarType ONE = SCT::one();
1239 std::vector<int> qcs( nq );
1240 for (
int i=0; i<nq; i++) {
1241 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1248 for (
int i=0; i<nq; i++) {
1251 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1258 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1261 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X );
1267 OPT::Apply( *(this->_Op), X, *MX);
1271 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1272 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1274 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1277 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX );
1284 for (
int j = 1; j < max_ortho_steps_; ++j) {
1286 for (
int i=0; i<nq; i++) {
1291 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1298 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1301 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X );
1307 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1311 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX );
1313 else if (xc <= qcs[i]) {
1315 OPT::Apply( *(this->_Op), X, *MX);
1326 template<
class ScalarType,
class MV,
class OP>
1328 ICGSOrthoManager<ScalarType, MV, OP>::blkOrtho ( MV &X,
Teuchos::RCP<MV> MX,
1333 int xc = MVT::GetNumberVecs( X );
1334 bool dep_flg =
false;
1335 const ScalarType ONE = SCT::one();
1337 std::vector<int> qcs( nq );
1338 for (
int i=0; i<nq; i++) {
1339 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1345 std::vector<ScalarType> oldDot( xc );
1347 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1350 MVT::MvDot( X, *MX, oldDot );
1355 for (
int i=0; i<nq; i++) {
1358 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1365 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1368 MVT::MvTimesMatAddMv( -ONE, *Q[i], *C[i], ONE, X );
1373 OPT::Apply( *(this->_Op), X, *MX);
1377 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1378 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1380 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1383 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *C[i], ONE, *MX );
1390 for (
int j = 1; j < max_ortho_steps_; ++j) {
1392 for (
int i=0; i<nq; i++) {
1397 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1404 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1407 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, X );
1413 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1417 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX );
1419 else if (xc <= qcs[i]) {
1421 OPT::Apply( *(this->_Op), X, *MX);
1428 std::vector<ScalarType> newDot(xc);
1430 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1433 MVT::MvDot( X, *MX, newDot );
1437 for (
int i=0; i<xc; i++){
1438 if (SCT::magnitude(newDot[i]) < SCT::magnitude(oldDot[i] * blk_tol_)) {
1447 template<
class ScalarType,
class MV,
class OP>
1449 ICGSOrthoManager<ScalarType, MV, OP>::blkOrthoSing ( MV &X,
Teuchos::RCP<MV> MX,
1456 const ScalarType ONE = SCT::one();
1457 const ScalarType ZERO = SCT::zero();
1460 int xc = MVT::GetNumberVecs( X );
1461 std::vector<int> indX( 1 );
1462 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
1464 std::vector<int> qcs( nq );
1465 for (
int i=0; i<nq; i++) {
1466 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1475 for (
int j=0; j<xc; j++) {
1477 bool dep_flg =
false;
1481 std::vector<int> index( j );
1482 for (
int ind=0; ind<j; ind++) {
1485 lastQ = MVT::CloneView( X, index );
1488 Q.push_back( lastQ );
1490 qcs.push_back( MVT::GetNumberVecs( *lastQ ) );
1495 Xj = MVT::CloneViewNonConst( X, indX );
1497 MXj = MVT::CloneViewNonConst( *MX, indX );
1505 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1508 MVT::MvDot( *Xj, *MXj, oldDot );
1513 for (
int i=0; i<Q.size(); i++) {
1520 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1526 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1530 MVT::MvTimesMatAddMv( -ONE, *Q[i], tempC, ONE, *Xj );
1535 OPT::Apply( *(this->_Op), *Xj, *MXj);
1539 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1540 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1542 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1545 MVT::MvTimesMatAddMv( -ONE, *MQ[i], tempC, ONE, *MXj );
1552 for (
int num_ortho_steps=1; num_ortho_steps < max_ortho_steps_; ++num_ortho_steps) {
1554 for (
int i=0; i<Q.size(); i++) {
1560 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1567 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1570 MVT::MvTimesMatAddMv( -ONE, *Q[i], C2, ONE, *Xj );
1577 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1580 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MXj );
1582 else if (xc <= qcs[i]) {
1584 OPT::Apply( *(this->_Op), *Xj, *MXj);
1592 if (SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*sing_tol_)) {
1598 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1600 MVT::MvScale( *Xj, ONE/diag );
1603 MVT::MvScale( *MXj, ONE/diag );
1613 MVT::MvRandom( *tempXj );
1615 tempMXj = MVT::Clone( X, 1 );
1616 OPT::Apply( *(this->_Op), *tempXj, *tempMXj );
1622 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1625 MVT::MvDot( *tempXj, *tempMXj, oldDot );
1628 for (
int num_orth=0; num_orth<max_ortho_steps_; num_orth++) {
1630 for (
int i=0; i<Q.size(); i++) {
1635 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1641 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1644 MVT::MvTimesMatAddMv( -ONE, *Q[i], product, ONE, *tempXj );
1650 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1654 MVT::MvTimesMatAddMv( -ONE, *MQ[i], product, ONE, *tempMXj );
1656 else if (xc <= qcs[i]) {
1658 OPT::Apply( *(this->_Op), *tempXj, *tempMXj);
1667 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1670 MVT::MvDot( *tempXj, *tempMXj, newDot );
1674 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) {
1675 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1681 MVT::MvAddMv( ONE/diag, *tempXj, ZERO, *tempXj, *Xj );
1683 MVT::MvAddMv( ONE/diag, *tempMXj, ZERO, *tempMXj, *MXj );
1703 template<
class ScalarType,
class MV,
class OP>
1707 using Teuchos::parameterList;
1710 RCP<ParameterList> params = parameterList (
"ICGS");
1715 "Maximum number of orthogonalization passes (includes the "
1716 "first). Default is 2, since \"twice is enough\" for Krylov "
1719 "Block reorthogonalization threshold.");
1721 "Singular block detection threshold.");
1726 template<
class ScalarType,
class MV,
class OP>
1732 RCP<ParameterList> params = getICGSDefaultParameters<ScalarType, MV, OP>();
1734 params->set (
"maxNumOrthogPasses",
1736 params->set (
"blkTol",
1738 params->set (
"singTol",
1746 #endif // BELOS_ICGS_ORTHOMANAGER_HPP