NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
LOCA::Epetra::Group Class Reference

Extension of the NOX::Epetra::Group to LOCA. More...

#include <LOCA_Epetra_Group.H>

Inheritance diagram for LOCA::Epetra::Group:
Inheritance graph
[legend]
Collaboration diagram for LOCA::Epetra::Group:
Collaboration graph
[legend]

Public Member Functions

 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data, Teuchos::ParameterList &printingParams, const Teuchos::RCP< LOCA::Epetra::Interface::Required > &i, NOX::Epetra::Vector &initialGuess, const LOCA::ParameterVector &p)
 Constructor with NO linear system (VERY LIMITED). More...
 
 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data, Teuchos::ParameterList &printingParams, const Teuchos::RCP< LOCA::Epetra::Interface::Required > &i, NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< NOX::Epetra::LinearSystem > &linSys, const LOCA::ParameterVector &p)
 Standard Constructor enabling most LOCA support. More...
 
 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data, Teuchos::ParameterList &printingParams, const Teuchos::RCP< LOCA::Epetra::Interface::TimeDependent > &i, NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< NOX::Epetra::LinearSystem > &linSys, const Teuchos::RCP< NOX::Epetra::LinearSystem > &shiftedLinSys, const LOCA::ParameterVector &p)
 Constructor with time-dependent interface and shifted linear system. More...
 
 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data, Teuchos::ParameterList &printingParams, const Teuchos::RCP< LOCA::Epetra::Interface::TimeDependentMatrixFree > &i, NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< NOX::Epetra::LinearSystem > &linSys, const Teuchos::RCP< NOX::Epetra::LinearSystem > &shiftedLinSys, const LOCA::ParameterVector &p)
 Constructor with time-dependent matrix-free interface and shifted linear system. More...
 
 Group (const Group &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor. If type is DeepCopy, takes ownership of valid shared Jacobian and shared preconditioning matrix.
 
virtual ~Group ()
 Destructor.
 
virtual Groupoperator= (const Group &source)
 Assignment operator.
 
virtual void setFreeEnergyInterface (const Teuchos::RCP< LOCA::Epetra::Interface::FreeEnergy > &iFE)
 Method to inject an interface for calucatiuong the free energy.
 
void declareSeparateMatrixMemory (bool separateMem=true)
 Method for calling code to guarantee to LOCA that separate matrix.
 
virtual Teuchos::RCP
< NOX::Epetra::Interface::Required
getUserInterface ()
 Return the userInterface.
 
virtual void printSolution (const NOX::Epetra::Vector &x, const double conParam) const
 Call the user interface print() routine, any vector.
 
void setScaleVector (const NOX::Abstract::Vector &s)
 Sets the scale vector.
 
void setJacobianOperatorForSolve (const Teuchos::RCP< const Epetra_Operator > &op) const
 Sets the Jacobian operator.
 
virtual Teuchos::RCP< const
NOX::Epetra::LinearSystem
getComplexLinearSystem () const
 Return the Linear System.
 
virtual Teuchos::RCP
< NOX::Epetra::LinearSystem
getComplexLinearSystem ()
 Return the Linear System.
 
virtual void getComplexMaps (Teuchos::RCP< const Epetra_BlockMap > &baseMap, Teuchos::RCP< const Epetra_BlockMap > &globalMap) const
 
Overloaded NOX::Epetra::Group methods.
virtual NOX::Abstract::Groupoperator= (const NOX::Abstract::Group &source)
 Assignment operator.
 
virtual NOX::Abstract::Groupoperator= (const NOX::Epetra::Group &source)
 Assignment operator.
 
virtual Teuchos::RCP
< NOX::Abstract::Group
clone (NOX::CopyType type=NOX::DeepCopy) const
 Cloning function.
 
virtual
NOX::Abstract::Group::ReturnType 
computeF ()
 Overloaded computeF() More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeJacobian ()
 Overloaded computeJacobian() More...
 
Implementation of LOCA::Abstract::TransposeSolveGroup methods.
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianTransposeInverse (Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Solve Jacobian-tranpose system. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianTransposeInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Solve Jacobian-tranpose system with multiple right-hand sides. More...
 
Implementation of LOCA::MultiContinuation::AbstractGroup virtual methods.
virtual void copy (const NOX::Abstract::Group &source)
 Copy.
 
virtual void setParams (const ParameterVector &p)
 Set the parameters.
 
virtual void setParam (int paramID, double val)
 Set parameter indexed by paramID.
 
virtual void setParam (std::string paramID, double val)
 Set parameter indexed by paramID.
 
const LOCA::ParameterVectorgetParams () const
 Return a const reference to the ParameterVector owned by the group.
 
virtual double getParam (int paramID) const
 Return copy of parameter indexed by paramID.
 
virtual double getParam (std::string paramID) const
 Return copy of parameter indexed by paramID.
 
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 conParam) const
 Call the user interface print() routine, solution vector.
 
virtual void printSolution (const NOX::Abstract::Vector &x, const double conParam) const
 Call the user interface print() routine, any vector.
 
virtual void scaleVector (NOX::Abstract::Vector &x) const
 Scales a vector using scaling vector. More...
 
Implementation of LOCA::Homotopy::AbstractGroup virtual methods.
virtual
NOX::Abstract::Group::ReturnType 
augmentJacobianForHomotopy (double a, double b)
 Replace Jacobian $ J$ by $ aJ+bI$ where $ I$ is the identity matrix.
 
Implementation of LOCA::TimeDependent::AbstractGroup virtual methods.
virtual
NOX::Abstract::Group::ReturnType 
computeShiftedMatrix (double alpha, double beta)
 Compute the shifted matrix.
 
virtual
NOX::Abstract::Group::ReturnType 
applyShiftedMatrix (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Multiply the shifted matrix by a vector.
 
virtual
NOX::Abstract::Group::ReturnType 
applyShiftedMatrixMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Multiply the shifted matrix by a multi-vector.
 
virtual
NOX::Abstract::Group::ReturnType 
applyShiftedMatrixInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Apply the inverse of the shifted matrix by a multi-vector, as needed by the shift-and-invert and generalized Cayley transformations.
 
virtual
NOX::Abstract::Group::ReturnType 
computeSecondShiftedMatrix (double alpha, double beta)
 Compute the second shifted matrix (uses different memory then Shifted matrix)
 
virtual
NOX::Abstract::Group::ReturnType 
applySecondShiftedMatrix (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Multiply the second shifted matrix by a vector.
 
virtual
NOX::Abstract::Group::ReturnType 
applySecondShiftedMatrixMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Multiply the second shifted matrix by a multi-vector.
 
Implementation of LOCA::Hopf::MooreSpence::AbstractGroup virtual methods.
virtual bool isComplex () const
 Is $ J+i\omega B$ valid.
 
virtual
NOX::Abstract::Group::ReturnType 
computeComplex (double frequency)
 Compute $ J+i\omega B$. 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
 Compute $(J+i\omega B)(y+iz)$.
 
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
 Compute $(J+i\omega B)(y+iz)$.
 
virtual
NOX::Abstract::Group::ReturnType 
applyComplexInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input_real, const NOX::Abstract::MultiVector &input_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag) const
 Solve $(J+i\omega B)(y+iz) = a+ib$.
 
Implementation of LOCA::Hopf::MinimallyAugmented::AbstractGroup virtual methods.
virtual
NOX::Abstract::Group::ReturnType 
applyComplexTranspose (const NOX::Abstract::Vector &input_real, const NOX::Abstract::Vector &input_imag, NOX::Abstract::Vector &result_real, NOX::Abstract::Vector &result_imag) const
 
virtual
NOX::Abstract::Group::ReturnType 
applyComplexTransposeMultiVector (const NOX::Abstract::MultiVector &input_real, const NOX::Abstract::MultiVector &input_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag) const
 
virtual
NOX::Abstract::Group::ReturnType 
applyComplexTransposeInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input_real, const NOX::Abstract::MultiVector &input_imag, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag) const
 Solve $(J+i\omega B)^H (x + iy) = a+ib$.
 
Implementation of LOCA::PhseTransition::AbstractGroup virtual methods.
virtual double computeFreeEnergy ()
 Computes the free energy at the current solution and parameter values.
 
- Public Member Functions inherited from NOX::Epetra::Group
 Group (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess)
 Constructor with NO linear system (VERY LIMITED). More...
 
 Group (Teuchos::ParameterList &printingParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const NOX::Epetra::Vector &initialGuess, const Teuchos::RCP< NOX::Epetra::LinearSystem > &linSys)
 Standard Constructor.
 
 Group (const NOX::Epetra::Group &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor. If type is DeepCopy, takes ownership of valid shared linear system.
 
virtual Teuchos::RCP
< NOX::Epetra::Interface::Required
getRequiredInterface ()
 Return the userInterface.
 
virtual Teuchos::RCP< const
NOX::Epetra::LinearSystem
getLinearSystem () const
 Return the Linear System.
 
virtual Teuchos::RCP
< NOX::Epetra::LinearSystem
getLinearSystem ()
 Return the Linear System.
 
virtual
NOX::Abstract::Group::ReturnType 
computeJacobianConditionNumber (int maxIters, double tolerance, int krylovSubspaceSize=100, bool printOutput=false)
 
virtual double getJacobianConditionNumber () const
 Returns the condition number of the Jacobian matrix.
 
virtual void disableLinearResidualComputation (const bool disableChoice)
 Sets option to disable linear resid computation. If disabled, this saves on a MatVec per Newton but disallows inexact Newton methods.
 
virtual void setX (const NOX::Epetra::Vector &y)
 
virtual void setX (const NOX::Abstract::Vector &y)
 Set the solution vector x to y. More...
 
virtual void computeX (const Group &grp, const NOX::Epetra::Vector &d, double step)
 
virtual void computeX (const NOX::Abstract::Group &grp, const NOX::Abstract::Vector &d, double step)
 Compute x = grp.x + step * d. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeGradient ()
 Compute and store gradient. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeNewton (Teuchos::ParameterList &params)
 Compute the Newton direction, using parameters for the linear solve. More...
 
virtual bool isF () const
 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 bool isNormNewtonSolveResidual () const
 Returns true if the value of the Norm of the linear model for a full Newton step ||Js + f|| is valid with respect to the current solution vector.
 
virtual bool isPreconditioner () const
 Returns true if an explicitly constructed preconditioner exists (i.e. one that is computed and saved for further use in multiple calls to applyRightPreconditioner).
 
virtual bool isConditionNumber () const
 Returns true if the condition number has been computed.
 
virtual const
NOX::Abstract::Vector
getX () const
 Return solution vector.
 
virtual const
NOX::Abstract::Vector
getF () const
 Return F(x)
 
virtual double getNormF () const
 Return 2-norm of F(x). More...
 
virtual const
NOX::Abstract::Vector
getGradient () const
 Return gradient.
 
virtual const
NOX::Abstract::Vector
getNewton () const
 Return Newton direction.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getXPtr () const
 Return RCP to solution vector.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getFPtr () const
 Return RCP to F(x)
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getGradientPtr () const
 Return RCP to gradient.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getNewtonPtr () const
 Return RCP to Newton direction.
 
virtual
NOX::Abstract::Group::ReturnType 
getNormLastLinearSolveResidual (double &residual) const
 Returns the 2-norm of the residual of the linear model used in the Newton solve computation, ||Js+f||. This does not account for line search adjustments to the step length!
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobian (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const
 
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::Epetra::Vector &input, NOX::Epetra::Vector &result) const
 
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 
applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList &params, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const
 
virtual
NOX::Abstract::Group::ReturnType 
applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Apply right preconditiong to the given input vector. More...
 
- Public Member Functions inherited from NOX::Abstract::Group
 Group ()
 Constructor. More...
 
virtual Teuchos::RCP
< NOX::Abstract::Group
getNestedGroup ()
 Return an internally stored group from this group. More...
 
virtual Teuchos::RCP< const
NOX::Abstract::Group
getNestedGroup () const
 Return an internally stored group from this group. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianInverse (Teuchos::ParameterList &params, 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 
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 &params, 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 &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 applyRightPreconditioning for multiple right-hand sides More...
 
virtual const
NOX::Abstract::Vector
getScaledX () const
 
- Public Member Functions inherited from LOCA::Abstract::Group
 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data)
 Constructor.
 
 Group (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::DerivUtils > &deriv)
 Constructor.
 
 Group (const Group &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual void setParamsMulti (const std::vector< int > &paramIDs, const NOX::Abstract::MultiVector::DenseMatrix &vals)
 Set parameters indexed by (integer) paramIDs.
 
virtual void notifyCompletedStep ()
 
- Public Member Functions inherited from LOCA::Homotopy::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
virtual ~AbstractGroup ()
 Destructor.
 
- Public Member Functions inherited from LOCA::MultiContinuation::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
- Public Member Functions inherited from LOCA::TurningPoint::MinimallyAugmented::FiniteDifferenceGroup
 FiniteDifferenceGroup ()
 Constructor.
 
 FiniteDifferenceGroup (const FiniteDifferenceGroup &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual ~FiniteDifferenceGroup ()
 Destructor.
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtJnDp (const std::vector< int > &paramIDs, const NOX::Abstract::Vector &w, const NOX::Abstract::Vector &nullVector, NOX::Abstract::MultiVector::DenseMatrix &result, bool isValid)
 Computes the derivative $\partial w^TJn/\partial p$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtJDp (const std::vector< int > &paramIDs, const NOX::Abstract::Vector &w, NOX::Abstract::MultiVector &result, bool isValid)
 Computes the derivative $\partial w^TJ/\partial p$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtJnDx (const NOX::Abstract::Vector &w, const NOX::Abstract::Vector &nullVector, NOX::Abstract::Vector &result)
 Computes the derivative $\frac{\partial w^TJn}{\partial x}$. More...
 
- Public Member Functions inherited from LOCA::TurningPoint::MinimallyAugmented::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
- Public Member Functions inherited from LOCA::TurningPoint::MooreSpence::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
- Public Member Functions inherited from LOCA::TurningPoint::MooreSpence::FiniteDifferenceGroup
 FiniteDifferenceGroup ()
 Constructor.
 
 FiniteDifferenceGroup (const FiniteDifferenceGroup &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual
NOX::Abstract::Group::ReturnType 
computeDJnDpMulti (const std::vector< int > &paramIDs, const NOX::Abstract::Vector &nullVector, NOX::Abstract::MultiVector &result, bool isValid)
 Computes the derivative $\partial Jn/\partial p$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeDJnDxaMulti (const NOX::Abstract::Vector &nullVector, const NOX::Abstract::MultiVector &aVector, NOX::Abstract::MultiVector &result)
 Computes the directional derivative $\frac{\partial Jn}{\partial x} a$ for the given direction $ a$. More...
 
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)
 Computes the directional derivative $\frac{\partial Jn}{\partial x} a$ for the given direction $ a$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtJnDxMulti (const NOX::Abstract::MultiVector &w, const NOX::Abstract::Vector &nullVector, NOX::Abstract::MultiVector &result)
 Computes the derivative $\frac{\partial w^TJn}{\partial x}$. More...
 
- Public Member Functions inherited from LOCA::MultiContinuation::FiniteDifferenceGroup
 FiniteDifferenceGroup ()
 Constructor.
 
 FiniteDifferenceGroup (const FiniteDifferenceGroup &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual void setDerivUtils (const Teuchos::RCP< LOCA::DerivUtils > &deriv)
 Set the LOCA::DerivUtils object.
 
virtual
NOX::Abstract::Group::ReturnType 
computeDfDpMulti (const std::vector< int > &paramIDs, NOX::Abstract::MultiVector &dfdp, bool isValidF)
 
- Public Member Functions inherited from LOCA::Pitchfork::MinimallyAugmented::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
virtual ~AbstractGroup ()
 Destructor.
 
- Public Member Functions inherited from LOCA::Pitchfork::MooreSpence::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
virtual double innerProduct (const NOX::Abstract::Vector &a, const NOX::Abstract::Vector &b) const
 Compute the inner product of a and b. More...
 
virtual void innerProduct (const NOX::Abstract::MultiVector &a, const NOX::Abstract::MultiVector &b, NOX::Abstract::MultiVector::DenseMatrix &c) const
 Compute the inner product of a and b. More...
 
- Public Member Functions inherited from LOCA::Hopf::MinimallyAugmented::FiniteDifferenceGroup
 FiniteDifferenceGroup ()
 Constructor.
 
 FiniteDifferenceGroup (const FiniteDifferenceGroup &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual ~FiniteDifferenceGroup ()
 Destructor.
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtCeDp (const std::vector< int > &paramIDs, const NOX::Abstract::Vector &w1, const NOX::Abstract::Vector &w2, const NOX::Abstract::Vector &y, const NOX::Abstract::Vector &x, double omega, NOX::Abstract::MultiVector::DenseMatrix &result_real, NOX::Abstract::MultiVector::DenseMatrix &result_imag, bool isValid)
 Computes the derivative $\partial w^TCe/\partial p$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeDwtCeDx (const NOX::Abstract::Vector &w1, const NOX::Abstract::Vector &w2, const NOX::Abstract::Vector &y, const NOX::Abstract::Vector &z, double omega, NOX::Abstract::Vector &result_real, NOX::Abstract::Vector &result_imag)
 Computes the derivative $\frac{\partial w^TCe}{\partial x}$. More...
 
- Public Member Functions inherited from LOCA::Hopf::MinimallyAugmented::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
virtual ~AbstractGroup ()
 Destructor.
 
- Public Member Functions inherited from LOCA::Hopf::MooreSpence::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
- Public Member Functions inherited from LOCA::TimeDependent::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
- Public Member Functions inherited from LOCA::Hopf::MooreSpence::FiniteDifferenceGroup
 FiniteDifferenceGroup ()
 Constructor.
 
 FiniteDifferenceGroup (const FiniteDifferenceGroup &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual
NOX::Abstract::Group::ReturnType 
computeDCeDp (const std::vector< int > &paramIDs, const NOX::Abstract::Vector &yVector, const NOX::Abstract::Vector &zVector, double w, NOX::Abstract::MultiVector &result_real, NOX::Abstract::MultiVector &result_imag, bool isValid)
 Computes the derivative $\frac{\partial (J+i\omega B)(y+iz)}{\partial p}$ where $ p$ is the parameter indexed by paramIDs. More...
 
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)
 Computes the directional derivative $\frac{\partial (J+i\omega B)(y+iz)}{\partial x} a$ for the given direction $ a$. More...
 
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)
 Computes the directional derivative $\frac{\partial (J+i\omega B)(y+iz)}{\partial x} a$ for the given direction $ a$. The arguments Ce_real and Ce_imag hold the real and imaginary components of $(J+i\omega B)(y+iz)$. More...
 
- Public Member Functions inherited from LOCA::PhaseTransition::AbstractGroup
 AbstractGroup ()
 Default constructor.
 
virtual ~AbstractGroup ()
 Destructor.
 
- Public Member Functions inherited from LOCA::Abstract::TransposeSolveGroup
 TransposeSolveGroup ()
 Constructor.
 
virtual ~TransposeSolveGroup ()
 Destructor.
 

Protected Member Functions

virtual void resetIsValid ()
 resets the isValid flags to false
 
- Protected Member Functions inherited from NOX::Epetra::Group
virtual void logLastLinearSolveStats (NOX::SolverStats &stats) const
 Adds statistics from last linear solve to the SovlerStats object.
 
virtual bool computeNormNewtonSolveResidual ()
 Computes the 2-norm of the residual of the linear model used in the Newton solve computation, ||Js+f||.
 

Protected Attributes

Teuchos::RCP< LOCA::GlobalDataglobalData
 Global data.
 
Teuchos::ParameterListprintParams
 Printing parameters.
 
LOCA::ParameterVector params
 Parameter vector.
 
Teuchos::RCP
< LOCA::Epetra::Interface::Required
userInterface
 Reference to the user supplied interface functions.
 
Teuchos::RCP
< LOCA::Epetra::Interface::TimeDependent
userInterfaceTime
 Interface for shifted-matrix.
 
Teuchos::RCP
< LOCA::Epetra::Interface::TimeDependentMatrixFree
userInterfaceTimeMF
 Interface for matrix-free shifted-matrix.
 
Teuchos::RCP
< LOCA::Epetra::Interface::FreeEnergy
userInterfaceFreeEnergy
 Interface for free enerfgy calculation for phase transitions.
 
Teuchos::RCP
< NOX::SharedObject
< NOX::Epetra::LinearSystem,
LOCA::Epetra::Group > > 
shiftedSharedLinearSystem
 Shared shifted linear system.
 
bool isValidShiftedPrec
 Is preconditioner for shifted matrix valid.
 
double alpha_
 $\alpha$ for matrix-free shifted matrix
 
double beta_
 $\beta$ for matrix-free shifted matrix
 
Teuchos::RCP< Epetra_VectortmpVectorPtr2
 Extra vector needed for intermediate calculations of LOCA routines. More...
 
Teuchos::RCP
< NOX::Abstract::Vector
scaleVecPtr
 Stores a pointer to the scale vector.
 
Teuchos::RCP
< LOCA::Epetra::TransposeLinearSystem::AbstractStrategy
tls_strategy
 Stores transpose linear solver strategy.
 
Teuchos::RCP
< NOX::SharedObject
< NOX::Epetra::LinearSystem,
LOCA::Epetra::Group > > 
complexSharedLinearSystem
 Shared complex system.
 
Teuchos::RCP
< EpetraExt::BlockCrsMatrix > 
complexMatrix
 Complex matrix.
 
Teuchos::RCP
< EpetraExt::BlockVector > 
complexVec
 Complex vector.
 
bool isValidComplex
 Is complex matrix valid.
 
bool isValidComplexPrec
 Is complex matrix preconditioner valid.
 
bool separateMatrixMemoryDeclared
 Is complex matrix valid.
 
- Protected Attributes inherited from NOX::Epetra::Group
const NOX::Utils utils
 Printing Utilities object.
 
double normNewtonSolveResidual
 2-Norm of the Newton solve residual: ||Js+f||
 
double conditionNumber
 condition number of Jacobian
 
Teuchos::RCP
< AztecOOConditionNumber > 
azConditionNumberPtr
 Pointer to the condition number object.
 
bool linearResidCompDisabled
 
Teuchos::RCP
< NOX::Epetra::Interface::Required
userInterfacePtr
 Reference to the user supplied interface functions.
 
bool linearSolveConverged
 
int numIterations
 
double achievedTol
 
Teuchos::RCP< NOX::Epetra::VectorxVectorPtr
 Solution vector pointer.
 
NOX::Epetra::VectorxVector
 Solution vector.
 
Teuchos::RCP< NOX::Epetra::VectorRHSVectorPtr
 Right-hand-side vector (function evaluation).
 
NOX::Epetra::VectorRHSVector
 Right-hand-side vector pointer (function evaluation).
 
Teuchos::RCP< NOX::Epetra::VectorgradVectorPtr
 Gradient vector pointer(steepest descent vector).
 
NOX::Epetra::VectorgradVector
 Gradient vector (steepest descent vector).
 
Teuchos::RCP< NOX::Epetra::VectorNewtonVectorPtr
 Newton direction vector pointer.
 
NOX::Epetra::VectorNewtonVector
 Newton direction vector.
 
Teuchos::RCP< Epetra_VectortmpVectorPtr
 An extra temporary vector, only allocated if needed.
 
bool isValidRHS
 
bool isValidJacobian
 
bool isValidGrad
 
bool isValidNewton
 
bool isValidNormNewtonSolveResidual
 
bool isValidPreconditioner
 
bool isValidSolverJacOp
 
bool isValidConditionNumber
 
Teuchos::RCP
< NOX::SharedObject
< NOX::Epetra::LinearSystem,
NOX::Epetra::Group > > 
sharedLinearSystemPtr
 Pointer to shared Jacobian matrix.
 
NOX::SharedObject
< NOX::Epetra::LinearSystem,
NOX::Epetra::Group > & 
sharedLinearSystem
 Reference to shared Jacobian matrix.
 
- Protected Attributes inherited from LOCA::Abstract::Group
Teuchos::RCP< LOCA::GlobalDataglobalData
 Global data.
 
- Protected Attributes inherited from LOCA::MultiContinuation::FiniteDifferenceGroup
Teuchos::RCP< DerivUtilsderivPtr
 Pointer to current DerivUtils derivative computation object.
 

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...
 
- Public Attributes inherited from NOX::Epetra::Group
virtual
NOX::Abstract::Group::ReturnType 
const
 Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result. More...
 

Detailed Description

Extension of the NOX::Epetra::Group to LOCA.

This class extends the NOX::Epetra::Group to LOCA enabling continuation and bifurcation capabilities using Epetra. It is derived from the NOX::Epetra::Group (basic Epetra support), the LOCA::Abstract::Group (brings in all LOCA abstract base classes), and the LOCA::Abstract::TransposeSolveGroup (for applyJacobianTransposeInverse() methods). It stores a parameter vector for setting/retrieving parameter values and overloads the computeF() and computeJacobian() methods of the NOX::Epetra::Group parent class to set the entire contents of the parameter vector in the problem interface before calling the NOX::Epetra::Group computeF() and computeJacobian().

Since it is derived from the LOCA::Abstract::Group (which is in turn derived from all FiniteDifference groups), it uses the finite-difference implementations for all parameter derivatives and second derivatives. However this behavior can be modified by calling the setDerivUtils() method of the LOCA::MultiContinuation::FiniteDifferenceGroup parent class.

This class provides complete support for all continuation and bifurcation methods including shift-invert and Cayley methods for computing eigenvalues and Hopf bifurcations. However this support is only enabled by calling the appropriate constructor described below.

Constructor & Destructor Documentation

LOCA::Epetra::Group::Group ( const Teuchos::RCP< LOCA::GlobalData > &  global_data,
Teuchos::ParameterList printingParams,
const Teuchos::RCP< LOCA::Epetra::Interface::Required > &  i,
NOX::Epetra::Vector initialGuess,
const LOCA::ParameterVector p 
)

Constructor with NO linear system (VERY LIMITED).

WARNING: If this constructor is used, then methods that require a Jacobian or preconditioning will not be available. You will be limited to simple algorithms like nonlinear-CG with no preconditioning.

LOCA::Epetra::Group::Group ( const Teuchos::RCP< LOCA::GlobalData > &  global_data,
Teuchos::ParameterList printingParams,
const Teuchos::RCP< LOCA::Epetra::Interface::Required > &  i,
NOX::Epetra::Vector initialGuess,
const Teuchos::RCP< NOX::Epetra::LinearSystem > &  linSys,
const LOCA::ParameterVector p 
)

Standard Constructor enabling most LOCA support.

This is the most commonly used constructor and provides support for all LOCA algorithms except shift-invert and Cayley transformations and Hopf bifurcations.

LOCA::Epetra::Group::Group ( const Teuchos::RCP< LOCA::GlobalData > &  global_data,
Teuchos::ParameterList printingParams,
const Teuchos::RCP< LOCA::Epetra::Interface::TimeDependent > &  i,
NOX::Epetra::Vector initialGuess,
const Teuchos::RCP< NOX::Epetra::LinearSystem > &  linSys,
const Teuchos::RCP< NOX::Epetra::LinearSystem > &  shiftedLinSys,
const LOCA::ParameterVector p 
)

Constructor with time-dependent interface and shifted linear system.

Use this constructor to enable shift-invert and Cayley transformations or Hopf bifurcations. It requires another interface to compute the shifted matrix $\alpha J + \beta M$ where $ J$ is the Jacobian matrix and $ M$ is the mass matrix, and a linear system object to solve this system. Setting linSys = shiftedLinSys is a valid option for passing the shifted solver, but this will cause the shifted matrix to overwrite the Jacobian possibly resulting in more matrix fills. See declareSeparateMatrixMemory() method below to assert separate memory.

References Teuchos::rcp(), and shiftedSharedLinearSystem.

Constructor with time-dependent matrix-free interface and shifted linear system.

This constructor may also be used for shift-invert and Cayley transformations, but should be only be used for a matrix-free method for solving the shifted system.

References Teuchos::rcp(), and shiftedSharedLinearSystem.

Member Function Documentation

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::applyComplexTranspose ( const NOX::Abstract::Vector input_real,
const NOX::Abstract::Vector input_imag,
NOX::Abstract::Vector result_real,
NOX::Abstract::Vector result_imag 
) const
virtual

Computes conjugate-tranpose matrix vector product $ (J+i\omega B)^H (x + iy) $.

Reimplemented from LOCA::Abstract::Group.

References NOX::Abstract::Group::BadDependency, NOX::Epetra::Vector::getEpetraVector(), and NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::applyComplexTransposeMultiVector ( const NOX::Abstract::MultiVector input_real,
const NOX::Abstract::MultiVector input_imag,
NOX::Abstract::MultiVector result_real,
NOX::Abstract::MultiVector result_imag 
) const
virtual

Computes conjugate-tranpose matrix vector product $ (J+i\omega B)^H (x + iy) $.

Reimplemented from LOCA::Abstract::Group.

References NOX::Abstract::MultiVector::numVectors(), and NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::applyJacobianTransposeInverse ( Teuchos::ParameterList params,
const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Solve Jacobian-tranpose system.

In addition to all regular linear solver parameters, this method references the following additional parameters:

  • "Transpose Solver Method" – [string] (default: "Transpose Preconditioner") Method for preconditioning the transpose linear system (LOCA::Epetra::TransposeLinearSystem::Factory). Available choices are:
    • "Transpose Preconditioner" – Use the transpose of the preconditioner for the original system.
    • "Left Preconditioning" – Use the transpose of the preconditioner, and apply using left preconditioning.
    • "Explicit Transpose" – Form the transpose of the matrix and compute the preconditioner. This method is available only if Trilinos is configured with EpetraExt support (–enable-epetraext).

Implements LOCA::Abstract::TransposeSolveGroup.

References LOCA::Epetra::TransposeLinearSystem::Factory::create(), LOCA::Epetra::TransposeLinearSystem::AbstractStrategy::createJacobianTranspose(), NOX::Epetra::LinearSystem::destroyPreconditioner(), NOX::Epetra::LinearSystem::getJacobianOperator(), NOX::Abstract::Group::NotConverged, NOX::Abstract::Group::Ok, Teuchos::rcp(), NOX::Epetra::LinearSystem::setJacobianOperatorForSolve(), and Epetra_Operator::SetUseTranspose().

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::applyJacobianTransposeInverseMultiVector ( Teuchos::ParameterList params,
const NOX::Abstract::MultiVector input,
NOX::Abstract::MultiVector result 
) const
virtual

Solve Jacobian-tranpose system with multiple right-hand sides.

In addition to all regular linear solver parameters, this method references the following additional parameters:

  • "Transpose Solver Method" – [string] (default: "Transpose Preconditioner") Method for preconditioning the transpose linear system (LOCA::Epetra::TransposeLinearSystem::Factory). Available choices are:
    • "Transpose Preconditioner" – Use the transpose of the preconditioner for the original system.
    • "Left Preconditioning" – Use the transpose of the preconditioner, and apply using left preconditioning.
    • "Explicit Transpose" – Form the transpose of the matrix and compute the preconditioner. This method is available only if Trilinos is configured with EpetraExt support (–enable-epetraext).

Implements LOCA::Abstract::TransposeSolveGroup.

References LOCA::Epetra::TransposeLinearSystem::Factory::create(), LOCA::Epetra::TransposeLinearSystem::AbstractStrategy::createJacobianTranspose(), NOX::Epetra::LinearSystem::destroyPreconditioner(), NOX::Epetra::LinearSystem::getJacobianOperator(), NOX::Abstract::Group::NotConverged, NOX::Abstract::MultiVector::numVectors(), NOX::Abstract::Group::Ok, Teuchos::rcp(), NOX::Epetra::LinearSystem::setJacobianOperatorForSolve(), and Epetra_Operator::SetUseTranspose().

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::computeComplex ( double  frequency)
virtual
NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::computeF ( )
virtual

Overloaded computeF()

Calls LOCA::Epetra::Interface::setParams before evalulating F.

Reimplemented from NOX::Epetra::Group.

References NOX::Epetra::Group::computeF(), and NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType LOCA::Epetra::Group::computeJacobian ( )
virtual

Overloaded computeJacobian()

Calls LOCA::Epetra::Interface::setParams before evalulating J.

Reimplemented from NOX::Epetra::Group.

References NOX::Epetra::Group::computeJacobian(), and NOX::Abstract::Group::Ok.

double LOCA::Epetra::Group::computeScaledDotProduct ( const NOX::Abstract::Vector a,
const NOX::Abstract::Vector b 
) const
virtual

Compute a scaled dot product.

The implementation here uses the scaling vector $ s$ if one is supplied:

\[ \sum_{i=1}^n a_i*b_i*s_i*s_i. \]

If the scaling vector is not provided, the standard dot product is used.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

References as(), NOX::Abstract::Vector::clone(), NOX::DeepCopy, NOX::Abstract::Vector::innerProduct(), NOX::Abstract::Vector::length(), and NOX::Abstract::Vector::scale().

void LOCA::Epetra::Group::postProcessContinuationStep ( LOCA::Abstract::Iterator::StepStatus  stepStatus)
virtual

Perform any postprocessing after a continuation step finishes.

The stepStatus argument indicates whether the step was successful. The implementation here is to call the corresponding method in the interface.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

void LOCA::Epetra::Group::preProcessContinuationStep ( LOCA::Abstract::Iterator::StepStatus  stepStatus)
virtual

Perform any preprocessing before a continuation step starts.

The stepStatus argument indicates whether the previous step was successful. The implementation here is to call the corresponding method in the interface.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

void LOCA::Epetra::Group::projectToDraw ( const NOX::Abstract::Vector x,
double *  px 
) const
virtual

Projects solution to a few scalars for multiparameter continuation.

This method is called every time a solution is saved by the multiparameter continuation code MF for later visualization and should project the solution vector down to a few scalars. The array px will be preallocated to the proper length given by projectToDrawDimension().

The implementation here is to call the corresponding method in the interface.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

int LOCA::Epetra::Group::projectToDrawDimension ( ) const
virtual

Returns the dimension of the project to draw array.

The implementation here is to call the corresponding method in the interface.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

void LOCA::Epetra::Group::scaleVector ( NOX::Abstract::Vector x) const
virtual

Scales a vector using scaling vector.

The implementation here uses the scaling vector $ s$ if one is supplied:

\[ x_i = a_i*s_i. \]

If the scaling vector is not provided, the vector is rescaled by the square root of its length.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

References NOX::Abstract::Vector::length(), and NOX::Abstract::Vector::scale().

Member Data Documentation

Teuchos::RCP<Epetra_Vector> LOCA::Epetra::Group::tmpVectorPtr2
protected

Extra vector needed for intermediate calculations of LOCA routines.

NOTE: there already is a tmpVectorPtr in the NOX::Epetra::Group. This is a second temporary vector if that one extra isn't enough.


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