NOX
Development
|
A group representing the minimally augemented turning point equations. More...
#include <LOCA_TurningPoint_MinimallyAugmented_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 > &tpParams, const Teuchos::RCP< LOCA::TurningPoint::MinimallyAugmented::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. | |
Teuchos::RCP< const NOX::Abstract::Vector > | getLeftNullVec () const |
Returns left null vector. | |
Teuchos::RCP< const NOX::Abstract::Vector > | getRightNullVec () const |
Returns right null vector v. | |
Teuchos::RCP< const NOX::Abstract::Vector > | getAVec () const |
Returns "A" vector. | |
Teuchos::RCP< const NOX::Abstract::Vector > | getBVec () const |
Returns "B". | |
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 turning point equation residual . More... | |
virtual NOX::Abstract::Group::ReturnType | computeJacobian () |
Compute the blocks of the Jacobian derivative of . More... | |
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 |
Computes the extended Jacobian transpose vector product. | |
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. | |
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 for extended 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. | |
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. | |
Implementation of LOCA::MultiContinuation::AbstractGroup | |
virtual methods | |
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 (std::string paramID, double val) |
Set parameter indexed by paramID. | |
virtual void | setParam (int 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... | |
Implementation of | |
LOCA::BorderedSystem::AbstractGroup virtual methods | |
virtual int | getBorderedWidth () const |
Return the total width of the bordered rows/columns. | |
virtual Teuchos::RCP< const NOX::Abstract::Group > | getUnborderedGroup () const |
Get bottom-level unbordered group. | |
virtual bool | isCombinedAZero () const |
Indicates whether combined A block is zero. | |
virtual bool | isCombinedBZero () const |
Indicates whether combined B block is zero. | |
virtual bool | isCombinedCZero () const |
Indicates whether combined C block is zero. | |
virtual void | extractSolutionComponent (const NOX::Abstract::MultiVector &v, NOX::Abstract::MultiVector &v_x) const |
virtual void | extractParameterComponent (bool use_transpose, const NOX::Abstract::MultiVector &v, NOX::Abstract::MultiVector::DenseMatrix &v_p) const |
virtual void | loadNestedComponents (const NOX::Abstract::MultiVector &v_x, const NOX::Abstract::MultiVector::DenseMatrix &v_p, NOX::Abstract::MultiVector &v) const |
virtual void | fillA (NOX::Abstract::MultiVector &A) const |
Fill the combined A block as described above. | |
virtual void | fillB (NOX::Abstract::MultiVector &B) const |
Fill the combined B block as described above. | |
virtual void | fillC (NOX::Abstract::MultiVector::DenseMatrix &C) const |
Fill the combined C block as described above. | |
Implementation of LOCA::Abstract::TransposeSolveGroup | |
virtual methods | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTransposeInverse (Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
Solve Jacobian-tranpose system. | |
virtual NOX::Abstract::Group::ReturnType | applyJacobianTransposeInverseMultiVector (Teuchos::ParameterList ¶ms, 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< 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... | |
Public Member Functions inherited from LOCA::BorderedSystem::AbstractGroup | |
AbstractGroup () | |
Constructor. | |
virtual | ~AbstractGroup () |
Destructor. | |
Public Member Functions inherited from LOCA::Abstract::TransposeSolveGroup | |
TransposeSolveGroup () | |
Constructor. | |
virtual | ~TransposeSolveGroup () |
Destructor. | |
Protected Member Functions | |
void | setBifParam (double param) |
Set bifurcation parameter. | |
Protected Attributes | |
Teuchos::RCP< LOCA::GlobalData > | globalData |
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::MinimallyAugmented::AbstractGroup > | grpPtr |
Pointer to base group that defines . | |
Teuchos::RCP < LOCA::TurningPoint::MinimallyAugmented::Constraint > | constraint |
Pointer to the constraint equation. | |
Teuchos::RCP < LOCA::MultiContinuation::ConstrainedGroup > | conGroup |
Pointer to constrained group implementation. | |
int | bifParamID |
Stores the bifurcation parameter index. | |
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 minimally augemented turning point equations.
The LOCA::TurningPoint::MinimallyAugmented::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:
where , is the solution vector, is the bifurcation parameter, is the Jacobian of , and is a measure of the singularity of and is defined via
for any vectors and in . Using these relationships, it is easy to show
The group stores an underlying group of type LOCA::TurningPoint::MinimallyAugmented::AbstractGroup to represent the equations and to manipulate the underlying Jacobian . This interface defines methods for computing the derivatives and as well.
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.
The class is intialized via the tpParams
parameter list argument to the constructor. The parameters this class recognizes are:
|
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.
|
virtual |
Compute the turning point equation residual .
This method fills the extended residual
Implements NOX::Abstract::Group.
|
virtual |
Compute the blocks of the Jacobian derivative of .
This method computes the , , and blocks of the extended Jacobian:
Reimplemented from NOX::Abstract::Group.
|
virtual |
Given the vector v
, extract the parameter components of all of the nested subvectors in v
down to the solution component for the unbordered group.
Implements LOCA::BorderedSystem::AbstractGroup.
|
virtual |
Given the vector v
, extract the underlying solution component corresponding to the unbordered group.
Implements LOCA::BorderedSystem::AbstractGroup.
|
virtual |
Given the solution component v_x
and combined parameter components v_p
, distribute these components through the nested sub-vectors in v
.
Implements LOCA::BorderedSystem::AbstractGroup.
|
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::Extended::Vector::getScalar(), LOCA::MultiContinuation::ExtendedVector::getXVec(), and NOX::Utils::StepperDetails.