| 
    NOX
    Development
    
   | 
 
Interface to underlying groups for Hopf point calculations using the Moore-Spence formulation. More...
#include <LOCA_Hopf_MooreSpence_AbstractGroup.H>


Public Member Functions | |
| AbstractGroup () | |
| Default constructor.  | |
| virtual | ~AbstractGroup () | 
| Destructor.  | |
Pure virtual methods  | |
These methods must be defined by any concrete implementation  | |
| virtual bool | isComplex () const =0 | 
Is   valid.  | |
| virtual  NOX::Abstract::Group::ReturnType  | computeComplex (double frequency)=0 | 
Compute  .  More... | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyComplex (const NOX::Abstract::Vector &input_real, const NOX::Abstract::Vector &input_imag, NOX::Abstract::Vector &result_real, NOX::Abstract::Vector &result_imag) const =0 | 
Compute  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyComplexMultiVector (const NOX::Abstract::MultiVector &input_real, const NOX::Abstract::MultiVector &input_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag) const =0 | 
Compute  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyComplexInverseMultiVector (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input_real, const NOX::Abstract::MultiVector &input_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag) const =0 | 
Solve  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDCeDp (const std::vector< int > ¶mIDs, const NOX::Abstract::Vector &yVector, const NOX::Abstract::Vector &zVector, double w, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag, bool isValid)=0 | 
Computes the derivative   where   is the parameter indexed by paramIDs.  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDCeDxa (const NOX::Abstract::Vector &yVector, const NOX::Abstract::Vector &zVector, double w, const NOX::Abstract::MultiVector &aVector, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag)=0 | 
Computes the directional derivative   for the given direction  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDCeDxa (const NOX::Abstract::Vector &yVector, const NOX::Abstract::Vector &zVector, double w, const NOX::Abstract::MultiVector &aVector, const NOX::Abstract::Vector &Ce_real, const NOX::Abstract::Vector &Ce_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag)=0 | 
Computes the directional derivative   for the given direction  . The arguments Ce_real and Ce_imag hold the real and imaginary components of  .  | |
  Public Member Functions inherited from LOCA::TurningPoint::MooreSpence::AbstractGroup | |
| AbstractGroup () | |
| Default constructor.  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDJnDpMulti (const std::vector< int > ¶mIDs, const NOX::Abstract::Vector &nullVector, NOX::Abstract::MultiVector &result, bool isValid)=0 | 
Computes the derivative  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDJnDxaMulti (const NOX::Abstract::Vector &nullVector, const NOX::Abstract::MultiVector &aVector, NOX::Abstract::MultiVector &result)=0 | 
Computes the directional derivative   for the given direction  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDJnDxaMulti (const NOX::Abstract::Vector &nullVector, const NOX::Abstract::Vector &JnVector, const NOX::Abstract::MultiVector &aVector, NOX::Abstract::MultiVector &result)=0 | 
Computes the directional derivative   for the given direction  .  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeDwtJnDxMulti (const NOX::Abstract::MultiVector &w, const NOX::Abstract::Vector &nullVector, NOX::Abstract::MultiVector &result)=0 | 
Computes the derivative  .  | |
  Public Member Functions inherited from LOCA::MultiContinuation::AbstractGroup | |
| AbstractGroup () | |
| Default constructor.  | |
| virtual void | copy (const NOX::Abstract::Group &source)=0 | 
| Copy the group (replaces operator = )  | |
| virtual void | setParamsMulti (const std::vector< int > ¶mIDs, const NOX::Abstract::MultiVector::DenseMatrix &vals)=0 | 
| Set parameters indexed by (integer) paramIDs.  | |
| virtual void | setParams (const LOCA::ParameterVector &p)=0 | 
| Set the parameter vector in the group to p (pVector = p).  | |
| virtual void | setParam (int paramID, double val)=0 | 
| Set parameter indexed by (integer) paramID.  | |
| virtual void | setParam (std::string paramID, double val)=0 | 
| Set parameter indexed by (std::string) paramID.  | |
| 
virtual const  LOCA::ParameterVector &  | getParams () const =0 | 
| Return a const reference to the ParameterVector owned by the group.  | |
| virtual double | getParam (int paramID) const =0 | 
| Return copy of parameter indexed by (integer) paramID.  | |
| virtual double | getParam (std::string paramID) const =0 | 
| Return copy of parameter indexed by (std::string) paramID.  | |
| virtual  NOX::Abstract::Group::ReturnType  | computeDfDpMulti (const std::vector< int > ¶mIDs, NOX::Abstract::MultiVector &dfdp, bool isValidF)=0 | 
| virtual void | preProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) | 
| Perform any preprocessing before a continuation step starts.  More... | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) | 
| Perform any postprocessing after a continuation step finishes.  More... | |
| virtual void | projectToDraw (const NOX::Abstract::Vector &x, double *px) const | 
| Projects solution to a few scalars for multiparameter continuation.  More... | |
| virtual int | projectToDrawDimension () const | 
| Returns the dimension of the project to draw array.  More... | |
| virtual double | computeScaledDotProduct (const NOX::Abstract::Vector &a, const NOX::Abstract::Vector &b) const | 
| Compute a scaled dot product.  More... | |
| virtual void | printSolution (const double) const | 
| Function to print out solution and parameter after successful step.  More... | |
| virtual void | printSolution (const NOX::Abstract::Vector &, const double) const | 
| Function to print out a vector and parameter after successful step.  More... | |
| virtual void | scaleVector (NOX::Abstract::Vector &x) const | 
| Scales a vector using scaling vector.  More... | |
  Public Member Functions inherited from NOX::Abstract::Group | |
| Group () | |
| Constructor.  More... | |
| virtual | ~Group () | 
| Destructor.  | |
| virtual NOX::Abstract::Group & | operator= (const NOX::Abstract::Group &source)=0 | 
| Copies the source group into this group.  More... | |
| virtual void | setX (const NOX::Abstract::Vector &y)=0 | 
| Set the solution vector x to y.  More... | |
| virtual void | computeX (const NOX::Abstract::Group &grp, const NOX::Abstract::Vector &d, double step)=0 | 
| Compute x = grp.x + step * d.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | computeF ()=0 | 
| Compute and store F(x).  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | computeJacobian () | 
| Compute and store Jacobian.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | computeGradient () | 
| Compute and store gradient.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | computeNewton (Teuchos::ParameterList ¶ms) | 
| Compute the Newton direction, using parameters for the linear solve.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobian (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const | 
| Applies Jacobian to the given input vector and puts the answer in the result.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const | 
| Applies Jacobian-Transpose to the given input vector and puts the answer in the result.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobianInverse (Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const | 
| Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const | 
| Apply right preconditiong to the given input vector.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobianMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const | 
| applyJacobian for multiple right-hand sides  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobianTransposeMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const | 
| applyJacobianTranspose for multiple right-hand sides  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyJacobianInverseMultiVector (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const | 
| applyJacobianInverse for multiple right-hand sides  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applyRightPreconditioningMultiVector (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const | 
| applyRightPreconditioning for multiple right-hand sides  More... | |
| virtual bool | isF () const =0 | 
| Return true if F is valid.  | |
| virtual bool | isJacobian () const | 
| Return true if the Jacobian is valid.  More... | |
| virtual bool | isGradient () const | 
| Return true if the gradient is valid.  More... | |
| virtual bool | isNewton () const | 
| Return true if the Newton direction is valid.  More... | |
| 
virtual const  NOX::Abstract::Vector &  | getX () const =0 | 
| Return solution vector.  | |
| virtual const  NOX::Abstract::Vector &  | getScaledX () const | 
| 
virtual const  NOX::Abstract::Vector &  | getF () const =0 | 
| Return F(x)  | |
| virtual double | getNormF () const =0 | 
| Return 2-norm of F(x).  More... | |
| 
virtual const  NOX::Abstract::Vector &  | getGradient () const =0 | 
| Return gradient.  | |
| 
virtual const  NOX::Abstract::Vector &  | getNewton () const =0 | 
| Return Newton direction.  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::Vector >  | getXPtr () const =0 | 
| Return RCP to solution vector.  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::Vector >  | getFPtr () const =0 | 
| Return RCP to F(x)  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::Vector >  | getGradientPtr () const =0 | 
| Return RCP to gradient.  | |
| 
virtual Teuchos::RCP< const  NOX::Abstract::Vector >  | getNewtonPtr () const =0 | 
| Return RCP to Newton direction.  | |
| virtual void | logLastLinearSolveStats (NOX::SolverStats &stats) const | 
| Adds statistics from last linear solve to the SovlerStats object.  | |
| virtual  NOX::Abstract::Group::ReturnType  | getNormLastLinearSolveResidual (double &residual) const | 
| Return the norm of the last linear solve residual as the result of either a call to computeNewton() or applyJacobianInverse().  More... | |
| virtual Teuchos::RCP < NOX::Abstract::Group >  | clone (NOX::CopyType type=NOX::DeepCopy) const =0 | 
| Create a new Group of the same derived type as this one by cloning this one, and return a ref count pointer to the new group.  More... | |
  Public Member Functions inherited from LOCA::TimeDependent::AbstractGroup | |
| AbstractGroup () | |
| Default constructor.  | |
| virtual  NOX::Abstract::Group::ReturnType  | computeSecondShiftedMatrix (double alpha, double beta)=0 | 
| Compute the second shifted matrix. Can avoid recomputing if two are stored.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applySecondShiftedMatrix (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const =0 | 
| Multiply the shifted matrix by a vector.  More... | |
| virtual  NOX::Abstract::Group::ReturnType  | applySecondShiftedMatrixMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const =0 | 
| Multiply the shifted matrix by a multi-vector.  More... | |
| 
virtual  NOX::Abstract::Group::ReturnType  | computeShiftedMatrix (double alpha, double beta)=0 | 
| Compute the shifted matrix.  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyShiftedMatrix (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const =0 | 
| Multiply the shifted matrix by a vector.  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyShiftedMatrixMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const =0 | 
| Multiply the shifted matrix by a multi-vector.  | |
| 
virtual  NOX::Abstract::Group::ReturnType  | applyShiftedMatrixInverseMultiVector (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const =0 | 
| Apply the inverse of the shifted matrix by a multi-vector, as needed by the shift-and-invert and generalized Cayley transformations.  | |
Additional Inherited Members | |
  Public Types inherited from NOX::Abstract::Group | |
| enum | ReturnType {  Ok, NotDefined, BadDependency, NotConverged, Failed }  | 
| The computation of, say, the Newton direction in computeNewton() may fail in many different ways, so we have included a variety of return codes to describe the failures. Of course, we also have a code for success.  More... | |
Interface to underlying groups for Hopf point calculations using the Moore-Spence formulation.
This abstract class provides the required interface for underlying groups to locate Hopf bifurcations using the bordering algorithm for the Moore-Spence Hopf fomulation (see LOCA::Hopf::MooreSpence::ExtendedGroup for a description of the governing equations).
This class is derived from the LOCA::TurningPoint::MooreSpence::AbstractGroup and LOCa::TimeDependent::AbstractGroup and declares several pure virtual methods for applying, solving, and computing derivatives of the complex matrix 
 where 
 is the Jacobian matrix, 
 is the mass matrix, and 
 is a (real) scalar. 
      
  | 
  pure virtual | 
Compute 
. 
The argument frequency stores 
. 
Implemented in LOCA::Epetra::Group, LOCA::LAPACK::Group, and LOCA::Abstract::Group.
Referenced by LOCA::DerivUtils::computeDCeDp(), LOCA::DerivUtils::computeDCeDxa(), LOCA::DerivUtils::computeDwtCeDp(), and LOCA::DerivUtils::computeDwtCeDx().
 1.8.5