47#ifndef BELOS_IMGS_ORTHOMANAGER_HPP
48#define BELOS_IMGS_ORTHOMANAGER_HPP
63#include "Teuchos_SerialDenseMatrix.hpp"
64#include "Teuchos_SerialDenseVector.hpp"
66#include "Teuchos_as.hpp"
67#ifdef BELOS_TEUCHOS_TIME_MONITOR
68#include "Teuchos_TimeMonitor.hpp"
74 template<
class ScalarType,
class MV,
class OP>
78 template<
class ScalarType,
class MV,
class OP>
81 template<
class ScalarType,
class MV,
class OP>
86 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
MagnitudeType;
87 typedef typename Teuchos::ScalarTraits<MagnitudeType>
MGT;
88 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
98 Teuchos::RCP<const OP> Op = Teuchos::null,
108#ifdef BELOS_TEUCHOS_TIME_MONITOR
109 std::stringstream ss;
112 std::string orthoLabel = ss.str() +
": Orthogonalization";
113 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel);
115 std::string updateLabel = ss.str() +
": Ortho (Update)";
116 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel);
118 std::string normLabel = ss.str() +
": Ortho (Norm)";
119 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel);
121 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
122 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel);
128 const std::string& label =
"Belos",
129 Teuchos::RCP<const OP> Op = Teuchos::null) :
138#ifdef BELOS_TEUCHOS_TIME_MONITOR
139 std::stringstream ss;
142 std::string orthoLabel = ss.str() +
": Orthogonalization";
143 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel);
145 std::string updateLabel = ss.str() +
": Ortho (Update)";
146 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel);
148 std::string normLabel = ss.str() +
": Ortho (Norm)";
149 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel);
151 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
152 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel);
165 using Teuchos::Exceptions::InvalidParameterName;
166 using Teuchos::ParameterList;
167 using Teuchos::parameterList;
171 RCP<ParameterList> params;
172 if (plist.is_null()) {
173 params = parameterList (*defaultParams);
188 int maxNumOrthogPasses;
193 maxNumOrthogPasses = params->get<
int> (
"maxNumOrthogPasses");
194 }
catch (InvalidParameterName&) {
195 maxNumOrthogPasses = defaultParams->get<
int> (
"maxNumOrthogPasses");
196 params->set (
"maxNumOrthogPasses", maxNumOrthogPasses);
208 }
catch (InvalidParameterName&) {
213 params->remove (
"depTol");
214 }
catch (InvalidParameterName&) {
217 params->set (
"blkTol", blkTol);
222 }
catch (InvalidParameterName&) {
224 params->set (
"singTol", singTol);
231 this->setMyParamList (params);
234 Teuchos::RCP<const Teuchos::ParameterList>
294 void project ( MV &X, Teuchos::RCP<MV> MX,
295 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
296 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const;
302 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
303 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const {
333 int normalize ( MV &X, Teuchos::RCP<MV> MX,
334 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B)
const;
339 int normalize ( MV &X, Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B )
const {
388 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
389 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
390 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const;
400 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
409 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
415 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
424 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
425 orthogError(
const MV &X1, Teuchos::RCP<const MV> MX1,
const MV &X2)
const;
434 void setLabel(
const std::string& label);
472#ifdef BELOS_TEUCHOS_TIME_MONITOR
473 Teuchos::RCP<Teuchos::Time> timerOrtho_, timerUpdate_, timerNorm_, timerInnerProd_;
480 int findBasis(MV &X, Teuchos::RCP<MV> MX,
481 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > C,
482 bool completeBasis,
int howMany = -1 )
const;
485 bool blkOrtho1 ( MV &X, Teuchos::RCP<MV> MX,
486 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
487 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const;
490 bool blkOrtho ( MV &X, Teuchos::RCP<MV> MX,
491 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
492 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const;
508 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
509 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
510 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ)
const;
514 template<
class ScalarType,
class MV,
class OP>
517 template<
class ScalarType,
class MV,
class OP>
520 = 10*Teuchos::ScalarTraits<typename IMGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::squareroot(
523 template<
class ScalarType,
class MV,
class OP>
526 = 10*Teuchos::ScalarTraits<typename IMGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::eps();
528 template<
class ScalarType,
class MV,
class OP>
531 template<
class ScalarType,
class MV,
class OP>
534 = Teuchos::ScalarTraits<typename IMGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero();
536 template<
class ScalarType,
class MV,
class OP>
539 = Teuchos::ScalarTraits<typename IMGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero();
543 template<
class ScalarType,
class MV,
class OP>
548#ifdef BELOS_TEUCHOS_TIME_MONITOR
549 std::stringstream ss;
552 std::string orthoLabel = ss.str() +
": Orthogonalization";
553 timerOrtho_ = Teuchos::TimeMonitor::getNewCounter(orthoLabel);
555 std::string updateLabel = ss.str() +
": Ortho (Update)";
556 timerUpdate_ = Teuchos::TimeMonitor::getNewCounter(updateLabel);
558 std::string normLabel = ss.str() +
": Ortho (Norm)";
559 timerNorm_ = Teuchos::TimeMonitor::getNewCounter(normLabel);
561 std::string ipLabel = ss.str() +
": Ortho (Inner Product)";
562 timerInnerProd_ = Teuchos::TimeMonitor::getNewCounter(ipLabel);
569 template<
class ScalarType,
class MV,
class OP>
570 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
572 const ScalarType ONE = SCT::one();
574 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(rank,rank);
576 for (
int i=0; i<rank; i++) {
579 return xTx.normFrobenius();
584 template<
class ScalarType,
class MV,
class OP>
585 typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
589 Teuchos::SerialDenseMatrix<int,ScalarType> xTx(r2,r1);
591 return xTx.normFrobenius();
596 template<
class ScalarType,
class MV,
class OP>
601 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
602 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
603 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const
605 using Teuchos::Array;
607 using Teuchos::is_null;
610 using Teuchos::SerialDenseMatrix;
611 typedef SerialDenseMatrix< int, ScalarType > serial_dense_matrix_type;
612 typedef typename Array< RCP< const MV > >::size_type size_type;
614#ifdef BELOS_TEUCHOS_TIME_MONITOR
615 Teuchos::TimeMonitor orthotimer(*timerOrtho_);
618 ScalarType ONE = SCT::one();
630 B = rcp (
new serial_dense_matrix_type (xc, xc));
640 for (size_type k = 0; k < nq; ++k)
643 const int numCols = xc;
646 C[k] = rcp (
new serial_dense_matrix_type (numRows, numCols));
647 else if (C[k]->numRows() != numRows || C[k]->numCols() != numCols)
649 int err = C[k]->reshape (numRows, numCols);
650 TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error,
651 "IMGS orthogonalization: failed to reshape "
652 "C[" << k <<
"] (the array of block "
653 "coefficients resulting from projecting X "
654 "against Q[1:" << nq <<
"]).");
660 if (MX == Teuchos::null) {
668 MX = Teuchos::rcp( &X,
false );
675 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
"Belos::IMGSOrthoManager::projectAndNormalize(): X must be non-empty" );
678 for (
int i=0; i<nq; i++) {
683 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
684 "Belos::IMGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of B" );
686 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
687 "Belos::IMGSOrthoManager::projectAndNormalize(): MVT returned negative dimensions for X,MX" );
689 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr, std::invalid_argument,
690 "Belos::IMGSOrthoManager::projectAndNormalize(): Size of X must be consistant with size of MX" );
696 bool dep_flg =
false;
699 Teuchos::RCP<MV> tmpX, tmpMX;
712 if ( B == Teuchos::null ) {
713 B = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
715 std::vector<ScalarType> diag(xc);
717#ifdef BELOS_TEUCHOS_TIME_MONITOR
718 Teuchos::TimeMonitor normTimer( *timerNorm_ );
722 (*B)(0,0) = SCT::squareroot(SCT::magnitude(diag[0]));
724 if (SCT::magnitude((*B)(0,0)) > ZERO) {
769 TEUCHOS_TEST_FOR_EXCEPTION( rank > xc || rank < 0, std::logic_error,
770 "Belos::IMGSOrthoManager::projectAndNormalize(): Debug error in rank variable." );
780 template<
class ScalarType,
class MV,
class OP>
782 MV &X, Teuchos::RCP<MV> MX,
783 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B )
const {
785#ifdef BELOS_TEUCHOS_TIME_MONITOR
786 Teuchos::TimeMonitor orthotimer(*timerOrtho_);
796 template<
class ScalarType,
class MV,
class OP>
798 MV &X, Teuchos::RCP<MV> MX,
799 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
800 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const {
816#ifdef BELOS_TEUCHOS_TIME_MONITOR
817 Teuchos::TimeMonitor orthotimer(*timerOrtho_);
823 std::vector<int> qcs(nq);
825 if (nq == 0 || xc == 0 || xr == 0) {
837 if (MX == Teuchos::null) {
845 MX = Teuchos::rcp( &X,
false );
851 TEUCHOS_TEST_FOR_EXCEPTION( xc<0 || xr<0 || mxc<0 || mxr<0, std::invalid_argument,
852 "Belos::IMGSOrthoManager::project(): MVT returned negative dimensions for X,MX" );
854 TEUCHOS_TEST_FOR_EXCEPTION( xc!=mxc || xr!=mxr || xr!=qr, std::invalid_argument,
855 "Belos::IMGSOrthoManager::project(): Size of X not consistant with MX,Q" );
858 for (
int i=0; i<nq; i++) {
860 "Belos::IMGSOrthoManager::project(): Q lengths not mutually consistant" );
862 TEUCHOS_TEST_FOR_EXCEPTION( qr < qcs[i], std::invalid_argument,
863 "Belos::IMGSOrthoManager::project(): Q has less rows than columns" );
866 if ( C[i] == Teuchos::null ) {
867 C[i] = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(qcs[i],xc) );
870 TEUCHOS_TEST_FOR_EXCEPTION( C[i]->numRows() != qcs[i] || C[i]->numCols() != xc , std::invalid_argument,
871 "Belos::IMGSOrthoManager::project(): Size of Q not consistant with size of C" );
883 template<
class ScalarType,
class MV,
class OP>
885 MV &X, Teuchos::RCP<MV> MX,
886 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
887 bool completeBasis,
int howMany )
const {
904 const ScalarType ONE = SCT::one();
921 if (MX == Teuchos::null) {
931 if ( B == Teuchos::null ) {
932 B = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<int,ScalarType>(xc,xc) );
939 TEUCHOS_TEST_FOR_EXCEPTION( xc == 0 || xr == 0, std::invalid_argument,
940 "Belos::IMGSOrthoManager::findBasis(): X must be non-empty" );
941 TEUCHOS_TEST_FOR_EXCEPTION( B->numRows() != xc || B->numCols() != xc, std::invalid_argument,
942 "Belos::IMGSOrthoManager::findBasis(): Size of X not consistant with size of B" );
943 TEUCHOS_TEST_FOR_EXCEPTION( xc != mxc || xr != mxr, std::invalid_argument,
944 "Belos::IMGSOrthoManager::findBasis(): Size of X not consistant with size of MX" );
945 TEUCHOS_TEST_FOR_EXCEPTION( xc > xr, std::invalid_argument,
946 "Belos::IMGSOrthoManager::findBasis(): Size of X not feasible for normalization" );
947 TEUCHOS_TEST_FOR_EXCEPTION( howMany < 0 || howMany > xc, std::invalid_argument,
948 "Belos::IMGSOrthoManager::findBasis(): Invalid howMany parameter" );
953 int xstart = xc - howMany;
955 for (
int j = xstart; j < xc; j++) {
964 std::vector<int> index(1);
967 Teuchos::RCP<MV> MXj;
977 Teuchos::RCP<MV> oldMXj;
983 Teuchos::SerialDenseVector<int,ScalarType> product(numX);
984 Teuchos::SerialDenseVector<int,ScalarType> P2(1);
985 Teuchos::RCP<const MV> prevX, prevMX;
987 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
992#ifdef BELOS_TEUCHOS_TIME_MONITOR
993 Teuchos::TimeMonitor normTimer( *timerNorm_ );
998 TEUCHOS_TEST_FOR_EXCEPTION( SCT::real(oldDot[0]) < ZERO,
OrthoError,
999 "Belos::IMGSOrthoManager::findBasis(): Negative definiteness discovered in inner product" );
1002 for (
int ii=0; ii<numX; ii++) {
1014#ifdef BELOS_TEUCHOS_TIME_MONITOR
1015 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1023#ifdef BELOS_TEUCHOS_TIME_MONITOR
1024 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1034#ifdef BELOS_TEUCHOS_TIME_MONITOR
1035 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1042 product[ii] = P2[0];
1044 product[ii] += P2[0];
1052#ifdef BELOS_TEUCHOS_TIME_MONITOR
1053 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1058 newDot[0] = oldDot[0];
1062 if (completeBasis) {
1066 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(
sing_tol_*oldDot[0]) ) {
1071 std::cout <<
"Belos::IMGSOrthoManager::findBasis() --> Random for column " << numX << std::endl;
1074 Teuchos::RCP<MV> tempXj =
MVT::Clone( X, 1 );
1075 Teuchos::RCP<MV> tempMXj;
1085#ifdef BELOS_TEUCHOS_TIME_MONITOR
1086 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1092 for (
int ii=0; ii<numX; ii++) {
1102#ifdef BELOS_TEUCHOS_TIME_MONITOR
1103 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1108#ifdef BELOS_TEUCHOS_TIME_MONITOR
1109 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1114#ifdef BELOS_TEUCHOS_TIME_MONITOR
1115 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1122 product[ii] = P2[0];
1124 product[ii] += P2[0];
1130#ifdef BELOS_TEUCHOS_TIME_MONITOR
1131 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1136 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*
sing_tol_) ) {
1152 if ( SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*
blk_tol_) ) {
1161 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1162 if (SCT::magnitude(diag) > ZERO) {
1180 for (
int i=0; i<numX; i++) {
1181 (*B)(i,j) = product(i);
1192 template<
class ScalarType,
class MV,
class OP>
1195 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
1196 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const
1200 const ScalarType ONE = SCT::one();
1202 std::vector<int> qcs( nq );
1203 for (
int i=0; i<nq; i++) {
1208 std::vector<int> index(1);
1209 Teuchos::RCP<const MV> tempQ;
1211 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq);
1213 for (
int i=0; i<nq; i++) {
1216 for (
int ii=0; ii<qcs[i]; ii++) {
1220 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], 1, 1, ii, 0 );
1224#ifdef BELOS_TEUCHOS_TIME_MONITOR
1225 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1231#ifdef BELOS_TEUCHOS_TIME_MONITOR
1232 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1255 for (
int i=0; i<nq; i++) {
1257 Teuchos::SerialDenseMatrix<int,ScalarType> C2(qcs[i],1);
1260 for (
int ii=0; ii<qcs[i]; ii++) {
1264 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], 1, 1, ii, 0 );
1265 Teuchos::SerialDenseMatrix<int,ScalarType> tempC2( Teuchos::View, C2, 1, 1, ii, 0 );
1269#ifdef BELOS_TEUCHOS_TIME_MONITOR
1270 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1276#ifdef BELOS_TEUCHOS_TIME_MONITOR
1277 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1290 else if (xc <= qcs[i]) {
1303 template<
class ScalarType,
class MV,
class OP>
1306 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
1307 Teuchos::ArrayView<Teuchos::RCP<const MV> > Q)
const
1311 bool dep_flg =
false;
1312 const ScalarType ONE = SCT::one();
1314 std::vector<int> qcs( nq );
1315 for (
int i=0; i<nq; i++) {
1322 std::vector<ScalarType> oldDot( xc );
1324#ifdef BELOS_TEUCHOS_TIME_MONITOR
1325 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1330 std::vector<int> index(1);
1331 Teuchos::Array<Teuchos::RCP<MV> > MQ(nq);
1332 Teuchos::RCP<const MV> tempQ;
1335 for (
int i=0; i<nq; i++) {
1338 for (
int ii=0; ii<qcs[i]; ii++) {
1342 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], 1, xc, ii, 0 );
1346#ifdef BELOS_TEUCHOS_TIME_MONITOR
1347 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1353#ifdef BELOS_TEUCHOS_TIME_MONITOR
1354 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1377 for (
int i=0; i<nq; i++) {
1378 Teuchos::SerialDenseMatrix<int,ScalarType> C2(qcs[i],xc);
1381 for (
int ii=0; ii<qcs[i]; ii++) {
1385 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], 1, xc, ii, 0 );
1386 Teuchos::SerialDenseMatrix<int,ScalarType> tempC2( Teuchos::View, C2, 1, xc, ii, 0 );
1390#ifdef BELOS_TEUCHOS_TIME_MONITOR
1391 Teuchos::TimeMonitor innerProdTimer( *timerInnerProd_ );
1397#ifdef BELOS_TEUCHOS_TIME_MONITOR
1398 Teuchos::TimeMonitor updateTimer( *timerUpdate_ );
1410 else if (xc <= qcs[i]) {
1419 std::vector<ScalarType> newDot(xc);
1421#ifdef BELOS_TEUCHOS_TIME_MONITOR
1422 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1428 for (
int i=0; i<xc; i++){
1429 if (SCT::magnitude(newDot[i]) < SCT::magnitude(oldDot[i] *
blk_tol_)) {
1438 template<
class ScalarType,
class MV,
class OP>
1441 Teuchos::Array<Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > > C,
1442 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B,
1443 Teuchos::ArrayView<Teuchos::RCP<const MV> > QQ)
const
1445 Teuchos::Array<Teuchos::RCP<const MV> > Q (QQ);
1447 const ScalarType ONE = SCT::one();
1448 const ScalarType ZERO = SCT::zero();
1452 std::vector<int> indX( 1 );
1453 std::vector<ScalarType> oldDot( 1 ), newDot( 1 );
1455 std::vector<int> qcs( nq );
1456 for (
int i=0; i<nq; i++) {
1461 Teuchos::RCP<const MV> lastQ;
1462 Teuchos::RCP<MV> Xj, MXj;
1463 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > lastC;
1466 for (
int j=0; j<xc; j++) {
1468 bool dep_flg =
false;
1472 std::vector<int> index( j );
1473 for (
int ind=0; ind<j; ind++) {
1479 Q.push_back( lastQ );
1496#ifdef BELOS_TEUCHOS_TIME_MONITOR
1497 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1502 Teuchos::Array<Teuchos::RCP<MV> > MQ(Q.size());
1503 Teuchos::RCP<const MV> tempQ;
1506 for (
int i=0; i<Q.size(); i++) {
1509 for (
int ii=0; ii<qcs[i]; ii++) {
1514 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], 1, 1, ii, j );
1532 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j );
1539 for (
int num_ortho_steps=1; num_ortho_steps <
max_ortho_steps_; ++num_ortho_steps) {
1541 for (
int i=0; i<Q.size(); i++) {
1542 Teuchos::SerialDenseMatrix<int,ScalarType> C2( qcs[i], 1 );
1545 for (
int ii=0; ii<qcs[i]; ii++) {
1550 Teuchos::SerialDenseMatrix<int,ScalarType> tempC2( Teuchos::View, C2, 1, 1, ii );
1558 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, *C[i], qcs[i], 1, 0, j );
1567 else if (xc <= qcs[i]) {
1578#ifdef BELOS_TEUCHOS_TIME_MONITOR
1579 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1585 if (SCT::magnitude(newDot[0]) < SCT::magnitude(oldDot[0]*
sing_tol_)) {
1591 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1604 Teuchos::RCP<MV> tempXj =
MVT::Clone( X, 1 );
1605 Teuchos::RCP<MV> tempMXj;
1615#ifdef BELOS_TEUCHOS_TIME_MONITOR
1616 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1623 for (
int i=0; i<Q.size(); i++) {
1624 Teuchos::SerialDenseMatrix<int,ScalarType> product( qcs[i], 1 );
1627 for (
int ii=0; ii<qcs[i]; ii++) {
1631 Teuchos::SerialDenseMatrix<int,ScalarType> tempC( Teuchos::View, product, 1, 1, ii );
1645 else if (xc <= qcs[i]) {
1655#ifdef BELOS_TEUCHOS_TIME_MONITOR
1656 Teuchos::TimeMonitor normTimer( *timerNorm_ );
1662 if ( SCT::magnitude(newDot[0]) >= SCT::magnitude(oldDot[0]*
sing_tol_) ) {
1663 ScalarType diag = SCT::squareroot(SCT::magnitude(newDot[0]));
1671 MVT::MvAddMv( ONE/diag, *tempMXj, ZERO, *tempMXj, *MXj );
1691 template<
class ScalarType,
class MV,
class OP>
1694 using Teuchos::ParameterList;
1695 using Teuchos::parameterList;
1698 RCP<ParameterList> params = parameterList (
"IMGS");
1703 "Maximum number of orthogonalization passes (includes the "
1704 "first). Default is 2, since \"twice is enough\" for Krylov "
1707 "Block reorthogonalization threshold.");
1709 "Singular block detection threshold.");
1714 template<
class ScalarType,
class MV,
class OP>
1717 using Teuchos::ParameterList;
1722 params->set (
"maxNumOrthogPasses",
1724 params->set (
"blkTol",
1726 params->set (
"singTol",
Belos header file which uses auto-configuration information to include necessary C++ headers.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
MagnitudeType sing_tol_
Singular block detection threshold.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &plist)
int max_ortho_steps_
Max number of (re)orthogonalization steps, including the first.
MultiVecTraits< ScalarType, MV > MVT
static const MagnitudeType sing_tol_fast_
Singular block detection threshold (fast).
const std::string & getLabel() const
This method returns the label being used by the timers in the orthogonalization manager.
Teuchos::RCP< Teuchos::ParameterList > defaultParams_
Default parameter list.
Teuchos::ScalarTraits< ScalarType > SCT
static const MagnitudeType sing_tol_default_
Singular block detection threshold (default).
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.
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.
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.
~IMGSOrthoManager()
Destructor.
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 ...
MagnitudeType blk_tol_
Block reorthogonalization tolerance.
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 ,...
void setLabel(const std::string &label)
This method sets the label used by the timers in the orthogonalization manager.
static const int max_ortho_steps_fast_
Max number of (re)orthogonalization steps, including the first (fast).
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.
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.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
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.
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 blk_tol_fast_
Block reorthogonalization threshold (fast).
static const int max_ortho_steps_default_
Max number of (re)orthogonalization steps, including the first (default).
MagnitudeType getBlkTol() const
Return parameter for block re-orthogonalization threshhold.
MagnitudeType getSingTol() const
Return parameter for singular block detection.
void setSingTol(const MagnitudeType sing_tol)
Set parameter for singular block detection.
std::string label_
Label for timers.
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.
void setBlkTol(const MagnitudeType blk_tol)
Set parameter for block re-orthogonalization threshhold.
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
OperatorTraits< ScalarType, MV, OP > OPT
static const MagnitudeType blk_tol_default_
Block reorthogonalization threshold (default).
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 .
Teuchos::ScalarTraits< MagnitudeType > MGT
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.
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...
MatOrthoManager(Teuchos::RCP< const OP > Op=Teuchos::null)
Default constructor.
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.
Teuchos::RCP< const OP > _Op
Traits class which defines basic operations on multivectors.
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
static void MvScale(MV &mv, const ScalarType alpha)
Scale each element of the vectors in mv with alpha.
static void MvTimesMatAddMv(const ScalarType alpha, const MV &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, MV &mv)
Update mv with .
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
static void MvDot(const MV &mv, const MV &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
static void Assign(const MV &A, MV &mv)
mv := A
Class which defines basic traits for the operator type.
static void Apply(const OP &Op, const MV &x, MV &y, ETrans trans=NOTRANS)
Apply Op to x, putting the result into y.
Exception thrown to signal error in an orthogonalization manager method.
Teuchos::RCP< Teuchos::ParameterList > getIMGSFastParameters()
"Fast" but possibly unsafe or less accurate parameters.
Teuchos::RCP< Teuchos::ParameterList > getIMGSDefaultParameters()
"Default" parameters for robustness and accuracy.