| 
    Anasazi Version of the Day
    
   | 
 
An implementation of the Anasazi::GenOrthoManager that performs orthogonalization using iterated classical Gram-Schmidt. More...
#include <AnasaziICGSOrthoManager.hpp>
  
Public Member Functions | |
Constructor/Destructor  | |
| ICGSOrthoManager (Teuchos::RCP< const OP > Op=Teuchos::null, int numIters=2, typename Teuchos::ScalarTraits< ScalarType >::magnitudeType eps=0.0, typename Teuchos::ScalarTraits< ScalarType >::magnitudeType tol=0.20) | |
| Constructor specifying the operator defining the inner product as well as the number of orthogonalization iterations.   | |
| ~ICGSOrthoManager () | |
| Destructor.   | |
Methods implementing Anasazi::GenOrthoManager  | |
| void | projectGen (MV &S, Teuchos::Array< Teuchos::RCP< const MV > > X, Teuchos::Array< Teuchos::RCP< const MV > > Y, bool isBiOrtho, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< MV > MS=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MX=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null)), Teuchos::Array< Teuchos::RCP< const MV > > MY=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const | 
| Applies a series of generic projectors.   | |
| int | projectAndNormalizeGen (MV &S, Teuchos::Array< Teuchos::RCP< const MV > > X, Teuchos::Array< Teuchos::RCP< const MV > > Y, bool isBiOrtho, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null, Teuchos::RCP< MV > MS=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MX=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null)), Teuchos::Array< Teuchos::RCP< const MV > > MY=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const | 
| Applies a series of generic projectors and returns an orthonormal basis for the residual data.   | |
Methods implementing Anasazi::MatOrthoManager  | |
| void | projectMat (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const | 
Given a list of mutually orthogonal and internally orthonormal bases Q, this method projects a multivector X 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().   | |
| int | normalizeMat (MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null) const | 
This method takes a multivector X and attempts to compute an orthonormal basis for  | |
| int | projectAndNormalizeMat (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null, Teuchos::RCP< MV > MX=Teuchos::null, Teuchos::Array< Teuchos::RCP< const MV > > MQ=Teuchos::tuple(Teuchos::RCP< const MV >(Teuchos::null))) const | 
Given a set of bases Q[i] and a multivector X, this method computes an orthonormal basis for  | |
Error methods  | |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthonormErrorMat (const MV &X, Teuchos::RCP< const MV > MX=Teuchos::null) 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.   | |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthogErrorMat (const MV &X1, const MV &X2, Teuchos::RCP< const MV > MX1, Teuchos::RCP< const MV > MX2) 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.   | |
  Public Member Functions inherited from Anasazi::GenOrthoManager< ScalarType, MV, OP > | |
| GenOrthoManager (Teuchos::RCP< const OP > Op=Teuchos::null) | |
| Default constructor.   | |
| virtual | ~GenOrthoManager () | 
| Destructor.   | |
  Public Member Functions inherited from Anasazi::MatOrthoManager< ScalarType, MV, OP > | |
| MatOrthoManager (Teuchos::RCP< const OP > Op=Teuchos::null) | |
| Default constructor.   | |
| virtual | ~MatOrthoManager () | 
| Destructor.   | |
| virtual void | setOp (Teuchos::RCP< const OP > Op) | 
| Set operator used for inner product.   | |
| virtual Teuchos::RCP< const OP > | getOp () const | 
| Get operator used for inner product.   | |
| int | getOpCounter () const | 
| Retrieve operator counter.   | |
| void | resetOpCounter () | 
| Reset the operator counter to zero.   | |
| void | innerProdMat (const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z, Teuchos::RCP< const MV > MX=Teuchos::null, Teuchos::RCP< const MV > MY=Teuchos::null) const | 
| Provides a matrix-based inner product.   | |
| void | normMat (const MV &X, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, Teuchos::RCP< const MV > MX=Teuchos::null) const | 
| Provides the norm induced by the matrix-based inner product.   | |
| void | innerProd (const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z) const | 
| Implements the interface OrthoManager::innerProd().   | |
| void | norm (const MV &X, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec) const | 
| Implements the interface OrthoManager::norm().   | |
| void | project (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null))) const | 
| Implements the interface OrthoManager::project().   | |
| int | normalize (MV &X, Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null) const | 
| Implements the interface OrthoManager::normalize().   | |
| int | projectAndNormalize (MV &X, Teuchos::Array< Teuchos::RCP< const MV > > Q, Teuchos::Array< Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > > C=Teuchos::tuple(Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)), Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B=Teuchos::null) const | 
| Implements the interface OrthoManager::projectAndNormalize().   | |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthonormError (const MV &X) const | 
| Implements the interface OrthoManager::orthonormError().   | |
| Teuchos::ScalarTraits< ScalarType >::magnitudeType | orthogError (const MV &X1, const MV &X2) const | 
| Implements the interface OrthoManager::orthogError().   | |
  Public Member Functions inherited from Anasazi::OrthoManager< ScalarType, MV > | |
| OrthoManager () | |
| Default constructor.   | |
| virtual | ~OrthoManager () | 
| Destructor.   | |
Accessor routines | |
| void | setNumIters (int numIters) | 
| Set parameter for re-orthogonalization threshold.   | |
| int | getNumIters () const | 
| Return parameter for re-orthogonalization threshold.   | |
Additional Inherited Members | |
  Protected Attributes inherited from Anasazi::MatOrthoManager< ScalarType, MV, OP > | 
An implementation of the Anasazi::GenOrthoManager that performs orthogonalization using iterated classical Gram-Schmidt.
Definition at line 39 of file AnasaziICGSOrthoManager.hpp.
| Anasazi::ICGSOrthoManager< ScalarType, MV, OP >::ICGSOrthoManager | ( | Teuchos::RCP< const OP > | Op = Teuchos::null,  | 
        
| int | numIters = 2,  | 
        ||
| typename Teuchos::ScalarTraits< ScalarType >::magnitudeType | eps = 0.0,  | 
        ||
| typename Teuchos::ScalarTraits< ScalarType >::magnitudeType | tol = 0.20  | 
        ||
| ) | 
Constructor specifying the operator defining the inner product as well as the number of orthogonalization iterations.
Definition at line 379 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  inline | 
Destructor.
Definition at line 58 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  virtual | 
Applies a series of generic projectors.
Given a list of bases X[i] and Y[i] (a projection pair), this method takes a multivector S and applies the projectors  
![\[
     P_{X[i],Y[i]} S = S - X[i] \langle Y[i], X[i] \rangle^{-1} \langle Y[i], S \rangle\ .
  \]](form_34.png)
 This operation projects S onto the space orthogonal to the Y[i], along the range of the X[i]. The inner product specified by 
The method also returns the coefficients C[i] associated with each projection pair, so that  
![\[
     S_{in} = S_{out} + \sum_i X[i] C[i]
  \]](form_36.png)
and therefore
![\[
     C[i] = \langle Y[i], X[i] \rangle^{-1} \langle Y[i], S \rangle\ .
  \]](form_37.png)
Lastly, for reasons of efficiency, the user must specify whether the projection pairs are bi-orthonormal with respect to innerProd(), i.e., whether ![$\langle Y[i], X[i] \rangle = I$](form_38.png)

S and the projection pairs under the inner product operator getOp().
projectGen() is implemented to apply the projectors via an iterated Classical Gram-Schmidt, where the iteration is performed getNumIters() number of times.
| S | [in/out] The multivector to be modified. On output, the columns of S will be orthogonal to each Y[i], satisfying  
 
  | 
| X | [in] Multivectors for bases under which  | 
| Y | [in] Multivectors for bases to which  | 
| isBiortho | [in] A flag specifying whether the bases X[i] and Y[i] are biorthonormal, i.e,. whether  | 
| C | [out] Coefficients for reconstructing X[i]. If C[i] is a non-null pointer and C[i] matches the dimensions of S and X[i], then the coefficients computed during the orthogonalization routine will be stored in the matrix C[i].If C[i] points to a Teuchos::SerialDenseMatrix with size inconsistent with S and , then a std::invalid_argument exception will be thrown.Otherwise, if C.size() < i or C[i] is a null pointer, the caller will not have access to the computed coefficients C[i]. | 
| MS | [in/out] If specified by the user, on input MS is required to be the image of S under the operator getOp(). On output, MS will be updated to reflect the changes in S. | 
| MX | [in] If specified by the user, on MX[i] is required to be the image of X[i] under the operator getOp().  | 
| MY | [in] If specified by the user, on MY[i] is required to be the image of Y[i] under the operator getOp(). | 
X[i] != Teuchos::null or Y[i] != Teuchos::null, then X[i] and Y[i] are required to have the same number of columns, and each should have the same number of rows as S. i != j, ![$\langle Y[i], X[j] \rangle == 0$](form_44.png)
biOrtho == true, ![$\langle Y[i], X[i]\rangle == I$](form_45.png)
biOrtho == false, then ![$\langle Y[i], X[i]\rangle$](form_46.png)
X[i] and Y[i] have 
S has 
C[i] if specified must be 
Implements Anasazi::GenOrthoManager< ScalarType, MV, OP >.
Definition at line 515 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  virtual | 
Applies a series of generic projectors and returns an orthonormal basis for the residual data.
Given a list of bases X[i] and Y[i] (a projection pair), this method takes a multivector S and applies the projectors  
![\[
     P_{X[i],Y[i]} S = S - X[i] \langle Y[i], X[i] \rangle^{-1} \langle Y[i], S \rangle\ .
  \]](form_34.png)
 These operation project S onto the space orthogonal to the range of the Y[i], along the range of X[i]. The inner product specified by 
The method returns in S an orthonormal basis for the residual  
![\[
     \left( \prod_{i} P_{X[i],Y[i]} \right) S_{in} = S_{out} B\ ,
  \]](form_51.png)
 where B contains the (not necessarily triangular) coefficients of the residual with respect to the new basis.
The method also returns the coefficients C[i] and B associated with each projection pair, so that  
![\[
     S_{in} = S_{out} B + \sum_i X[i] C[i]
  \]](form_52.png)
and
![\[
     C[i] = \langle Y[i], X[i] \rangle^{-1} \langle Y[i], S \rangle\ .
  \]](form_37.png)
Lastly, for reasons of efficiency, the user must specify whether the projection pairs are bi-orthonormal with respect to innerProd(), i.e., whether ![$\langle Y[i], X[i] \rangle = I$](form_38.png)
S and the projection pairs under the inner product operator getOp(). 
| S | [in/out] The multivector to be modified. On output, the columns of S will be orthogonal to each Y[i], satisfying  
 
 m is the number of rows in S, n is the number of columns in S, and rank is the value returned from the method. | 
| X | [in] Multivectors for bases under which  | 
| Y | [in] Multivectors for bases to which  | 
| isBiortho | [in] A flag specifying whether the bases X[i] and Y[i] are biorthonormal, i.e,. whether  | 
| C | [out] Coefficients for reconstructing X[i]. 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] points to a Teuchos::SerialDenseMatrix with size inconsistent with S and , then a std::invalid_argument exception will be thrown.Otherwise, if C.size() < i or C[i] is a null pointer, the caller will not have access to the computed coefficients C[i]. | 
| B | [out] The coefficients of the original S with respect to the computed basis. If B is a non-null pointer and B matches the dimensions of B, then the coefficients computed during the orthogonalization routine will be stored in B, similar to calling innerProd( Sout, Sin, B ); 
void innerProd(const MV &X, const MV &Y, Teuchos::SerialDenseMatrix< int, ScalarType > &Z) const Implements the interface OrthoManager::innerProd(). Definition AnasaziMatOrthoManager.hpp:315 B points to a Teuchos::SerialDenseMatrix with size inconsistent with S, then a std::invalid_argument exception will be thrown.Otherwise, if B is null, the caller will not have access to the computed coefficients.The normalization uses classical Gram-Schmidt iteration, so that B is an upper triangular matrix with positive diagonal elements. | 
| MS | [in/out] If specified by the user, on input MS is required to be the image of S under the operator getOp(). On output, MS will be updated to reflect the changes in S. | 
| MX | [in] If specified by the user, on MX[i] is required to be the image of X[i] under the operator getOp().  | 
| MY | [in] If specified by the user, on MY[i] is required to be the image of Y[i] under the operator getOp(). | 
X[i] != Teuchos::null or Y[i] != Teuchos::null, then X[i] and Y[i] are required to have the same number of columns, and each should have the same number of rows as S. i != j, ![$\langle Y[i], X[j] \rangle == 0$](form_44.png)
biOrtho == true, ![$\langle Y[i], X[i]\rangle == I$](form_45.png)
biOrtho == false, then ![$\langle Y[i], X[i]\rangle$](form_46.png)
X[i] and Y[i] have 
S has 
C[i] if specified must be 
S has 
B if specified must be 
Implements Anasazi::GenOrthoManager< ScalarType, MV, OP >.
Definition at line 860 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  virtual | 
Given a list of mutually orthogonal and internally orthonormal bases Q, this method projects a multivector X 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(). 
This method calls projectGen() as follows:
See projectGen() for argument requirements.
Implements Anasazi::MatOrthoManager< ScalarType, MV, OP >.
Definition at line 434 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  virtual | 
This method takes a multivector X and attempts to compute an orthonormal basis for 
This method calls projectAndNormalizeGen() as follows:
See projectAndNormalizeGen() for argument requirements.
Implements Anasazi::MatOrthoManager< ScalarType, MV, OP >.
Definition at line 450 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  virtual | 
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])$](form_16.png)
This method calls projectAndNormalizeGen() as follows:
See projectAndNormalizeGen() for argument requirements.
Implements Anasazi::MatOrthoManager< ScalarType, MV, OP >.
Definition at line 499 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  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 Anasazi::MatOrthoManager< ScalarType, MV, OP >.
Definition at line 405 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  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 Anasazi::MatOrthoManager< ScalarType, MV, OP >.
Definition at line 422 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  inline | 
Set parameter for re-orthogonalization threshold.
Definition at line 350 of file AnasaziICGSOrthoManager.hpp.
      
  | 
  inline | 
Return parameter for re-orthogonalization threshold.
Definition at line 357 of file AnasaziICGSOrthoManager.hpp.