Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | List of all members
Belos::ICGSOrthoManager< ScalarType, MV, OP > Class Template Reference

An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps of classical Gram-Schmidt. More...

#include <BelosICGSOrthoManager.hpp>

Inheritance diagram for Belos::ICGSOrthoManager< ScalarType, MV, OP >:
Inheritance graph
[legend]

Public Member Functions

Teuchos::RCP< const
Teuchos::ParameterList
getFastParameters () const
 "Fast" but possibly unsafe or less accurate parameters. More...
 
Constructor/Destructor
 ICGSOrthoManager (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. More...
 
 ICGSOrthoManager (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. More...
 
 ~ICGSOrthoManager ()
 Destructor. More...
 
Implementation of Teuchos::ParameterListAcceptorDefaultBase interface
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &plist)
 
Teuchos::RCP< const
Teuchos::ParameterList
getValidParameters () const
 
Accessor routines
void setBlkTol (const MagnitudeType blk_tol)
 Set parameter for block re-orthogonalization threshhold. More...
 
void setSingTol (const MagnitudeType sing_tol)
 Set parameter for singular block detection. More...
 
MagnitudeType getBlkTol () const
 Return parameter for block re-orthogonalization threshhold. More...
 
MagnitudeType getSingTol () const
 Return parameter for singular block detection. More...
 
Error methods
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...
 
Label methods
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...
 
- Public Member Functions inherited from Belos::MatOrthoManager< ScalarType, MV, OP >
 MatOrthoManager (Teuchos::RCP< const OP > Op=Teuchos::null)
 Default constructor. More...
 
virtual ~MatOrthoManager ()
 Destructor. More...
 
void setOp (Teuchos::RCP< const OP > Op)
 Set operator. More...
 
Teuchos::RCP< const OP > getOp () const
 Get operator. More...
 
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. More...
 
void innerProd (const MV &X, const MV &Y, Teuchos::RCP< const MV > MY, Teuchos::SerialDenseMatrix< int, ScalarType > &Z) const
 Provides the inner product defining the orthogonality concepts, using the provided operator. The method has the options of exploiting a caller-provided MX. More...
 
void norm (const MV &X, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec) const
 Provides the norm induced by innerProd(). More...
 
void norm (const MV &X, Teuchos::RCP< const MV > MX, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec) const
 Compute norm of each column of X. More...
 
int projectAndNormalize (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 $colspan(X) - \sum_i colspan(Q[i])$. More...
 
- Public Member Functions inherited from Belos::OrthoManager< ScalarType, MV >
 OrthoManager ()
 Default constructor. More...
 
virtual ~OrthoManager ()
 Destructor. More...
 
int projectAndNormalize (MV &X, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
 Project X against the Q[i] and normalize X. More...
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
RCP< ParameterListgetNonconstParameterList ()
 
RCP< ParameterListunsetParameterList ()
 
RCP< const ParameterListgetParameterList () const
 

Static Public Attributes

Default orthogonalization constants
static const int max_ortho_steps_default_ = 2
 Max number of (re)orthogonalization steps, including the first (default). More...
 
static const MagnitudeType blk_tol_default_
 Block reorthogonalization threshold (default). More...
 
static const MagnitudeType sing_tol_default_ = 10*Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::eps()
 Singular block detection threshold (default). More...
 
static const int max_ortho_steps_fast_ = 1
 Max number of (re)orthogonalization steps, including the first (fast). More...
 
static const MagnitudeType blk_tol_fast_ = Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero()
 Block reorthogonalization threshold (fast). More...
 
static const MagnitudeType sing_tol_fast_ = Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero()
 Singular block detection threshold (fast). More...
 

Orthogonalization methods

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 $colspan(X)$, 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 $colspan(X) - \sum_i colspan(Q[i])$. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Belos::MatOrthoManager< ScalarType, MV, OP >
virtual int projectAndNormalizeImpl (MV &X, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B, Teuchos::ArrayView< Teuchos::RCP< const MV > > Q) const
 
- Protected Member Functions inherited from Belos::OrthoManager< ScalarType, MV >
- Protected Member Functions inherited from Teuchos::ParameterListAcceptorDefaultBase
void setMyParamList (const RCP< ParameterList > &paramList)
 
RCP< ParameterListgetMyNonconstParamList ()
 
RCP< const ParameterListgetMyParamList () const
 
- Protected Attributes inherited from Belos::MatOrthoManager< ScalarType, MV, OP >
Teuchos::RCP< const OP > _Op
 
bool _hasOp
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Belos::ICGSOrthoManager< ScalarType, MV, OP >

An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps of classical Gram-Schmidt.

Author
Chris Baker, Ulrich Hetmaniuk, Rich Lehoucq, and Heidi Thornquist

Definition at line 80 of file BelosICGSOrthoManager.hpp.

Constructor & Destructor Documentation

template<class ScalarType, class MV, class OP>
Belos::ICGSOrthoManager< ScalarType, MV, OP >::ICGSOrthoManager ( 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_ 
)
inline

Constructor specifying re-orthogonalization tolerance.

Definition at line 95 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Belos::ICGSOrthoManager< ScalarType, MV, OP >::ICGSOrthoManager ( const Teuchos::RCP< Teuchos::ParameterList > &  plist,
const std::string &  label = "Belos",
Teuchos::RCP< const OP >  Op = Teuchos::null 
)
inline

Constructor that takes a list of parameters.

Definition at line 125 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Belos::ICGSOrthoManager< ScalarType, MV, OP >::~ICGSOrthoManager ( )
inline

Destructor.

Definition at line 155 of file BelosICGSOrthoManager.hpp.

Member Function Documentation

template<class ScalarType, class MV, class OP>
void Belos::ICGSOrthoManager< ScalarType, MV, OP >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  plist)
inlinevirtual

Implements Teuchos::ParameterListAcceptor.

Definition at line 162 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const Teuchos::ParameterList> Belos::ICGSOrthoManager< ScalarType, MV, OP >::getValidParameters ( ) const
inlinevirtual

Reimplemented from Teuchos::ParameterListAcceptor.

Definition at line 234 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::RCP<const Teuchos::ParameterList> Belos::ICGSOrthoManager< ScalarType, MV, OP >::getFastParameters ( ) const
inline

"Fast" but possibly unsafe or less accurate parameters.

Use this parameter list when you care more about speed than accuracy of the orthogonalization.

Definition at line 250 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
void Belos::ICGSOrthoManager< ScalarType, MV, OP >::setBlkTol ( const MagnitudeType  blk_tol)
inline

Set parameter for block re-orthogonalization threshhold.

Definition at line 272 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
void Belos::ICGSOrthoManager< ScalarType, MV, OP >::setSingTol ( const MagnitudeType  sing_tol)
inline

Set parameter for singular block detection.

Definition at line 286 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::getBlkTol ( ) const
inline

Return parameter for block re-orthogonalization threshhold.

Definition at line 300 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::getSingTol ( ) const
inline

Return parameter for singular block detection.

Definition at line 303 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
void Belos::ICGSOrthoManager< ScalarType, 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
virtual

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 classical Gram-Schmidt. The algebraically equivalent projection matrix is $P_Q = I - Q Q^H Op$, 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 $ MX != 0$: On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If $ MX == 0$ or $ Op == 0$: 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.

Implements Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 843 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
void Belos::ICGSOrthoManager< ScalarType, MV, OP >::project ( MV &  X,
Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > >  C,
Teuchos::ArrayView< Teuchos::RCP< const MV > >  Q 
) const
inlinevirtual

This method calls project(X,Teuchos::null,C,Q); see documentation for that function.

Reimplemented from Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 345 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
int Belos::ICGSOrthoManager< ScalarType, MV, OP >::normalize ( MV &  X,
Teuchos::RCP< MV >  MX,
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >  B 
) const
virtual

This method takes a multivector X and attempts to compute an orthonormal basis for $colspan(X)$, with respect to innerProd().

The method uses classical 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 $ MX != 0$: On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If $ MX == 0$ or $ Op == 0$: 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.

Implements Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 826 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
int Belos::ICGSOrthoManager< ScalarType, MV, OP >::normalize ( MV &  X,
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >  B 
) const
inlinevirtual

This method calls normalize(X,Teuchos::null,B); see documentation for that function.

Reimplemented from Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 383 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
int Belos::ICGSOrthoManager< ScalarType, 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 $colspan(X) - \sum_i colspan(Q[i])$.

This routine returns an integer rank stating the rank of the computed basis. If the subspace $colspan(X) - \sum_i colspan(Q[i])$ 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 $ MX != 0$: On input, this is expected to be consistent with X. On output, this is updated consistent with updates to X. If $ MX == 0$ or $ Op == 0$: 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.

Implements Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 644 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::orthonormError ( const MV &  X) const
inlinevirtual

This method computes the error in orthonormality of a multivector, measured as the Frobenius norm of the difference innerProd(X,Y) - I.

Reimplemented from Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 446 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::ScalarTraits< ScalarType >::magnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::orthonormError ( const MV &  X,
Teuchos::RCP< const MV >  MX 
) const
virtual

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.

Implements Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 616 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::orthogError ( const MV &  X1,
const MV &  X2 
) const
inlinevirtual

This method computes the error in orthogonality of two multivectors, measured as the Frobenius norm of innerProd(X,Y).

Reimplemented from Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 461 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::ScalarTraits< ScalarType >::magnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::orthogError ( const MV &  X1,
Teuchos::RCP< const MV >  MX1,
const MV &  X2 
) const
virtual

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.

Implements Belos::MatOrthoManager< ScalarType, MV, OP >.

Definition at line 631 of file BelosICGSOrthoManager.hpp.

template<class ScalarType , class MV , class OP >
void Belos::ICGSOrthoManager< ScalarType, MV, OP >::setLabel ( const std::string &  label)
virtual

This method sets the label used by the timers in the orthogonalization manager.

Implements Belos::OrthoManager< ScalarType, MV >.

Definition at line 589 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const std::string& Belos::ICGSOrthoManager< ScalarType, MV, OP >::getLabel ( ) const
inlinevirtual

This method returns the label being used by the timers in the orthogonalization manager.

Implements Belos::OrthoManager< ScalarType, MV >.

Definition at line 483 of file BelosICGSOrthoManager.hpp.

Member Data Documentation

template<class ScalarType, class MV, class OP>
const int Belos::ICGSOrthoManager< ScalarType, MV, OP >::max_ortho_steps_default_ = 2
static

Max number of (re)orthogonalization steps, including the first (default).

Definition at line 491 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const ICGSOrthoManager< ScalarType, MV, OP >::MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::blk_tol_default_
static
Initial value:
Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::eps() )

Block reorthogonalization threshold (default).

Definition at line 493 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const ICGSOrthoManager< ScalarType, MV, OP >::MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::sing_tol_default_ = 10*Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::eps()
static

Singular block detection threshold (default).

Definition at line 495 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const int Belos::ICGSOrthoManager< ScalarType, MV, OP >::max_ortho_steps_fast_ = 1
static

Max number of (re)orthogonalization steps, including the first (fast).

Definition at line 498 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const ICGSOrthoManager< ScalarType, MV, OP >::MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::blk_tol_fast_ = Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero()
static

Block reorthogonalization threshold (fast).

Definition at line 500 of file BelosICGSOrthoManager.hpp.

template<class ScalarType, class MV, class OP>
const ICGSOrthoManager< ScalarType, MV, OP >::MagnitudeType Belos::ICGSOrthoManager< ScalarType, MV, OP >::sing_tol_fast_ = Teuchos::ScalarTraits<typename ICGSOrthoManager<ScalarType,MV,OP>::MagnitudeType>::zero()
static

Singular block detection threshold (fast).

Definition at line 502 of file BelosICGSOrthoManager.hpp.


The documentation for this class was generated from the following file:

Generated on Mon Nov 4 2024 09:26:15 for Belos by doxygen 1.8.5