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

A group representing the Moore-Spence turning point equations. More...

#include <LOCA_TurningPoint_MooreSpence_ExtendedGroup.H>

Inheritance diagram for LOCA::TurningPoint::MooreSpence::ExtendedGroup:
Inheritance graph
[legend]
Collaboration diagram for LOCA::TurningPoint::MooreSpence::ExtendedGroup:
Collaboration graph
[legend]

Public Member Functions

 ExtendedGroup (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &tpParams, const Teuchos::RCP< LOCA::TurningPoint::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.
 
double getBifParam () const
 Get bifurcation parameter.
 
double lTransNorm (const NOX::Abstract::Vector &n) const
 Defines null vector normalization $l^Tn$ equation.
 
void lTransNorm (const NOX::Abstract::MultiVector &n, NOX::Abstract::MultiVector::DenseMatrix &result) const
 null vector normalization for multivectors More...
 
Teuchos::RCP
< NOX::Abstract::Vector
getLengthVector () const
 Return length normalization vector.
 
Implementation of NOX::Abstract::Group virtual methods
virtual NOX::Abstract::Groupoperator= (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 turning point equation residual $G$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeJacobian ()
 Compute the blocks of the Jacobian derivative of $G$. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeGradient ()
 Gradient computation is not defined for this group.
 
virtual
NOX::Abstract::Group::ReturnType 
computeNewton (Teuchos::ParameterList &params)
 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. More...
 
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 &params, 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 &params, 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 $G$ 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 $z$.
 
virtual const
NOX::Abstract::Vector
getF () const
 Return extended equation residual $G(z)$.
 
virtual double getNormF () const
 Return 2-norm of $G(z)$.
 
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 extended solution vector $z$.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getFPtr () const
 Return extended equation residual $G(z)$.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getGradientPtr () const
 Vector returned is not valid.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getNewtonPtr () const
 Return 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.
 
Implementation of LOCA::MultiContinuation::AbstractGroup

virtual methods

virtual void copy (const NOX::Abstract::Group &source)
 Assignment operator.
 
virtual void setParamsMulti (const std::vector< int > &paramIDs, 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 ParameterVectorgetParams () 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 > &paramIDs, 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...
 
Implementation of LOCA::Abstract::TransposeSolveGroup

virtual methods

virtual
NOX::Abstract::Group::ReturnType 
applyJacobianTransposeInverse (Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Solve Jacobian-tranpose system.
 
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.
 
- Public Member Functions inherited from LOCA::Extended::MultiAbstractGroup
 MultiAbstractGroup ()
 Default constructor.
 
virtual ~MultiAbstractGroup ()
 Destructor.
 
virtual Teuchos::RCP
< NOX::Abstract::Group
getNestedGroup ()
 Override from NOX::Abstract::Group base class. Calls getUnderlyingGroup() from this class.
 
virtual Teuchos::RCP< const
NOX::Abstract::Group
getNestedGroup () const
 Override from NOX::Abstract::Group base class. Calls getUnderlyingGroup() from this class.
 
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 &params, 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 &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
 
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...
 
- Public Member Functions inherited from LOCA::Abstract::TransposeSolveGroup
 TransposeSolveGroup ()
 Constructor.
 
virtual ~TransposeSolveGroup ()
 Destructor.
 

Protected Types

enum  NullVectorScaling { NVS_None, NVS_OrderOne, NVS_OrderN }
 Enumerated type determining type of scaling. More...
 

Protected Member Functions

void setBifParam (double param)
 Set bifurcation parameter.
 
void setupViews ()
 Sets up multivector views.
 
void init (bool perturbSoln=false, double perturbSize=0.0)
 Initializes group.
 
void scaleNullVector (NOX::Abstract::Vector &a)
 Scale null vector.
 

Protected Attributes

Teuchos::RCP< LOCA::GlobalDataglobalData
 Pointer LOCA global data object.
 
Teuchos::RCP
< LOCA::Parameter::SublistParser
parsedParams
 Parsed top-level parameters.
 
Teuchos::RCP
< Teuchos::ParameterList
turningPointParams
 Bifurcation parameter list.
 
Teuchos::RCP
< LOCA::TurningPoint::MooreSpence::AbstractGroup
grpPtr
 Pointer to base group that defines $F$.
 
LOCA::TurningPoint::MooreSpence::ExtendedMultiVector xMultiVec
 Stores the extended solution vector.
 
LOCA::TurningPoint::MooreSpence::ExtendedMultiVector fMultiVec
 Stores the extended residual vector and df/dp.
 
LOCA::TurningPoint::MooreSpence::ExtendedMultiVector newtonMultiVec
 Stores the extended Newton vector.
 
Teuchos::RCP
< NOX::Abstract::MultiVector
lengthMultiVec
 Stores the length normalization vector.
 
Teuchos::RCP
< LOCA::TurningPoint::MooreSpence::ExtendedVector
xVec
 Stores view of first column of xMultiVec.
 
Teuchos::RCP
< LOCA::TurningPoint::MooreSpence::ExtendedVector
fVec
 Stores view of first column of fMultiVec.
 
Teuchos::RCP
< LOCA::TurningPoint::MooreSpence::ExtendedMultiVector
ffMultiVec
 Stores view of first column of fMultiVec as a multivec.
 
Teuchos::RCP
< LOCA::TurningPoint::MooreSpence::ExtendedMultiVector
dfdpMultiVec
 Stores view of df/dp columns of fMultiVec.
 
Teuchos::RCP
< LOCA::TurningPoint::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
< LOCA::TurningPoint::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.
 
bool updateVectorsEveryContinuationStep
 Flag indicating whether to update $n$ every continuation step.
 
NullVectorScaling nullVecScaling
 Null vector scaling method.
 
bool multiplyMass
 Multiply null vectors by mass matrix.
 
Teuchos::RCP
< LOCA::TimeDependent::AbstractGroup
tdGrp
 Time dependent group interface for multiplying by mass matrix.
 
Teuchos::RCP
< NOX::Abstract::Vector
tmp_mass
 Temporary vector for mulplying null vectors by mass matrix.
 

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...
 

Detailed Description

A group representing the Moore-Spence turning point equations.

The LOCA::TurningPoint::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 turning point:

\[ G(z) = \left[ \begin{array}{c} F(x,p) \\ Jn \\ l^Tn-1 \end{array} \right] = 0 \]

where $z = [x, n, p]\in\Re^{2n+1}$, $x$ is the solution vector, $n$ is the null vector, $l$ is the length normalization vector and $J$ is the Jacobian of F.

The group stores an underlying group of type LOCA::TurningPoint::MooreSpence AbstractGroup to represent the equations $F(x,p) = 0$ and to manipulate the underlying Jacobian $J$. Note that the entire extended Jacobian $D_z G$ is not stored in memory, rather a block-elimination algorithm (bordering algorithm) is used to compute linear solves of the extended Jacobian (see LOCA::TurningPoint::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 turning point or to the LOCA::Stepper to compute a family of turning 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 turning point equations.

The class is intialized via the tpParams parameter list argument to the constructor. The parameters this class recognizes are:

Member Enumeration Documentation

Enumerated type determining type of scaling.

Enumerator
NVS_OrderOne 

No scaling.

NVS_OrderN 

Scale to O(1)

Scale to O(N) when N is the vector length

Member Function Documentation

NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobian ( const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
virtual

Computes the extended Jacobian vector product.

This method computes the extended Jacobian vector product

\[ \begin{bmatrix} J & 0 & \frac{\partial F}{\partial p} \\ \frac{\partial Jn}{\partial x} & J & \frac{\partial Jn}{\partial p} \\ 0 & l^T & 0 \end{bmatrix} \begin{bmatrix} a \\ b \\ c \end{bmatrix} = \begin{bmatrix} Ja + \frac{\partial F}{\partial p}c \\ \frac{\partial Jn}{\partial x}a + Jb + \frac{\partial Jn}{\partial p}c \\ l^T b \end{bmatrix} \]

using the applyJacobian and computeDJnDxa methods of the underlying group where $a$, $b$, and $c$ are the solution, null-vector, and paramter components of the given vector input. Vectors input and result must be of type LOCA::TurningPoint::MooreSpence::ExtendedVector, otherwise an error is thrown.

Reimplemented from NOX::Abstract::Group.

References NOX::Abstract::Vector::createMultiVector(), and NOX::DeepCopy.

NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianInverse ( Teuchos::ParameterList params,
const NOX::Abstract::Vector input,
NOX::Abstract::Vector result 
) const
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.

NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianInverseMultiVector ( Teuchos::ParameterList params,
const NOX::Abstract::MultiVector input,
NOX::Abstract::MultiVector result 
) const
virtual

Applies Jacobian inverse for extended system.

Uses a LOCA::TurningPoint::MooreSpence::SolverStrategy instantiated by the LOCA::TurningPoint::MooreSpence::SolverFactory to implement the solve.

Reimplemented from NOX::Abstract::Group.

NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::computeDfDpMulti ( const std::vector< int > &  paramIDs,
NOX::Abstract::MultiVector dfdp,
bool  isValidF 
)
virtual
NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::computeF ( )
virtual

Compute the turning point equation residual $G$.

This method fills the extended residual

\[ G(z) = \left[ \begin{array}{c} F(x,p) \\ Jn \\ l^Tn-1 \end{array} \right]. \]

The solution component residual $F(x,p)$ and the null-vector residual $Jn$ are calculated via the computeF and applyJacobian methods of the underlying group.

Implements NOX::Abstract::Group.

References NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType LOCA::TurningPoint::MooreSpence::ExtendedGroup::computeJacobian ( )
virtual

Compute the blocks of the Jacobian derivative of $G$.

This method computes the $J$, $\partial F/\partial p$, and $\partial Jn/\partial p$ blocks of the extended Jacobian:

\[ D_z G(z) = \begin{bmatrix} J & 0 & \frac{\partial F}{\partial p} \\ \frac{\partial Jn}{\partial x} & J & \frac{\partial Jn}{\partial p} \\ 0 & l^T & 0 \end{bmatrix} \]

by calling the computeJacobian, computeDfDp, and computeDJnDp methods of the underlying group. The second derivative matrix $\partial Jn/\partial x$ is not calculated since only its action on vectors is needed for linear solves using the bordering algorithm.

Reimplemented from NOX::Abstract::Group.

References NOX::Abstract::Group::Ok, and Teuchos::rcp().

void LOCA::TurningPoint::MooreSpence::ExtendedGroup::lTransNorm ( const NOX::Abstract::MultiVector n,
NOX::Abstract::MultiVector::DenseMatrix result 
) const

null vector normalization for multivectors

Note: result should have 1 row and n.numVectors() columns.

References NOX::Abstract::MultiVector::multiply().

void LOCA::TurningPoint::MooreSpence::ExtendedGroup::postProcessContinuationStep ( LOCA::Abstract::Iterator::StepStatus  stepStatus)
virtual

Perform any postprocessing after a continuation step finishes.

The stepStatus argument indicates whether the step was successful.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

References NOX::Utils::StepperDetails, and LOCA::Abstract::Iterator::Successful.

void LOCA::TurningPoint::MooreSpence::ExtendedGroup::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.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

References NOX::Utils::StepperDetails.

void LOCA::TurningPoint::MooreSpence::ExtendedGroup::printSolution ( const double  conParam) const
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.

void LOCA::TurningPoint::MooreSpence::ExtendedGroup::printSolution ( const NOX::Abstract::Vector x_,
const double  conParam 
) const
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::TurningPoint::MooreSpence::ExtendedVector::getBifParam(), LOCA::TurningPoint::MooreSpence::ExtendedVector::getNullVec(), LOCA::TurningPoint::MooreSpence::ExtendedVector::getXVec(), and NOX::Utils::StepperDetails.


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