NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Member Functions | List of all members
LOCA::MultiPredictor::AbstractStrategy Class Referenceabstract

Abstract interface class for predictor strategies. More...

#include <LOCA_MultiPredictor_AbstractStrategy.H>

Inheritance diagram for LOCA::MultiPredictor::AbstractStrategy:
Inheritance graph
[legend]

Public Member Functions

 AbstractStrategy ()
 Constructor.
 
virtual ~AbstractStrategy ()
 Destructor.
 
virtual
LOCA::MultiPredictor::AbstractStrategy
operator= (const LOCA::MultiPredictor::AbstractStrategy &source)=0
 Assignment operator.
 
virtual Teuchos::RCP
< LOCA::MultiPredictor::AbstractStrategy
clone (NOX::CopyType type=NOX::DeepCopy) const =0
 Clone function.
 
virtual
NOX::Abstract::Group::ReturnType 
compute (bool baseOnSecant, const std::vector< double > &stepSize, LOCA::MultiContinuation::ExtendedGroup &grp, const LOCA::MultiContinuation::ExtendedVector &prevXVec, const LOCA::MultiContinuation::ExtendedVector &xVec)=0
 Compute the predictor given the current and previous solution vectors. Set baseOnSecant to false if the predictor orientation should not be based on the secant vector (first or last steps of a continuation run). More...
 
virtual
NOX::Abstract::Group::ReturnType 
evaluate (const std::vector< double > &stepSize, const LOCA::MultiContinuation::ExtendedVector &xVec, LOCA::MultiContinuation::ExtendedMultiVector &result) const =0
 Evaluate predictor with step size stepSize. More...
 
virtual
NOX::Abstract::Group::ReturnType 
computeTangent (LOCA::MultiContinuation::ExtendedMultiVector &tangent)=0
 Compute tangent to predictor and store in tangent. More...
 
virtual bool isTangentScalable () const =0
 Is the tangent vector for this predictor scalable. More...
 

Protected Member Functions

virtual void setPredictorOrientation (bool baseOnSecant, const std::vector< double > &stepSize, const LOCA::MultiContinuation::ExtendedGroup &grp, const LOCA::MultiContinuation::ExtendedVector &prevXVec, const LOCA::MultiContinuation::ExtendedVector &xVec, LOCA::MultiContinuation::ExtendedVector &secant, LOCA::MultiContinuation::ExtendedMultiVector &tangent)
 Sets orientation of predictor based on parameter change from previous steps. More...
 

Detailed Description

Abstract interface class for predictor strategies.

AbstractStrategy defines an abstract interface for predictor strategies. It is used by the LOCA::Stepper and LOCA::MultiContinuation groups to compute predictors and tangent approximations to continuation curves and surfaces.

The interface defines several pure virtual methods that derived classes should implement for a particular predictor strategy. Note that predictor strategies are assumed to have a state, and therefore need to define copy constructors, assignment operators, and a clone function. Constructors of derived classes should be of the form:

class Derived : public AbstractStrategy {
public:
Derived(
const Teuchos::RCP<LOCA::GlobalData>& global_data,
const Teuchos::RCP<Teuchos::ParameterList>& predictorParams);
...
};

where global_data is the LOCA global data object, topParams is the parsed top-level parameter list, and predictorParams is a parameter list of predictor parameters.

This class and its children follow the Strategy pattern as defined in Erich Gamma, et al. "Design Patterns: Elements of Reusable Object-Oriented Software." Addison Wesley, Boston, MA, 1995.

Member Function Documentation

virtual NOX::Abstract::Group::ReturnType LOCA::MultiPredictor::AbstractStrategy::compute ( bool  baseOnSecant,
const std::vector< double > &  stepSize,
LOCA::MultiContinuation::ExtendedGroup grp,
const LOCA::MultiContinuation::ExtendedVector prevXVec,
const LOCA::MultiContinuation::ExtendedVector xVec 
)
pure virtual

Compute the predictor given the current and previous solution vectors. Set baseOnSecant to false if the predictor orientation should not be based on the secant vector (first or last steps of a continuation run).

As an example for a first-order predictor, this method should compute the approximate tangent to the continuation curve.

Implemented in LOCA::MultiPredictor::Secant, LOCA::MultiPredictor::Random, LOCA::MultiPredictor::Tangent, LOCA::MultiPredictor::Restart, and LOCA::MultiPredictor::Constant.

virtual NOX::Abstract::Group::ReturnType LOCA::MultiPredictor::AbstractStrategy::computeTangent ( LOCA::MultiContinuation::ExtendedMultiVector tangent)
pure virtual

Compute tangent to predictor and store in tangent.

For a first-order predictor, this is the predictor direction itself.

Implemented in LOCA::MultiPredictor::Secant, LOCA::MultiPredictor::Random, LOCA::MultiPredictor::Tangent, LOCA::MultiPredictor::Restart, and LOCA::MultiPredictor::Constant.

virtual NOX::Abstract::Group::ReturnType LOCA::MultiPredictor::AbstractStrategy::evaluate ( const std::vector< double > &  stepSize,
const LOCA::MultiContinuation::ExtendedVector xVec,
LOCA::MultiContinuation::ExtendedMultiVector result 
) const
pure virtual

Evaluate predictor with step size stepSize.

For a first-order predictor, this method should compute result[i] = xVec[i] + stepSize[i] * v[i] for each i, where v[i] is the ith predictor direction.

Implemented in LOCA::MultiPredictor::Secant, LOCA::MultiPredictor::Random, LOCA::MultiPredictor::Tangent, LOCA::MultiPredictor::Restart, and LOCA::MultiPredictor::Constant.

virtual bool LOCA::MultiPredictor::AbstractStrategy::isTangentScalable ( ) const
pure virtual

Is the tangent vector for this predictor scalable.

This method determines whether the approximate tangent computed by this strategy is appropriate for scaling.

Implemented in LOCA::MultiPredictor::Secant, LOCA::MultiPredictor::Random, LOCA::MultiPredictor::Tangent, LOCA::MultiPredictor::Restart, and LOCA::MultiPredictor::Constant.

void LOCA::MultiPredictor::AbstractStrategy::setPredictorOrientation ( bool  baseOnSecant,
const std::vector< double > &  stepSize,
const LOCA::MultiContinuation::ExtendedGroup grp,
const LOCA::MultiContinuation::ExtendedVector prevXVec,
const LOCA::MultiContinuation::ExtendedVector xVec,
LOCA::MultiContinuation::ExtendedVector secant,
LOCA::MultiContinuation::ExtendedMultiVector tangent 
)
protectedvirtual

Sets orientation of predictor based on parameter change from previous steps.

The implementation here looks at the sign of the scaled dot product between the secant vector and the predictor and changes the sign of the predictor if this scaled dot product is a different sign than the step size. If baseOnSecant is false, then the sign is chosen so the parameter component of the predictor is positive for cases when a secant vector is not available (first step in a continuation run) or may give incorrect information (last step of a continuation run).

References LOCA::MultiContinuation::ExtendedGroup::computeScaledDotProduct(), LOCA::Extended::MultiVector::getScalar(), LOCA::Extended::MultiVector::numVectors(), LOCA::Extended::MultiVector::scale(), and LOCA::Extended::Vector::update().


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