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::Homotopy::Group Class Reference

LOCA's Homotopy Algorithm. More...

#include <LOCA_Homotopy_Group.H>

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

Public Member Functions

 Group (Teuchos::ParameterList &locaSublist, const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Homotopy::AbstractGroup > &g, double scaleRandom=1.0, double scaleInitialGuess=0.0)
 Constructor to set the base group and generate the "%Stepper" sublist for homotopy continuation. More...
 
 Group (Teuchos::ParameterList &locaSublist, const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Homotopy::AbstractGroup > &g, const NOX::Abstract::Vector &randomVector)
 Constructor with a user supplied random vector.
 
 Group (const Group &source, NOX::CopyType type=NOX::DeepCopy)
 Copy constructor.
 
virtual ~Group ()
 Destructor.
 
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 homotopy residual $g$.
 
virtual
NOX::Abstract::Group::ReturnType 
computeJacobian ()
 Compute the Jacobian derivative of the homotopy residual $g$.
 
virtual
NOX::Abstract::Group::ReturnType 
computeGradient ()
 Compute gradient of homotopy residual $g$.
 
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 homotopy Jacobian vector product.
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Computes the homotopy Jacobian-transpose vector product.
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianInverse (Teuchos::ParameterList &params, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const
 Applies the inverse of the homotopy Jacobian matrix.
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Applies Jacobian for homotopy system.
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianTransposeMultiVector (const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Applies Jacobian-transpose for homotopy system.
 
virtual
NOX::Abstract::Group::ReturnType 
applyJacobianInverseMultiVector (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector &input, NOX::Abstract::MultiVector &result) const
 Applies Jacobian inverse for homotopy system.
 
virtual bool isF () const
 Return true if the homotopy residual $ g$ is valid.
 
virtual bool isJacobian () const
 Return true if the homotopy Jacobian is valid.
 
virtual bool isGradient () const
 Return true if the homotopy gradient is valid.
 
virtual bool isNewton () const
 Return true if the homotopy Newton direction is valid.
 
virtual const
NOX::Abstract::Vector
getX () const
 Return homotopy solution vector $ x$.
 
virtual const
NOX::Abstract::Vector
getF () const
 Return homotopy residual $ g$.
 
virtual double getNormF () const
 Return 2-norm of $ g$.
 
virtual const
NOX::Abstract::Vector
getGradient () const
 Return homotopy gradient.
 
virtual const
NOX::Abstract::Vector
getNewton () const
 Return homotopy Newton direction.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getXPtr () const
 Return homotopy solution vector $ x$.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getFPtr () const
 Return homotopy residual $ g$.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getGradientPtr () const
 Return homotopy gradient.
 
virtual Teuchos::RCP< const
NOX::Abstract::Vector
getNewtonPtr () const
 Return homotopy Newton direction.
 
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.
 
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 solution and continuation parameter after successful continuation step.
 
virtual void printSolution (const NOX::Abstract::Vector &x_, const double conParam) const
 Function to print out solution and continuation parameter after successful continuation step.
 
- 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 NOX::Abstract::Group
 Group ()
 Constructor. More...
 
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::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...
 

Protected Member Functions

void resetIsValidFlags ()
 Reset the isValid flags to false. More...
 
void setStepperParameters (Teuchos::ParameterList &params)
 Creates and sets the "Stepper" parameter sublist.
 

Protected Attributes

Teuchos::RCP< LOCA::GlobalDataglobalData
 Pointer LOCA global data object.
 
Teuchos::RCP
< LOCA::Homotopy::AbstractGroup
grpPtr
 Stores the underlying loca group.
 
Teuchos::RCP
< NOX::Abstract::Vector
gVecPtr
 Stores the homotopy residual vector, $ g $.
 
Teuchos::RCP
< NOX::Abstract::Vector
randomVecPtr
 Stores the random Vector, $ a $.
 
Teuchos::RCP
< NOX::Abstract::Vector
newtonVecPtr
 Stores the homotopy Newton vector, $ \frac{\partial r}{\partial x} $.
 
Teuchos::RCP
< NOX::Abstract::Vector
gradVecPtr
 
bool isValidF
 Is residual vector valid.
 
bool isValidJacobian
 Is Jacobian matrix valid.
 
bool isValidNewton
 Is Newton vector valid.
 
bool isValidGradient
 Is gradient vector valid.
 
LOCA::ParameterVector paramVec
 Copy of the ParameterVector for the underlying grpPtr. More...
 
double conParam
 Value of the homotopy continuation parameter. More...
 
int conParamID
 Continuatioin parameter ID number from the ParameterVector.
 
const std::string conParamLabel
 Contains the std::string used to identify the homotopy parameter in the ParameterVector object.
 
bool augmentJacForHomotopyNotImplemented
 Tracks whether the LOCA::Homotopy::Group method augmentJacobianForHomotopy is implemented. If not, the augmentation is applied during the applyJacobian assuming a matrix-free implementation.
 

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

LOCA's Homotopy Algorithm.

The HomotopyGroup is a concrete implementation of the LOCA::Continuation::AbstractGroup that modifies the set of nonlinear equations to be solved to allow for Homotopy to be applied to the system. This object should be used in conjunction with the LOCA::Stepper object to drive the continuation. This algorithm solves a system of nonlinear equations supplied by the user ( $ F(x) $) through continuation. An artificial parameter $ \lambda $ is used to control the continuation. The idea is to solve a simple equation starting at $ \lambda $ = 0 and, using the solution from the previous step, solve systems of equations that gets progressively closer to the true system of interest ( at $ \lambda $ = 1.0 we recover the original equations $ F(x) $). By constraining the definition of $ g(x, \lambda) $ and using artificial parameter contiuation, the continuation branch should be free of multiplicity and bifurcation phenomena.

The modified system of equations, $ g(x, \lambda) $, supplied by the HomotopyGroup is defined as:

\[ g(x, \lambda) = \lambda F(x) + (1.0 - \lambda)(x - a) \]

where $ x$ is the solution vector, $ \lambda $ is an artificial parameter, $ F(x) $ is the set of nonlinear equations the user supplies, $ g(x) $ is the corresponding set of homotopy equations that LOCA will solve, and $ a $ is a random vector.

This group requires the loca Stepper for continuation from $ \lambda $ = 0.0 (a simple set of equations to solve) to $ \lambda $ = 1.0 (the set of equations requested by the user, $ F(x) $). The Homotopy::Group will generate the Stepper parameter sublist in the parameter list that is passed in to the constructor. The user is free to modify this list (it sets default values) before passing it into the stepper object but should NOT change the starting and stopping values for the continuation parameter.

References:

Constructor & Destructor Documentation

LOCA::Homotopy::Group::Group ( Teuchos::ParameterList locaSublist,
const Teuchos::RCP< LOCA::GlobalData > &  global_data,
const Teuchos::RCP< LOCA::Homotopy::AbstractGroup > &  g,
double  scaleRandom = 1.0,
double  scaleInitialGuess = 0.0 
)

Constructor to set the base group and generate the "%Stepper" sublist for homotopy continuation.

The locaSublist variable is the "LOCA" sublist (of type Teuchos::ParameterList) that will be used in loca continuation runs.

The variables scalarRandomVector and scalarInitialGuess are used to give some control over the generation of the random vector. In certain instances we have seen the random vector force the solution to a set of variables that are unphysical and could break the function evaluations (cause them to return nan). For example, in heat transfer problems, the temperature could be the dependent variable. If the solution vector has an unphysical temperature ( the random vector could force the temperature to negative or near zero values for the solution at $ \lambda = 0$) then property evaluations could break. The random vector can be modified to keep the values near the initial guess based on values supplied to the constructor of the HomotopyGroup:

\[ a = abs(r) * \mbox{scalarRandom} + x_o * \mbox{scalarInitialGuess} \]

where $ r $ is the random vector generated by a call to NOX::Abstract::Vector::random(), $ \mbox{scalarRandom} $ is a scalar value, $ x_o $ is the initial guess to the solution vector, and $ \mbox{scalarInitialGuess} $ is a scalar value. The defualt values force the random vector to be calculated as:

\[ a = abs(r) \]

IMPORTANT: For homotopy to work correctly you should not change the starting and stopping parameter values (0.0 and 1.0 respectively) set in the "%Stepper" sublist.

References NOX::Abstract::Vector::abs(), LOCA::ParameterVector::addParameter(), conParam, conParamID, conParamLabel, LOCA::ParameterVector::getIndex(), grpPtr, paramVec, NOX::Abstract::Vector::random(), randomVecPtr, resetIsValidFlags(), setStepperParameters(), and NOX::Abstract::Vector::update().

Member Function Documentation

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

Compute $\partial F/\partial p$ for each parameter $ p$ 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 NOX::Abstract::MultiVector::scale(), NOX::Abstract::MultiVector::subView(), and NOX::Abstract::MultiVector::update().

void LOCA::Homotopy::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.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

void LOCA::Homotopy::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.

Reimplemented from LOCA::MultiContinuation::AbstractGroup.

void LOCA::Homotopy::Group::resetIsValidFlags ( )
protected

Reset the isValid flags to false.

  This is called when the solution vector or parameter vector

is changed.

Referenced by Group().

Member Data Documentation

double LOCA::Homotopy::Group::conParam
protected

Value of the homotopy continuation parameter.

Ranges from 0.0 (easy solution) to 1.0 (solution to the system of interest).

Referenced by copy(), and Group().

Teuchos::RCP<NOX::Abstract::Vector> LOCA::Homotopy::Group::gradVecPtr
protected

Stores the homotopy gradient vector if needed, $ \frac{\partial r}{\partial x} $.

Referenced by copy(), and Group().

LOCA::ParameterVector LOCA::Homotopy::Group::paramVec
protected

Copy of the ParameterVector for the underlying grpPtr.

We copy this and then add the homotopy parameter to the list.

Referenced by copy(), and Group().


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