15 #ifndef BELOS_IMGS_ORTHOMANAGER_HPP
16 #define BELOS_IMGS_ORTHOMANAGER_HPP
35 #ifdef BELOS_TEUCHOS_TIME_MONITOR
37 #endif // BELOS_TEUCHOS_TIME_MONITOR
42 template<
class ScalarType,
class MV,
class OP>
46 template<
class ScalarType,
class MV,
class OP>
49 template<
class ScalarType,
class MV,
class OP>
76 #ifdef BELOS_TEUCHOS_TIME_MONITOR
80 std::string orthoLabel = ss.str() +
": Orthogonalization";
83 std::string updateLabel = ss.str() +
": Ortho (Update)";
86 std::string normLabel = ss.str() +
": Ortho (Norm)";
89 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
96 const std::string& label =
"Belos",
106 #ifdef BELOS_TEUCHOS_TIME_MONITOR
107 std::stringstream ss;
110 std::string orthoLabel = ss.str() +
": Orthogonalization";
113 std::string updateLabel = ss.str() +
": Ortho (Update)";
116 std::string normLabel = ss.str() +
": Ortho (Norm)";
119 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
135 using Teuchos::parameterList;
139 RCP<ParameterList> params;
141 params = parameterList (*defaultParams);
156 int maxNumOrthogPasses;
161 maxNumOrthogPasses = params->
get<
int> (
"maxNumOrthogPasses");
162 }
catch (InvalidParameterName&) {
163 maxNumOrthogPasses = defaultParams->get<
int> (
"maxNumOrthogPasses");
164 params->set (
"maxNumOrthogPasses", maxNumOrthogPasses);
176 }
catch (InvalidParameterName&) {
181 params->remove (
"depTol");
182 }
catch (InvalidParameterName&) {
185 params->set (
"blkTol", blkTol);
190 }
catch (InvalidParameterName&) {
192 params->set (
"singTol", singTol);
206 defaultParams_ = Belos::getIMGSDefaultParameters<ScalarType, MV, OP>();
402 void setLabel(
const std::string& label);
440 #ifdef BELOS_TEUCHOS_TIME_MONITOR
442 #endif // BELOS_TEUCHOS_TIME_MONITOR
450 bool completeBasis,
int howMany = -1 )
const;
482 template<
class ScalarType,
class MV,
class OP>
485 template<
class ScalarType,
class MV,
class OP>
491 template<
class ScalarType,
class MV,
class OP>
496 template<
class ScalarType,
class MV,
class OP>
499 template<
class ScalarType,
class MV,
class OP>
504 template<
class ScalarType,
class MV,
class OP>
511 template<
class ScalarType,
class MV,
class OP>
514 if (label != label_) {
516 #ifdef BELOS_TEUCHOS_TIME_MONITOR
517 std::stringstream ss;
518 ss << label_ +
": IMGS[" << max_ortho_steps_ <<
"]";
520 std::string orthoLabel = ss.str() +
": Orthogonalization";
523 std::string updateLabel = ss.str() +
": Ortho (Update)";
526 std::string normLabel = ss.str() +
": Ortho (Norm)";
529 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
537 template<
class ScalarType,
class MV,
class OP>
540 const ScalarType ONE = SCT::one();
541 int rank = MVT::GetNumberVecs(X);
544 for (
int i=0; i<rank; i++) {
552 template<
class ScalarType,
class MV,
class OP>
555 int r1 = MVT::GetNumberVecs(X1);
556 int r2 = MVT::GetNumberVecs(X2);
564 template<
class ScalarType,
class MV,
class OP>
580 typedef typename Array< RCP< const MV > >::size_type size_type;
582 #ifdef BELOS_TEUCHOS_TIME_MONITOR
586 ScalarType ONE = SCT::one();
590 int xc = MVT::GetNumberVecs( X );
591 ptrdiff_t xr = MVT::GetGlobalLength( X );
598 B =
rcp (
new serial_dense_matrix_type (xc, xc));
608 for (size_type k = 0; k < nq; ++k)
610 const int numRows = MVT::GetNumberVecs (*Q[k]);
611 const int numCols = xc;
614 C[k] =
rcp (
new serial_dense_matrix_type (numRows, numCols));
615 else if (
C[k]->numRows() != numRows ||
C[k]->numCols() != numCols)
617 int err =
C[k]->reshape (numRows, numCols);
619 "IMGS orthogonalization: failed to reshape "
620 "C[" << k <<
"] (the array of block "
621 "coefficients resulting from projecting X "
622 "against Q[1:" << nq <<
"]).");
630 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
631 OPT::Apply(*(this->_Op),X,*MX);
639 int mxc = MVT::GetNumberVecs( *MX );
640 ptrdiff_t mxr = MVT::GetGlobalLength( *MX );
643 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
"Belos::IMGSOrthoManager::projectAndNormalize(): X must be non-empty" );
646 for (
int i=0; i<nq; i++) {
647 numbas += MVT::GetNumberVecs( *Q[i] );
652 "Belos::IMGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of B" );
655 "Belos::IMGSOrthoManager::projectAndNormalize(): MVT returned negative dimensions for X,MX" );
657 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument,
658 "Belos::IMGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of MX" );
664 bool dep_flg =
false;
668 tmpX = MVT::CloneCopy(X);
670 tmpMX = MVT::CloneCopy(*MX);
677 dep_flg = blkOrtho1( X, MX,
C, Q );
683 std::vector<ScalarType> diag(xc);
685 #ifdef BELOS_TEUCHOS_TIME_MONITOR
688 MVT::MvDot( X, *MX, diag );
690 (*B)(0,0) = SCT::squareroot(SCT::magnitude(diag[0]));
692 if (SCT::magnitude((*
B)(0,0)) > ZERO) {
694 MVT::MvScale( X, ONE/(*
B)(0,0) );
697 MVT::MvScale( *MX, ONE/(*
B)(0,0) );
704 dep_flg = blkOrtho( X, MX,
C, Q );
710 rank = blkOrthoSing( *tmpX, tmpMX,
C,
B, Q );
713 MVT::Assign( *tmpX, X );
715 MVT::Assign( *tmpMX, *MX );
720 rank = findBasis( X, MX,
B,
false );
725 rank = blkOrthoSing( *tmpX, tmpMX,
C,
B, Q );
728 MVT::Assign( *tmpX, X );
730 MVT::Assign( *tmpMX, *MX );
737 TEUCHOS_TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error,
738 "Belos::IMGSOrthoManager::projectAndNormalize(): Debug error in rank variable." );
748 template<
class ScalarType,
class MV,
class OP>
753 #ifdef BELOS_TEUCHOS_TIME_MONITOR
758 return findBasis(X, MX,
B,
true);
764 template<
class ScalarType,
class MV,
class OP>
784 #ifdef BELOS_TEUCHOS_TIME_MONITOR
788 int xc = MVT::GetNumberVecs( X );
789 ptrdiff_t xr = MVT::GetGlobalLength( X );
791 std::vector<int> qcs(nq);
793 if (nq == 0 || xc == 0 || xr == 0) {
796 ptrdiff_t qr = MVT::GetGlobalLength ( *Q[0] );
807 MX = MVT::Clone(X,MVT::GetNumberVecs(X));
808 OPT::Apply(*(this->_Op),X,*MX);
815 int mxc = MVT::GetNumberVecs( *MX );
816 ptrdiff_t mxr = MVT::GetGlobalLength( *MX );
820 "Belos::IMGSOrthoManager::project(): MVT returned negative dimensions for X,MX" );
823 "Belos::IMGSOrthoManager::project(): Size of X not consistant with MX,Q" );
826 for (
int i=0; i<nq; i++) {
828 "Belos::IMGSOrthoManager::project(): Q lengths not mutually consistant" );
829 qcs[i] = MVT::GetNumberVecs( *Q[i] );
831 "Belos::IMGSOrthoManager::project(): Q has less rows than columns" );
839 "Belos::IMGSOrthoManager::project(): Size of Q not consistant with size of C" );
844 blkOrtho( X, MX,
C, Q );
851 template<
class ScalarType,
class MV,
class OP>
855 bool completeBasis,
int howMany )
const {
872 const ScalarType ONE = SCT::one();
875 int xc = MVT::GetNumberVecs( X );
876 ptrdiff_t xr = MVT::GetGlobalLength( X );
891 MX = MVT::Clone(X,xc);
892 OPT::Apply(*(this->_Op),X,*MX);
903 int mxc = (this->_hasOp) ? MVT::GetNumberVecs( *MX ) : xc;
904 ptrdiff_t mxr = (this->_hasOp) ? MVT::GetGlobalLength( *MX ) : xr;
908 "Belos::IMGSOrthoManager::findBasis(): X must be non-empty" );
910 "Belos::IMGSOrthoManager::findBasis(): Size of X not consistant with size of B" );
912 "Belos::IMGSOrthoManager::findBasis(): Size of X not consistant with size of MX" );
913 TEUCHOS_TEST_FOR_EXCEPTION( xc > xr, std::invalid_argument,
914 "Belos::IMGSOrthoManager::findBasis(): Size of X not feasible for normalization" );
915 TEUCHOS_TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::invalid_argument,
916 "Belos::IMGSOrthoManager::findBasis(): Invalid howMany parameter" );
921 int xstart = xc - howMany;
923 for (
int j = xstart; j < xc; j++) {
932 std::vector<int> index(1);
936 if ((this->_hasOp)) {
938 MXj = MVT::CloneViewNonConst( *MX, index );
947 oldMXj = MVT::CloneCopy( *MXj );
955 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
960 #ifdef BELOS_TEUCHOS_TIME_MONITOR
963 MVT::MvDot( *Xj, *MXj, oldDot );
966 TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(oldDot[0]) < ZERO,
OrthoError,
967 "Belos::IMGSOrthoManager::findBasis(): Negative definiteness discovered in inner product" );
970 for (
int ii=0; ii<numX; ii++) {
973 prevX = MVT::CloneView( X, index );
975 prevMX = MVT::CloneView( *MX, index );
978 for (
int i=0; i<max_ortho_steps_; ++i) {
982 #ifdef BELOS_TEUCHOS_TIME_MONITOR
991 #ifdef BELOS_TEUCHOS_TIME_MONITOR
994 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *Xj );
1002 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1005 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *MXj );
1010 product[ii] = P2[0];
1012 product[ii] += P2[0];
1020 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1023 MVT::MvDot( *Xj, *oldMXj, newDot );
1026 newDot[0] = oldDot[0];
1030 if (completeBasis) {
1034 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(sing_tol_*oldDot[0]) ) {
1039 std::cout <<
"Belos::IMGSOrthoManager::findBasis() --> Random for column " << numX << std::endl;
1044 MVT::MvRandom( *tempXj );
1046 tempMXj = MVT::Clone( X, 1 );
1047 OPT::Apply( *(this->_Op), *tempXj, *tempMXj );
1053 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1056 MVT::MvDot( *tempXj, *tempMXj, oldDot );
1060 for (
int ii=0; ii<numX; ii++) {
1063 prevX = MVT::CloneView( X, index );
1065 prevMX = MVT::CloneView( *MX, index );
1068 for (
int num_orth=0; num_orth<max_ortho_steps_; num_orth++){
1070 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1076 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1079 MVT::MvTimesMatAddMv( -ONE, *prevX, P2, ONE, *tempXj );
1082 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1085 MVT::MvTimesMatAddMv( -ONE, *prevMX, P2, ONE, *tempMXj );
1090 product[ii] = P2[0];
1092 product[ii] += P2[0];
1098 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1101 MVT::MvDot( *tempXj, *tempMXj, newDot );
1104 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) {
1106 MVT::Assign( *tempXj, *Xj );
1108 MVT::Assign( *tempMXj, *MXj );
1120 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*blk_tol_) ) {
1129 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1130 if (SCT::magnitude(diag) > ZERO) {
1131 MVT::MvScale( *Xj, ONE/diag );
1134 MVT::MvScale( *MXj, ONE/diag );
1148 for (
int i=0; i<numX; i++) {
1149 (*B)(i,j) = product(i);
1160 template<
class ScalarType,
class MV,
class OP>
1167 int xc = MVT::GetNumberVecs( X );
1168 const ScalarType ONE = SCT::one();
1170 std::vector<int> qcs( nq );
1171 for (
int i=0; i<nq; i++) {
1172 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1176 std::vector<int> index(1);
1181 for (
int i=0; i<nq; i++) {
1184 for (
int ii=0; ii<qcs[i]; ii++) {
1187 tempQ = MVT::CloneView( *Q[i], index );
1192 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1199 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1202 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC, ONE, X );
1209 OPT::Apply( *(this->_Op), X, *MX);
1213 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1214 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1215 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *
C[i], ONE, *MX );
1221 for (
int j = 1; j < max_ortho_steps_; ++j) {
1223 for (
int i=0; i<nq; i++) {
1228 for (
int ii=0; ii<qcs[i]; ii++) {
1231 tempQ = MVT::CloneView( *Q[i], index );
1237 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1244 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1247 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC2, ONE, X );
1256 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX );
1258 else if (xc <= qcs[i]) {
1260 OPT::Apply( *(this->_Op), X, *MX);
1271 template<
class ScalarType,
class MV,
class OP>
1278 int xc = MVT::GetNumberVecs( X );
1279 bool dep_flg =
false;
1280 const ScalarType ONE = SCT::one();
1282 std::vector<int> qcs( nq );
1283 for (
int i=0; i<nq; i++) {
1284 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1290 std::vector<ScalarType> oldDot( xc );
1292 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1295 MVT::MvDot( X, *MX, oldDot );
1298 std::vector<int> index(1);
1303 for (
int i=0; i<nq; i++) {
1306 for (
int ii=0; ii<qcs[i]; ii++) {
1309 tempQ = MVT::CloneView( *Q[i], index );
1314 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1321 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1324 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC, ONE, X );
1331 OPT::Apply( *(this->_Op), X, *MX);
1335 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1336 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1337 MVT::MvTimesMatAddMv( -ONE, *MQ[i], *
C[i], ONE, *MX );
1343 for (
int j = 1; j < max_ortho_steps_; ++j) {
1345 for (
int i=0; i<nq; i++) {
1349 for (
int ii=0; ii<qcs[i]; ii++) {
1352 tempQ = MVT::CloneView( *Q[i], index );
1358 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1365 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1368 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC2, ONE, X );
1376 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MX );
1378 else if (xc <= qcs[i]) {
1380 OPT::Apply( *(this->_Op), X, *MX);
1387 std::vector<ScalarType> newDot(xc);
1389 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1392 MVT::MvDot( X, *MX, newDot );
1396 for (
int i=0; i<xc; i++){
1397 if (SCT::magnitude(newDot[i]) < SCT::magnitude(oldDot[i] * blk_tol_)) {
1406 template<
class ScalarType,
class MV,
class OP>
1415 const ScalarType ONE = SCT::one();
1416 const ScalarType ZERO = SCT::zero();
1419 int xc = MVT::GetNumberVecs( X );
1420 std::vector<int> indX( 1 );
1421 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
1423 std::vector<int> qcs( nq );
1424 for (
int i=0; i<nq; i++) {
1425 qcs[i] = MVT::GetNumberVecs( *Q[i] );
1434 for (
int j=0; j<xc; j++) {
1436 bool dep_flg =
false;
1440 std::vector<int> index( j );
1441 for (
int ind=0; ind<j; ind++) {
1444 lastQ = MVT::CloneView( X, index );
1449 qcs.push_back( MVT::GetNumberVecs( *lastQ ) );
1454 Xj = MVT::CloneViewNonConst( X, indX );
1456 MXj = MVT::CloneViewNonConst( *MX, indX );
1464 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1467 MVT::MvDot( *Xj, *MXj, oldDot );
1474 for (
int i=0; i<Q.
size(); i++) {
1477 for (
int ii=0; ii<qcs[i]; ii++) {
1480 tempQ = MVT::CloneView( *Q[i], indX );
1488 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC, ONE, *Xj );
1494 OPT::Apply( *(this->_Op), *Xj, *MXj);
1498 MQ[i] = MVT::Clone( *Q[i], qcs[i] );
1499 OPT::Apply( *(this->_Op), *Q[i], *MQ[i] );
1501 MVT::MvTimesMatAddMv( -ONE, *MQ[i], tempC, ONE, *MXj );
1507 for (
int num_ortho_steps=1; num_ortho_steps < max_ortho_steps_; ++num_ortho_steps) {
1509 for (
int i=0; i<Q.
size(); i++) {
1513 for (
int ii=0; ii<qcs[i]; ii++) {
1516 tempQ = MVT::CloneView( *Q[i], indX );
1522 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC2, ONE, *Xj );
1533 MVT::MvTimesMatAddMv( -ONE, *MQ[i], C2, ONE, *MXj );
1535 else if (xc <= qcs[i]) {
1537 OPT::Apply( *(this->_Op), *Xj, *MXj);
1546 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1549 MVT::MvDot( *Xj, *MXj, newDot );
1553 if (SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*sing_tol_)) {
1559 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1561 MVT::MvScale( *Xj, ONE/diag );
1564 MVT::MvScale( *MXj, ONE/diag );
1574 MVT::MvRandom( *tempXj );
1576 tempMXj = MVT::Clone( X, 1 );
1577 OPT::Apply( *(this->_Op), *tempXj, *tempMXj );
1583 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1586 MVT::MvDot( *tempXj, *tempMXj, oldDot );
1589 for (
int num_orth=0; num_orth<max_ortho_steps_; num_orth++) {
1591 for (
int i=0; i<Q.
size(); i++) {
1595 for (
int ii=0; ii<qcs[i]; ii++) {
1598 tempQ = MVT::CloneView( *Q[i], indX );
1603 MVT::MvTimesMatAddMv( -ONE, *tempQ, tempC, ONE, *tempXj );
1611 MVT::MvTimesMatAddMv( -ONE, *MQ[i], product, ONE, *tempMXj );
1613 else if (xc <= qcs[i]) {
1615 OPT::Apply( *(this->_Op), *tempXj, *tempMXj);
1623 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1626 MVT::MvDot( *tempXj, *tempMXj, newDot );
1630 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*sing_tol_) ) {
1631 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1637 MVT::MvAddMv( ONE/diag, *tempXj, ZERO, *tempXj, *Xj );
1639 MVT::MvAddMv( ONE/diag, *tempMXj, ZERO, *tempMXj, *MXj );
1659 template<
class ScalarType,
class MV,
class OP>
1663 using Teuchos::parameterList;
1666 RCP<ParameterList> params = parameterList (
"IMGS");
1671 "Maximum number of orthogonalization passes (includes the "
1672 "first). Default is 2, since \"twice is enough\" for Krylov "
1675 "Block reorthogonalization threshold.");
1677 "Singular block detection threshold.");
1682 template<
class ScalarType,
class MV,
class OP>
1688 RCP<ParameterList> params = getIMGSDefaultParameters<ScalarType, MV, OP>();
1690 params->set (
"maxNumOrthogPasses",
1692 params->set (
"blkTol",
1694 params->set (
"singTol",
1702 #endif // BELOS_IMGS_ORTHOMANAGER_HPP
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &plist)
bool is_null(const boost::shared_ptr< T > &p)
Teuchos::ScalarTraits< ScalarType >::magnitudeType orthogError(const MV &X1, const MV &X2) const
This method computes the error in orthogonality of two multivectors, measured as the Frobenius norm o...
MagnitudeType blk_tol_
Block reorthogonalization tolerance.
static const int max_ortho_steps_fast_
Max number of (re)orthogonalization steps, including the first (fast).
T & get(ParameterList &l, const std::string &name)
void setBlkTol(const MagnitudeType blk_tol)
Set parameter for block re-orthogonalization threshhold.
int findBasis(MV &X, Teuchos::RCP< MV > MX, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > C, bool completeBasis, int howMany=-1) const
Routine to find an orthonormal basis for X.
void project(MV &X, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
This method calls project(X,Teuchos::null,C,Q); see documentation for that function.
bool is_null(const std::shared_ptr< T > &p)
Teuchos::ScalarTraits< ScalarType >::magnitudeType orthonormError(const MV &X) const
This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of ...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool blkOrtho1(MV &X, Teuchos::RCP< MV > MX, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
Routine to compute the block orthogonalization.
Declaration of basic traits for the multivector type.
int blkOrthoSing(MV &X, Teuchos::RCP< MV > MX, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B, Teuchos::ArrayView< Teuchos::RCP< const MV > > QQ) const
Project X against QQ and normalize X, one vector at a time.
Class which defines basic traits for the operator type.
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
const std::string & getLabel() const
This method returns the label being used by the timers in the orthogonalization manager.
Traits class which defines basic operations on multivectors.
void innerProd(const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z) const
Provides the inner product defining the orthogonality concepts, using the provided operator...
IMGSOrthoManager(const std::string &label="Belos", Teuchos::RCP< const OP > Op=Teuchos::null, const int max_ortho_steps=max_ortho_steps_default_, const MagnitudeType blk_tol=blk_tol_default_, const MagnitudeType sing_tol=sing_tol_default_)
Constructor specifying re-orthogonalization tolerance.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
std::string label_
Label for timers.
bool is_null(const RCP< T > &p)
int normalize(MV &X, Teuchos::RCP< MV > MX, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B) const
This method takes a multivector X and attempts to compute an orthonormal basis for ...
MagnitudeType getBlkTol() const
Return parameter for block re-orthogonalization threshhold.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::ScalarTraits< MagnitudeType > MGT
Exception thrown to signal error in an orthogonalization manager method.
MagnitudeType sing_tol_
Singular block detection threshold.
void setMyParamList(const RCP< ParameterList > ¶mList)
int max_ortho_steps_
Max number of (re)orthogonalization steps, including the first.
MultiVecTraits< ScalarType, MV > MVT
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
void resize(size_type new_size, const value_type &x=value_type())
void setLabel(const std::string &label)
This method sets the label used by the timers in the orthogonalization manager.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
void push_back(const value_type &x)
void project(MV &X, Teuchos::RCP< MV > MX, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
Given a list of (mutually and internally) orthonormal bases Q, this method takes a multivector X and ...
static const MagnitudeType sing_tol_default_
Singular block detection threshold (default).
Teuchos::RCP< Teuchos::ParameterList > getIMGSDefaultParameters()
"Default" parameters for robustness and accuracy.
Teuchos::ScalarTraits< ScalarType > SCT
Teuchos::RCP< Teuchos::ParameterList > defaultParams_
Default parameter list.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
static const int max_ortho_steps_default_
Max number of (re)orthogonalization steps, including the first (default).
IMGSOrthoManager(const Teuchos::RCP< Teuchos::ParameterList > &plist, const std::string &label="Belos", Teuchos::RCP< const OP > Op=Teuchos::null)
Constructor that takes a list of parameters.
MagnitudeType getSingTol() const
Return parameter for singular block detection.
Class which defines basic traits for the operator type.
int normalize(MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B) const
This method calls normalize(X,Teuchos::null,B); see documentation for that function.
virtual int projectAndNormalizeWithMxImpl(MV &X, Teuchos::RCP< MV > MX, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
Given a set of bases Q[i] and a multivector X, this method computes an orthonormal basis for ...
~IMGSOrthoManager()
Destructor.
Belos's templated virtual class for providing routines for orthogonalization and orthonormzalition of...
void setSingTol(const MagnitudeType sing_tol)
Set parameter for singular block detection.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< Teuchos::ParameterList > getIMGSFastParameters()
"Fast" but possibly unsafe or less accurate parameters.
static const MagnitudeType blk_tol_default_
Block reorthogonalization threshold (default).
static const MagnitudeType sing_tol_fast_
Singular block detection threshold (fast).
bool blkOrtho(MV &X, Teuchos::RCP< MV > MX, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
Routine to compute the block orthogonalization.
static const MagnitudeType blk_tol_fast_
Block reorthogonalization threshold (fast).
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...