NOX
Development
|
A group representing the Moore-Spence Hopf equations. More...
#include <LOCA_Hopf_MooreSpence_ExtendedGroup.H>
Public Member Functions | |
ExtendedGroup (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &hopfParams, const Teuchos::RCP< LOCA::Hopf::MooreSpence::AbstractGroup > &g) | |
Constructor with initial data passed through parameter lists. | |
ExtendedGroup (const ExtendedGroup &source, NOX::CopyType type=NOX::DeepCopy) | |
Copy constructor. | |
virtual | ~ExtendedGroup () |
Destructor. | |
Implementation of NOX::Abstract::Group virtual methods | |
virtual NOX::Abstract::Group & | operator= (const NOX::Abstract::Group &source) |
Assignment operator. | |
virtual Teuchos::RCP < NOX::Abstract::Group > | clone (NOX::CopyType type=NOX::DeepCopy) const |
Cloning function. | |
virtual void | setX (const NOX::Abstract::Vector &y) |
Set the solution vector, x, to y. | |
virtual void | computeX (const NOX::Abstract::Group &g, const NOX::Abstract::Vector &d, double step) |
Compute this.x = grp.x + step * d. | |
virtual NOX::Abstract::Group::ReturnType | computeF () |
Compute the Hopf point equation residual . | |
virtual NOX::Abstract::Group::ReturnType | computeJacobian () |
Compute the blocks of the Jacobian derivative of . | |
virtual NOX::Abstract::Group::ReturnType | computeGradient () |
Gradient computation is not defined for this group. | |
virtual NOX::Abstract::Group::ReturnType | computeNewton (Teuchos::ParameterList ¶ms) |
Compute Newton direction using applyJacobianInverse(). | |
virtual NOX::Abstract::Group::ReturnType | applyJacobian (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Computes the extended Jacobian vector product. | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Jacobian transpose product is not defined by this group. | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianInverse (Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Applies the inverse of the extended Jacobian matrix using the bordering algorithm. More... | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
Applies Jacobian for extended system. | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTransposeMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
Jacobian transpose not defined for this system. | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianInverseMultiVector (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
Applies Jacobian inverse for extended system. More... | |
virtual bool | isF () const |
Return true if the extended residual is valid. | |
virtual bool | isJacobian () const |
Return true if the extended Jacobian is valid. | |
virtual bool | isGradient () const |
Always returns false. | |
virtual bool | isNewton () const |
Return true if the extended Newton direction is valid. | |
virtual const NOX::Abstract::Vector & | getX () const |
Return extended solution vector . | |
virtual const NOX::Abstract::Vector & | getF () const |
Return extended equation residual . | |
virtual double | getNormF () const |
Return 2-norm of . | |
virtual const NOX::Abstract::Vector & | getGradient () const |
Vector returned is not valid. | |
virtual const NOX::Abstract::Vector & | getNewton () const |
Return extended Newton direction. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getXPtr () const |
Return RCP to extended solution vector . | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getFPtr () const |
Return RCP to extended equation residual . | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getGradientPtr () const |
Vector returned is not valid. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getNewtonPtr () const |
Return RCP to extended Newton direction. | |
virtual double | getNormNewtonSolveResidual () const |
Return the norm of the Newton solve residual. | |
Implementation of LOCA::Extended::MultiAbstractGroup | |
virtual methods | |
virtual Teuchos::RCP< const LOCA::MultiContinuation::AbstractGroup > | getUnderlyingGroup () const |
Return underlying group. | |
virtual Teuchos::RCP < LOCA::MultiContinuation::AbstractGroup > | getUnderlyingGroup () |
Return underlying group. | |
Public Member Functions inherited from LOCA::Extended::MultiAbstractGroup | |
MultiAbstractGroup () | |
Default constructor. | |
virtual | ~MultiAbstractGroup () |
Destructor. | |
virtual Teuchos::RCP< const LOCA::MultiContinuation::AbstractGroup > | getBaseLevelUnderlyingGroup () const |
Return base-level underlying group. More... | |
virtual Teuchos::RCP < LOCA::MultiContinuation::AbstractGroup > | getBaseLevelUnderlyingGroup () |
Return base-level underlying group. More... | |
Public Member Functions inherited from NOX::Abstract::Group | |
Group () | |
Constructor. More... | |
virtual | ~Group () |
Destructor. | |
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 | applyRightPreconditioningMultiVector (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const |
applyRightPreconditioning for multiple right-hand sides More... | |
virtual const NOX::Abstract::Vector & | getScaledX () const |
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... | |
Public Member Functions inherited from LOCA::MultiContinuation::AbstractGroup | |
AbstractGroup () | |
Default constructor. | |
virtual | ~AbstractGroup () |
Destructor. | |
virtual double | computeScaledDotProduct (const NOX::Abstract::Vector &a, const NOX::Abstract::Vector &b) const |
Compute a scaled dot product. More... | |
virtual void | scaleVector (NOX::Abstract::Vector &x) const |
Scales a vector using scaling vector. More... | |
Implementation of LOCA::MultiContinuation::AbstractGroup | |
Teuchos::RCP< LOCA::GlobalData > | globalData |
Pointer LOCA global data object. | |
Teuchos::RCP < LOCA::Parameter::SublistParser > | parsedParams |
Parsed top-level parameters. | |
Teuchos::RCP < Teuchos::ParameterList > | hopfParams |
Bifurcation parameter list. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::AbstractGroup > | grpPtr |
Pointer to base group that defines . | |
LOCA::Hopf::MooreSpence::ExtendedMultiVector | xMultiVec |
Stores the extended solution vector. | |
LOCA::Hopf::MooreSpence::ExtendedMultiVector | fMultiVec |
Stores the extended residual vector and df/dp. | |
LOCA::Hopf::MooreSpence::ExtendedMultiVector | newtonMultiVec |
Stores the extended Newton vector. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | lengthMultiVec |
Stores the length normalization vector. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::ExtendedVector > | xVec |
Stores view of first column of xMultiVec. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::ExtendedVector > | fVec |
Stores view of first column of fMultiVec. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::ExtendedMultiVector > | ffMultiVec |
Stores view of first column of fMultiVec as a multivec. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::ExtendedMultiVector > | dfdpMultiVec |
Stores view of df/dp columns of fMultiVec. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::ExtendedVector > | newtonVec |
Stores view of first column of newtonMultiVec. | |
Teuchos::RCP < NOX::Abstract::Vector > | lengthVec |
Stores view of first column of lengthMultiVec. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | massTimesY |
Stores mass matrix times real component of eigenvector. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | minusMassTimesZ |
Stores negative of mass matrix times imaginary component of eigenvector. | |
Teuchos::RCP < LOCA::Hopf::MooreSpence::SolverStrategy > | solverStrategy |
Stores bordering strategy. | |
std::vector< int > | index_f |
Stores indices for getting f part of fMultiVec. | |
std::vector< int > | index_dfdp |
Stores indices for getting df/dp part of fMultiVec. | |
std::vector< int > | bifParamID |
Stores the bifurcation parameter index. | |
bool | isValidF |
Is residual vector valid. | |
bool | isValidJacobian |
Is Jacobian matrix valid. | |
bool | isValidNewton |
Is Newton vector valid. | |
virtual void | copy (const NOX::Abstract::Group &source) |
Assignment operator. | |
virtual void | setParamsMulti (const std::vector< int > ¶mIDs, const NOX::Abstract::MultiVector::DenseMatrix &vals) |
Set parameters indexed by (integer) paramIDs. | |
virtual void | setParams (const ParameterVector &p) |
Set the parameter vector in the group to p. | |
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. | |
virtual const ParameterVector & | getParams () const |
Return a const reference to the paramter vector 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 NOX::Abstract::Group::ReturnType | computeDfDpMulti (const std::vector< int > ¶mIDs, NOX::Abstract::MultiVector &dfdp, bool isValidF) |
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. | |
virtual int | projectToDrawDimension () const |
Returns the dimension of the project to draw array. | |
virtual void | printSolution (const double conParam) const |
Function to print out extended solution and continuation parameter after successful continuation step. More... | |
virtual void | printSolution (const NOX::Abstract::Vector &x_, const double conParam) const |
Function to print out extended solution and continuation parameter after successful continuation step. More... | |
double | getBifParam () const |
Get bifurcation parameter. | |
double | getFrequency () const |
Get bifurcation frequency. | |
double | lTransNorm (const NOX::Abstract::Vector &z) const |
Defines null vector normalization equation. | |
void | lTransNorm (const NOX::Abstract::MultiVector &z, NOX::Abstract::MultiVector::DenseMatrix &result) const |
null vector normalization for multivectors More... | |
void | setBifParam (double param) |
Set bifurcation parameter. | |
void | setupViews () |
Sets up multivector views. | |
void | init (bool perturbSoln=false, double perturbSize=0.0) |
Initializes group. | |
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... | |
A group representing the Moore-Spence Hopf equations.
The LOCA::Hopf::MooreSpence::ExtendedGroup is a concrete implementation of the NOX::Abstract::Group, LOCA::MultiContinuation::AbstractGroup and LOCA::Extended::MultiAbstractGroup that defines the following extended set of equations that are regular at a generic Hopf point:
where , is the solution vector, is the complex eigenvector of with corresponding eigenvalues , is the length normalization vector and is the Jacobian of F w.r.t .
The group stores an underlying group of type LOCA::Hopf::MooreSpence AbstractGroup to represent the equations and to manipulate the underlying Jacobian . Note that the entire extended Jacobian is not stored in memory, rather a block-elimination algorithm (bordering algorithm) is used to compute linear solves of the extended Jacobian (see LOCA::Hopf::MooreSpence::SolverFactory() for more details).
This class implements all of the NOX::Abstract::Group, LOCA::MultiContinuation::AbstractGroup, and LOCA::Extended::MultiAbstractGroup methods for this extended set of equations and therefore is a complete group which can be passed to most NOX solvers to locate a single Hopf point or to the LOCA::Stepper to compute a family of Hopf points in a second parameter.
However, Jacobian-tranpose operations and gradient calculations cannot be implemented efficiently and therefore gradient-base nonlinear solvers such as steepest descent and Trust region methods cannot be used to solve the extended Hopf point equations.
The class is intialized via the hopfParams
parameter list argument to the constructor. The parameters this class recognizes are:
|
virtual |
Applies the inverse of the extended Jacobian matrix using the bordering algorithm.
This method is a special case of applyJacobianInverseMultiVector() for a single right-hand-side.
Reimplemented from NOX::Abstract::Group.
References NOX::Abstract::Vector::createMultiVector(), and NOX::DeepCopy.
|
virtual |
Applies Jacobian inverse for extended system.
Uses a LOCA::Hopf::MooreSpence::SolverStrategy instantiated by the LOCA::Hopf::MooreSpence::SolverFactory to implement the solve.
Reimplemented from NOX::Abstract::Group.
|
virtual |
Compute for each parameter indexed by paramIDs. The first column of dfdp holds F, which is valid if isValidF is true. Otherwise F must be computed.
Implements LOCA::MultiContinuation::AbstractGroup.
References LOCA::Hopf::MooreSpence::ExtendedMultiVector::getImagEigenMultiVec(), LOCA::Hopf::MooreSpence::ExtendedMultiVector::getRealEigenMultiVec(), LOCA::Extended::MultiVector::getScalar(), LOCA::Hopf::MooreSpence::ExtendedMultiVector::getXMultiVec(), NOX::Abstract::MultiVector::numVectors(), and NOX::Abstract::Group::Ok.
void LOCA::Hopf::MooreSpence::ExtendedGroup::lTransNorm | ( | const NOX::Abstract::MultiVector & | z, |
NOX::Abstract::MultiVector::DenseMatrix & | result | ||
) | const |
null vector normalization for multivectors
Note: result should have 1 row and z.numVectors() columns.
References NOX::Abstract::MultiVector::multiply().
|
virtual |
Perform any postprocessing after a continuation step finishes.
The stepStatus
argument indicates whether the step was successful.
Reimplemented from LOCA::MultiContinuation::AbstractGroup.
|
virtual |
Perform any preprocessing before a continuation step starts.
The stepStatus
argument indicates whether the previous step was successful.
Reimplemented from LOCA::MultiContinuation::AbstractGroup.
|
virtual |
Function to print out extended solution and continuation parameter after successful continuation step.
This method prints the solution, null-vector, and parameter components of the extended solution vector using the printSolution method of the underlying group.
Reimplemented from LOCA::MultiContinuation::AbstractGroup.
References NOX::Utils::StepperDetails.
|
virtual |
Function to print out extended solution and continuation parameter after successful continuation step.
This method prints the solution, null-vector, and parameter components of the extended solution vector using the printSolution method of the underlying group.
Reimplemented from LOCA::MultiContinuation::AbstractGroup.
References LOCA::Hopf::MooreSpence::ExtendedVector::getBifParam(), LOCA::Hopf::MooreSpence::ExtendedVector::getFrequency(), LOCA::Hopf::MooreSpence::ExtendedVector::getImagEigenVec(), LOCA::Hopf::MooreSpence::ExtendedVector::getRealEigenVec(), LOCA::Hopf::MooreSpence::ExtendedVector::getXVec(), and NOX::Utils::StepperDetails.