| 
    NOX
    Development
    
   | 
 
Anasazi operator for computing generalized eigenvalues using Cayley transformations. More...
#include <LOCA_AnasaziOperator_Cayley2Matrix.H>


Public Member Functions | |
| Cayley2Matrix (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &eigenParams, const Teuchos::RCP< Teuchos::ParameterList > &solverParams, const Teuchos::RCP< LOCA::TimeDependent::AbstractGroup > &grp) | |
| Constructor.  More... | |
| virtual | ~Cayley2Matrix () | 
| Destructor.  | |
| virtual const std::string & | label () const | 
| Return name of this operator.  | |
| virtual void | apply (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &output) const | 
| Apply the operator.  More... | |
| virtual void | preProcessSeedVector (NOX::Abstract::MultiVector &ivec) | 
| PreProcess the random seed vector.  More... | |
| virtual void | beginPostProcessing () | 
| Begin PostProcessing of eigenvalues.  More... | |
| virtual void | transformEigenvalue (double &ev_r, double &ev_i) const | 
| Transform eigenvalue.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | rayleighQuotient (NOX::Abstract::Vector &evec_r, NOX::Abstract::Vector &evec_i, double &rq_r, double &rq_i) const | 
| Compute Rayleigh quotient.  More... | |
  Public Member Functions inherited from LOCA::AnasaziOperator::AbstractStrategy | |
| AbstractStrategy () | |
| Constructor.  | |
| virtual | ~AbstractStrategy () | 
| Destructor.  | |
Protected Attributes | |
| Teuchos::RCP< LOCA::GlobalData > | globalData | 
| Global data.  | |
| std::string | myLabel | 
| Name of this operator.  | |
| 
Teuchos::RCP < Teuchos::ParameterList >  | eigenParams | 
| Stores parameters relating to the operator.  | |
| 
Teuchos::RCP < Teuchos::ParameterList >  | solverParams | 
| Stores linear solver parameters.  | |
| 
Teuchos::RCP < LOCA::TimeDependent::AbstractGroup >  | grp | 
| Stores group representing Jacobian and Mass matrix.  | |
| 
Teuchos::RCP < NOX::Abstract::MultiVector >  | tmp_r | 
| Stores a temporary vector for computing Rayleigh quotients.  | |
| 
Teuchos::RCP < NOX::Abstract::MultiVector >  | tmp_i | 
| Stores a temporary vector for computing Rayleigh quotients.  | |
| double | sigma | 
Stores Cayley pole  .  | |
| double | mu | 
Stores Cayley zero  .  | |
Anasazi operator for computing generalized eigenvalues using Cayley transformations.
This class implements the LOCA::AnasaziOperator::AbstractStrategy interface for computing generalized eigenvalues 
 and eigenvectors 
 of the system 
 where 
 is the Jacobian matrix and 
 is the mass matrix. The eigenvalues are computed using a Cayley transformation, i.e. solving 
 where 
 is the Cayley pole and 
 is the Cayley zero.
The parameters used by this class supplied in the constructor are:
 as defined above (Default 0.0) 
 as defined above (Default 0.0) Also the grp argument to the constructor must be a child of LOCA::TimeDependent::AbstractGroup for the shift-invert operations.
| LOCA::AnasaziOperator::Cayley2Matrix::Cayley2Matrix | ( | const Teuchos::RCP< LOCA::GlobalData > & | global_data, | 
| const Teuchos::RCP< LOCA::Parameter::SublistParser > & | topParams, | ||
| const Teuchos::RCP< Teuchos::ParameterList > & | eigenParams, | ||
| const Teuchos::RCP< Teuchos::ParameterList > & | solverParams, | ||
| const Teuchos::RCP< LOCA::TimeDependent::AbstractGroup > & | grp | ||
| ) | 
Constructor.
Argument grp must be of type LOCA::TimeDependent::AbstractGroup. See class description for a list of eigenParams.
References eigenParams, Teuchos::ParameterList::get(), mu, and sigma.
      
  | 
  virtual | 
Apply the operator.
Applies the inverse of the shifted operator, i.e., solves
 for 
, where 
 and 
. 
Implements LOCA::AnasaziOperator::AbstractStrategy.
References NOX::Abstract::MultiVector::clone(), NOX::Abstract::MultiVector::numVectors(), NOX::Abstract::Group::Ok, and NOX::ShapeCopy.
      
  | 
  virtual | 
Begin PostProcessing of eigenvalues.
Compute Jacobian and mass matrix once, for use in subsequent repeated calls to rayleighQuotient
Reimplemented from LOCA::AnasaziOperator::AbstractStrategy.
      
  | 
  virtual | 
PreProcess the random seed vector.
Performs one backward Euler iteration on the random initial seed vector, to satisfy contraints
Reimplemented from LOCA::AnasaziOperator::AbstractStrategy.
References NOX::Abstract::MultiVector::clone(), NOX::Abstract::MultiVector::numVectors(), NOX::Abstract::Group::Ok, and NOX::ShapeCopy.
      
  | 
  virtual | 
Compute Rayleigh quotient.
Computes the Rayleigh quotient 
 for the eigenvector 
. 
Implements LOCA::AnasaziOperator::AbstractStrategy.
References NOX::Abstract::Vector::createMultiVector(), NOX::Abstract::Vector::innerProduct(), NOX::Abstract::Group::Ok, NOX::Abstract::Vector::scale(), and NOX::ShapeCopy.
      
  | 
  virtual | 
Transform eigenvalue.
Transforms the given eigenvalue to the eigenvalue of the Jacobian-mass matrix system by shifting and inverting it.
Implements LOCA::AnasaziOperator::AbstractStrategy.
 1.8.5