#include <Belos_IMGS_OrthoManager_MP_Vector.hpp>
Inherits MatOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >.
|
| 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. More...
|
| |
| 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. More...
|
| |
| 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. More...
|
| |
| 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 |
| |
|
| 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 projects it onto the space orthogonal to the individual Q[i], optionally returning the coefficients of X for the individual Q[i]. All of this is done with respect to the inner product innerProd(). More...
|
| |
| 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. More...
|
| |
| 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 , with respect to innerProd(). More...
|
| |
| 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. More...
|
| |
| 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 . More...
|
| |
|
| 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 the difference innerProd(X,Y) - I. More...
|
| |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthonormError (const MV &X, Teuchos::RCP< const MV > MX) const |
| | This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of the difference innerProd(X,Y) - I. The method has the option of exploiting a caller-provided MX. More...
|
| |
| 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 of innerProd(X,Y). More...
|
| |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthogError (const MV &X1, Teuchos::RCP< const MV > MX1, const MV &X2) const |
| | This method computes the error in orthogonality of two multivectors, measured as the Frobenius norm of innerProd(X,Y). The method has the option of exploiting a caller-provided MX. More...
|
| |
|
| void | setLabel (const std::string &label) |
| | This method sets the label used by the timers in the orthogonalization manager. More...
|
| |
| const std::string & | getLabel () const |
| | This method returns the label being used by the timers in the orthogonalization manager. More...
|
| |
template<class Storage, class MV, class OP>
class Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >
Definition at line 60 of file Belos_IMGS_OrthoManager_MP_Vector.hpp.
◆ ScalarType
template<class Storage , class MV , class OP >
◆ MagnitudeType
template<class Storage , class MV , class OP >
◆ MGT
template<class Storage , class MV , class OP >
◆ SCT
template<class Storage , class MV , class OP >
◆ MVT
template<class Storage , class MV , class OP >
◆ OPT
template<class Storage , class MV , class OP >
◆ IMGSOrthoManager() [1/2]
template<class Storage , class MV , class OP >
◆ IMGSOrthoManager() [2/2]
template<class Storage , class MV , class OP >
| Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::IMGSOrthoManager |
( |
const Teuchos::RCP< Teuchos::ParameterList > & |
plist, |
|
|
const std::string & |
label = "Belos", |
|
|
Teuchos::RCP< const OP > |
Op = Teuchos::null |
|
) |
| |
|
inline |
◆ ~IMGSOrthoManager()
template<class Storage , class MV , class OP >
◆ setParameterList()
template<class Storage , class MV , class OP >
| void Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::setParameterList |
( |
const Teuchos::RCP< Teuchos::ParameterList > & |
plist | ) |
|
|
inline |
◆ getValidParameters()
template<class Storage , class MV , class OP >
| Teuchos::RCP<const Teuchos::ParameterList> Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::getValidParameters |
( |
| ) |
const |
|
inline |
◆ setBlkTol()
template<class Storage , class MV , class OP >
◆ setSingTol()
template<class Storage , class MV , class OP >
◆ getBlkTol()
template<class Storage , class MV , class OP >
◆ getSingTol()
template<class Storage , class MV , class OP >
◆ project() [1/2]
template<class Storage , class MV , class OP >
| void Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::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 projects it onto the space orthogonal to the individual Q[i], optionally returning the coefficients of X for the individual Q[i]. All of this is done with respect to the inner product innerProd().
After calling this routine, X will be orthogonal to each of the Q[i].
The method uses either one or two steps of modified Gram-Schmidt. The algebraically equivalent projection matrix is
, if Op is the matrix specified for use in the inner product. Note, this is not an orthogonal projector.
- Parameters
-
| X | [in/out] The multivector to be modified. On output, X will be orthogonal to Q[i] with respect to innerProd(). |
| MX | [in/out] The image of X under the operator Op. If : On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If or : MX is not referenced. |
| C | [out] The coefficients of X in the *Q[i], with respect to innerProd(). If C[i] is a non-null pointer and *C[i] matches the dimensions of X and *Q[i], then the coefficients computed during the orthogonalization routine will be stored in the matrix *C[i]. If C[i] is a non-null pointer whose size does not match the dimensions of X and *Q[i], then a std::invalid_argument std::exception will be thrown. Otherwise, if C.size() < i or C[i] is a null pointer, then the orthogonalization manager will declare storage for the coefficients and the user will not have access to them. |
| Q | [in] A list of multivector bases specifying the subspaces to be orthogonalized against. Each Q[i] is assumed to have orthonormal columns, and the Q[i] are assumed to be mutually orthogonal. |
◆ project() [2/2]
template<class Storage , class MV , class OP >
| void Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::project |
( |
MV & |
X, |
|
|
Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > |
C, |
|
|
Teuchos::ArrayView< Teuchos::RCP< const MV > > |
Q |
|
) |
| const |
|
inline |
◆ normalize() [1/2]
template<class Storage , class MV , class OP >
| int Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::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
, with respect to innerProd().
The method uses modified Gram-Schmidt, so that the coefficient matrix B is upper triangular.
This routine returns an integer rank stating the rank of the computed basis. If X does not have full rank and the normalize() routine does not attempt to augment the subspace, then rank may be smaller than the number of columns in X. In this case, only the first rank columns of output X and first rank rows of B will be valid.
The method attempts to find a basis with dimension the same as the number of columns in X. It does this by augmenting linearly dependant vectors in X with random directions. A finite number of these attempts will be made; therefore, it is possible that the dimension of the computed basis is less than the number of vectors in X.
- Parameters
-
| X | [in/out] The multivector to the modified. On output, X will have some number of orthonormal columns (with respect to innerProd()). |
| MX | [in/out] The image of X under the operator Op. If : On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If or : MX is not referenced. |
| B | [out] The coefficients of the original X with respect to the computed basis. The first rows in B corresponding to the valid columns in X will be upper triangular. |
- Returns
- Rank of the basis computed by this method.
◆ normalize() [2/2]
template<class Storage , class MV , class OP >
◆ projectAndNormalizeWithMxImpl()
template<class Storage , class MV , class OP >
| virtual int Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::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 |
|
protectedvirtual |
Given a set of bases Q[i] and a multivector X, this method computes an orthonormal basis for
.
This routine returns an integer rank stating the rank of the computed basis. If the subspace
does not have dimension as large as the number of columns of X and the orthogonalization manager doe not attempt to augment the subspace, then rank may be smaller than the number of columns of X. In this case, only the first rank columns of output X and first rank rows of B will be valid.
The method attempts to find a basis with dimension the same as the number of columns in X. It does this by augmenting linearly dependant vectors with random directions. A finite number of these attempts will be made; therefore, it is possible that the dimension of the computed basis is less than the number of vectors in X.
- Parameters
-
| X | [in/out] The multivector to the modified. On output, the relevant rows of X will be orthogonal to the Q[i] and will have orthonormal columns (with respect to innerProd()). |
| MX | [in/out] The image of X under the operator Op. If : On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If or : MX is not referenced. |
| C | [out] The coefficients of the original X in the Q[i], with respect to innerProd(). If C[i] is a non-null pointer and *C[i] matches the dimensions of X and *Q[i], then the coefficients computed during the orthogonalization routine will be stored in the matrix C[i]. If C[i] is a non-null pointer whose size does not match the dimensions of X and *Q[i], then *C[i] will first be resized to the correct size. This will destroy the original contents of the matrix. (This is a change from previous behavior, in which a std::invalid_argument exception was thrown if *C[i] was of the wrong size.) Otherwise, if C.size() < i<> or C[i] is a null pointer, then the orthogonalization manager will declare storage for the coefficients and the user will not have access to them. |
| B | [out] The coefficients of the original X with respect to the computed basis. The first rows in B corresponding to the valid columns in X will be upper triangular. |
| Q | [in] A list of multivector bases specifying the subspaces to be orthogonalized against. Each Q[i] is assumed to have orthonormal columns, and the Q[i] are assumed to be mutually orthogonal. |
- Returns
Rank of the basis computed by this method.
◆ orthonormError() [1/2]
template<class Storage , class MV , class OP >
This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of the difference innerProd(X,Y) - I.
Definition at line 380 of file Belos_IMGS_OrthoManager_MP_Vector.hpp.
◆ orthonormError() [2/2]
template<class Storage , class MV , class OP >
| Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::orthonormError |
( |
const MV & |
X, |
|
|
Teuchos::RCP< const MV > |
MX |
|
) |
| const |
This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of the difference innerProd(X,Y) - I. The method has the option of exploiting a caller-provided MX.
◆ orthogError() [1/2]
template<class Storage , class MV , class OP >
◆ orthogError() [2/2]
template<class Storage , class MV , class OP >
| Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::orthogError |
( |
const MV & |
X1, |
|
|
Teuchos::RCP< const MV > |
MX1, |
|
|
const MV & |
X2 |
|
) |
| const |
This method computes the error in orthogonality of two multivectors, measured as the Frobenius norm of innerProd(X,Y). The method has the option of exploiting a caller-provided MX.
◆ setLabel()
template<class Storage , class MV , class OP >
This method sets the label used by the timers in the orthogonalization manager.
◆ getLabel()
template<class Storage , class MV , class OP >
◆ findBasis()
template<class Storage , class MV , class OP >
| int Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::findBasis |
( |
MV & |
X, |
|
|
Teuchos::RCP< MV > |
MX, |
|
|
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > |
C, |
|
|
bool |
completeBasis, |
|
|
int |
howMany = -1 |
|
) |
| const |
|
private |
Routine to find an orthonormal basis for X.
◆ blkOrtho1()
template<class Storage , class MV , class OP >
| bool Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::blkOrtho1 |
( |
MV & |
X, |
|
|
Teuchos::RCP< MV > |
MX, |
|
|
Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > |
C, |
|
|
Teuchos::ArrayView< Teuchos::RCP< const MV > > |
Q |
|
) |
| const |
|
private |
Routine to compute the block orthogonalization.
◆ blkOrtho()
template<class Storage , class MV , class OP >
| bool Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::blkOrtho |
( |
MV & |
X, |
|
|
Teuchos::RCP< MV > |
MX, |
|
|
Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > |
C, |
|
|
Teuchos::ArrayView< Teuchos::RCP< const MV > > |
Q |
|
) |
| const |
|
private |
Routine to compute the block orthogonalization.
◆ blkOrthoSing()
template<class Storage , class MV , class OP >
| int Belos::IMGSOrthoManager< Sacado::MP::Vector< Storage >, MV, OP >::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 |
|
private |
Project X against QQ and normalize X, one vector at a time
- Note
- QQ is called QQ, rather than Q, because we convert it internally from an ArrayView to an Array (named Q inside). This is because the C++ compiler doesn't know how to do type inference (Array has a constructor that takes an ArrayView input). This routine wants an Array rather than an ArrayView internally, because it likes to add (via push_back()) and remove (via resize()) elements to the Q array. Remember that Arrays can be passed by value, just like std::vector objects, so this routine can add whatever it likes to the Q array without changing it from the caller's perspective.
◆ max_ortho_steps_default_
template<class Storage , class MV , class OP >
◆ blk_tol_default_
template<class Storage , class MV , class OP >
◆ sing_tol_default_
template<class Storage , class MV , class OP >
◆ max_ortho_steps_fast_
template<class Storage , class MV , class OP >
◆ blk_tol_fast_
template<class Storage , class MV , class OP >
◆ sing_tol_fast_
template<class Storage , class MV , class OP >
◆ max_ortho_steps_
template<class Storage , class MV , class OP >
◆ blk_tol_
template<class Storage , class MV , class OP >
◆ sing_tol_
template<class Storage , class MV , class OP >
◆ label_
template<class Storage , class MV , class OP >
◆ defaultParams_
template<class Storage , class MV , class OP >
The documentation for this class was generated from the following file: