|
NOX
Development
|
Specialization of LOCA::MultiContinuation::ExtendedGroup to natural continuation. More...
#include <LOCA_MultiContinuation_NaturalGroup.H>
Public Member Functions | |
| NaturalGroup (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &continuationParams, const Teuchos::RCP< LOCA::MultiContinuation::AbstractGroup > &grp, const Teuchos::RCP< LOCA::MultiPredictor::AbstractStrategy > &pred, const std::vector< int > ¶mIDs) | |
| Constructor. More... | |
| NaturalGroup (const NaturalGroup &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor. | |
| virtual | ~NaturalGroup () |
| 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 |
| Clone function. | |
Implementation of LOCA::MultiContinuation::AbstractStrategy virtual methods | |
| virtual void | copy (const NOX::Abstract::Group &source) |
| Copy. | |
Public Member Functions inherited from LOCA::MultiContinuation::ExtendedGroup | |
| ExtendedGroup (const ExtendedGroup &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor. | |
| virtual | ~ExtendedGroup () |
| Destructor. | |
| virtual void | setX (const NOX::Abstract::Vector &y) |
| Set the solution vector to y. | |
| virtual void | computeX (const NOX::Abstract::Group &g, const NOX::Abstract::Vector &d, double step) |
| Compute and return solution vector, x, where this.x = grp.x + step * d. | |
|
virtual NOX::Abstract::Group::ReturnType | computeF () |
| Compute extended continuation equations. | |
|
virtual NOX::Abstract::Group::ReturnType | computeJacobian () |
| Compute extended continuation jacobian. | |
|
virtual NOX::Abstract::Group::ReturnType | computeGradient () |
| Gradient is not defined for this system. | |
|
virtual NOX::Abstract::Group::ReturnType | computeNewton (Teuchos::ParameterList ¶ms) |
| Compute Newton direction for extended continuation system. | |
|
virtual NOX::Abstract::Group::ReturnType | applyJacobian (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
| Applies Jacobian for extended system. | |
|
virtual NOX::Abstract::Group::ReturnType | applyJacobianTranspose (const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
| Jacobian transpose not defined for this system. | |
|
virtual NOX::Abstract::Group::ReturnType | applyJacobianInverse (Teuchos::ParameterList ¶ms, const NOX::Abstract::Vector &input, NOX::Abstract::Vector &result) const |
| Applies Jacobian inverse for extended system. | |
|
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. | |
| virtual bool | isF () const |
Return true if 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 residual. | |
| virtual double | getNormF () const |
| Return 2-norm of extended residual. | |
|
virtual const NOX::Abstract::Vector & | getGradient () const |
| Gradient is never 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. | |
|
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getFPtr () const |
| Return extended residual. | |
|
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getGradientPtr () const |
| Gradient is never valid. | |
|
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getNewtonPtr () const |
| Return extended Newton direction. | |
| virtual double | getNormNewtonSolveResidual () const |
| Returns 2-norm of extended Newton solve residual. | |
|
virtual Teuchos::RCP< const LOCA::MultiContinuation::AbstractGroup > | getUnderlyingGroup () const |
| Return underlying group. | |
|
virtual Teuchos::RCP < LOCA::MultiContinuation::AbstractGroup > | getUnderlyingGroup () |
| Return underlying group. | |
| virtual int | getNumParams () const |
| Returns number of parameters. | |
| 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 NOX::Abstract::Group::ReturnType | computePredictor () |
| Compute predictor directions. | |
| virtual bool | isPredictor () const |
| Is Predictor valid. | |
| virtual void | scaleTangent () |
| Scales tangent to predictor. | |
| virtual void | setPredictorTangentDirection (const LOCA::MultiContinuation::ExtendedVector &v, int i) |
| Sets tangent to predictor. More... | |
|
virtual const LOCA::MultiContinuation::ExtendedMultiVector & | getPredictorTangent () const |
| Returns tangent to predictor. | |
|
virtual const LOCA::MultiContinuation::ExtendedMultiVector & | getScaledPredictorTangent () const |
| Returns scaled tangent to predictor. | |
| virtual void | setPrevX (const NOX::Abstract::Vector &y) |
| Set the previous solution vector y. | |
|
virtual const LOCA::MultiContinuation::ExtendedVector & | getPrevX () const |
| Gets the previous solution vector. | |
| virtual void | setStepSize (double deltaS, int i=0) |
| Set step size for continuation constraint equation i. | |
| virtual double | getStepSize (int i=0) const |
| Get step size for continuation constraint equation i. | |
| virtual void | setContinuationParameter (double val, int i=0) |
| Sets the value for continuation parameter i. | |
| virtual double | getContinuationParameter (int i=0) const |
| Returns the value for continuation parameter i. | |
| virtual int | getContinuationParameterID (int i=0) const |
| Get the continuation parameter id for parameter i. | |
| virtual const std::vector< int > & | getContinuationParameterIDs () const |
| Get the continuation parameter ids. | |
| virtual std::string | getContinuationParameterName (int i=0) const |
| Get the continuation parameter id for parameter i. | |
| virtual double | getStepSizeScaleFactor (int i=0) const |
| Returns step size scale factor for constraint equation i. | |
| virtual void | printSolution () const |
| Prints the group. | |
| virtual double | computeScaledDotProduct (const NOX::Abstract::Vector &x, const NOX::Abstract::Vector &y) const |
| Computes a scaled dot product between two continuation vectors. | |
| virtual int | projectToDrawDimension () const |
| Returns dimension of project to draw array. | |
| virtual void | projectToDraw (const LOCA::MultiContinuation::ExtendedVector &x, double *px) const |
| Fills the project to draw array. | |
| 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. | |
Public Member Functions inherited from LOCA::MultiContinuation::AbstractStrategy | |
| AbstractStrategy () | |
| Constructor. | |
| virtual | ~AbstractStrategy () |
| Destructor. | |
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::BorderedSystem::AbstractGroup | |
| AbstractGroup () | |
| Constructor. | |
| virtual | ~AbstractGroup () |
| Destructor. | |
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... | |
Protected Member Functions inherited from LOCA::MultiContinuation::ExtendedGroup | |
| ExtendedGroup (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &continuationParams, const Teuchos::RCP< LOCA::MultiContinuation::AbstractGroup > &grp, const Teuchos::RCP< LOCA::MultiPredictor::AbstractStrategy > &pred, const std::vector< int > ¶mIDs) | |
| Constructor used by derived classes. | |
| virtual void | setConstraints (const Teuchos::RCP< LOCA::MultiContinuation::ConstraintInterface > &constraints, bool skip_dfdp) |
| Set constraint object. More... | |
Protected Attributes inherited from LOCA::MultiContinuation::ExtendedGroup | |
| Teuchos::RCP< LOCA::GlobalData > | globalData |
| Pointer LOCA global data object. | |
|
Teuchos::RCP < LOCA::Parameter::SublistParser > | parsedParams |
| Parsed top-level parameters. | |
|
Teuchos::RCP < Teuchos::ParameterList > | continuationParams |
| Continuation parameter list. | |
|
Teuchos::RCP < LOCA::MultiContinuation::AbstractGroup > | grpPtr |
| Pointer to underlying group. | |
|
Teuchos::RCP < LOCA::MultiPredictor::AbstractStrategy > | predictor |
| Pointer to predictor object. | |
|
Teuchos::RCP < LOCA::MultiContinuation::ConstrainedGroup > | conGroup |
| Pointer to constrained group implementation. | |
| int | numParams |
| Number of parameters. | |
| LOCA::MultiContinuation::ExtendedMultiVector | tangentMultiVec |
| Stores the tangent to the predictor. | |
| LOCA::MultiContinuation::ExtendedMultiVector | scaledTangentMultiVec |
| Stores the scaled tangent to the predictor. | |
| LOCA::MultiContinuation::ExtendedVector | prevXVec |
| Stores the previous extended solution vector. | |
| std::vector< int > | conParamIDs |
| integer id of continuation parameters | |
| std::vector< double > | stepSize |
| continuation step size | |
| std::vector< double > | stepSizeScaleFactor |
| step size scale factors | |
| bool | isValidPredictor |
| Is Predictor vector valid. | |
| bool | baseOnSecant |
| Flag indicating whether to base predictor direction on secant. | |
Specialization of LOCA::MultiContinuation::ExtendedGroup to natural continuation.
Natural continuation corresponds to a continuation equation
with
given by
where
is the parameter component of the predictor direction
. This corresponds geometrically to constraining the nonlinear solver steps used in calculating
to be orthogonal to the parameter axis. The natural constraint
is represented by a LOCA::MultiContinuation::NaturalConstraint object.
| LOCA::MultiContinuation::NaturalGroup::NaturalGroup | ( | const Teuchos::RCP< LOCA::GlobalData > & | global_data, |
| const Teuchos::RCP< LOCA::Parameter::SublistParser > & | topParams, | ||
| const Teuchos::RCP< Teuchos::ParameterList > & | continuationParams, | ||
| const Teuchos::RCP< LOCA::MultiContinuation::AbstractGroup > & | grp, | ||
| const Teuchos::RCP< LOCA::MultiPredictor::AbstractStrategy > & | pred, | ||
| const std::vector< int > & | paramIDs | ||
| ) |
Constructor.
| global_data | [in] Global data object |
| topParams | [in] Parsed top-level parameter list. |
| continuationParams | [in] Continuation parameters. |
| grp | [in] Group representing . |
| pred | [in] Predictor strategy. |
| paramIDs | [in] Parameter IDs of continuation parameters. |
References Teuchos::ParameterList::get(), LOCA::MultiContinuation::ExtendedGroup::globalData, Teuchos::rcp(), and LOCA::MultiContinuation::ExtendedGroup::setConstraints().
1.8.5