lineSearch Class Referenceabstract

Abstract base class for line search methods. More...

Inheritance diagram for lineSearch:
Collaboration diagram for lineSearch:

Public Member Functions

 TypeName ("lineSearch")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, lineSearch, dictionary,(const dictionary &dict, const Time &time, updateMethod &UpdateMethod),(dict, time, UpdateMethod))
 
 lineSearch (const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
 Construct from components. More...
 
virtual ~lineSearch ()=default
 Destructor. More...
 
virtual void setDeriv (const scalar deriv)
 Set directional derivative. More...
 
virtual void setNewDeriv (const scalar deriv)
 Set new directional derivative. More...
 
void setDirection (const scalarField &direction)
 Set direction. More...
 
void setNewMeritValue (const scalar value)
 Set new objective value. More...
 
void setOldMeritValue (const scalar value)
 Set old objective value. More...
 
virtual void reset ()
 Reset step to initial value. More...
 
label innerIter () const
 Get inner line search iteration. More...
 
label maxIters () const
 Get max number of iterations. More...
 
scalar step () const
 Get current step. More...
 
virtual bool converged ()=0
 Return the correction of the design variables. More...
 
virtual void updateStep ()=0
 Update the line search step based on the specific line search strategy, e.g. bisection, quadratic fit, etc. More...
 
virtual void updateStep (const scalar newStep)
 Update the step using the supplied value. More...
 
virtual void updateCorrection (scalarField &correction)
 Update the correction. More...
 
virtual bool loop ()
 Return true if lineSearch should continue and if so increment inner. More...
 
virtual bool computeGradient () const
 Does line search need to update the gradient? More...
 
virtual void postUpdate ()
 Execute steps at the end of the line search iterations. More...
 
virtual lineSearchoperator++ ()
 Increment iteration number and store merit value corresponding to the previous optimisation cycle. More...
 
virtual lineSearchoperator++ (int)
 Postfix increment. Necessary for compilation. More...
 

Static Public Member Functions

static autoPtr< lineSearchNew (const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
 Return a reference to the selected turbulence model. More...
 

Protected Member Functions

const dictionarycoeffsDict ()
 Optional coeffs dict. More...
 

Protected Attributes

const dictionary dict_
 Subdict within updateMethod. More...
 
IOdictionary lineSearchDict_
 IOdictionary under time/uniform for continuation. More...
 
scalar directionalDeriv_
 Directional derivative of the merit function. More...
 
scalarField direction_
 Update direction. More...
 
scalar oldMeritValue_
 Old merit value from this opt cycle. More...
 
scalar newMeritValue_
 New merit value from this opt cycle. More...
 
scalar prevMeritDeriv_
 Merit directional deriv from the previous opt cycle. More...
 
scalar initialStep_
 Correction multiplier at the first step of line search. More...
 
scalar minStep_
 Minimum allowed correction multiplier. More...
 
scalar step_
 Correction multiplier. More...
 
label iter_
 Optimisation cycle. More...
 
label innerIter_
 Inner line search iteration. More...
 
label maxIters_
 Maximum line search iterations. More...
 
bool extrapolateInitialStep_
 Whether to extrapolate the correction multiplier for this optimisation cycle based on the previous ones. More...
 
autoPtr< stepUpdatestepUpdate_
 Mechanism to update method if line search conditions are not set. More...
 
updateMethodupdateMethod_
 Reference to the update method related to the line search. More...
 

Detailed Description

Abstract base class for line search methods.

Source files

Definition at line 53 of file lineSearch.H.

Constructor & Destructor Documentation

◆ lineSearch()

lineSearch ( const dictionary dict,
const Time time,
updateMethod UpdateMethod 
)

Construct from components.

Definition at line 46 of file lineSearch.C.

◆ ~lineSearch()

virtual ~lineSearch ( )
virtualdefault

Destructor.

Member Function Documentation

◆ coeffsDict()

const Foam::dictionary & coeffsDict ( )
protected

Optional coeffs dict.

Definition at line 37 of file lineSearch.C.

References lineSearch::dict_, dictionary::optionalSubDict(), and Foam::type().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "lineSearch"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
lineSearch  ,
dictionary  ,
(const dictionary &dict, const Time &time, updateMethod &UpdateMethod)  ,
(dict, time, UpdateMethod)   
)

◆ New()

Foam::autoPtr< Foam::lineSearch > New ( const dictionary dict,
const Time time,
updateMethod UpdateMethod 
)
static

Return a reference to the selected turbulence model.

Definition at line 96 of file lineSearch.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::getOrDefault(), Foam::Info, and autoPtr< T >::reset().

Here is the call graph for this function:

◆ setDeriv()

void setDeriv ( const scalar  deriv)
virtual

Set directional derivative.

Definition at line 137 of file lineSearch.C.

◆ setNewDeriv()

void setNewDeriv ( const scalar  deriv)
virtual

Set new directional derivative.

Does nothing in base. Only used by methods that require the gradient information to be computed at the new point

Definition at line 144 of file lineSearch.C.

◆ setDirection()

void setDirection ( const scalarField direction)
inline

Set direction.

Definition at line 239 of file lineSearch.H.

References lineSearch::direction_.

◆ setNewMeritValue()

void setNewMeritValue ( const scalar  value)

Set new objective value.

Definition at line 150 of file lineSearch.C.

◆ setOldMeritValue()

void setOldMeritValue ( const scalar  value)

Set old objective value.

Definition at line 157 of file lineSearch.C.

◆ reset()

void reset ( )
virtual

Reset step to initial value.

Definition at line 164 of file lineSearch.C.

References Foam::clamp(), Foam::endl(), and Foam::Info.

Here is the call graph for this function:

◆ innerIter()

label innerIter ( ) const
inline

Get inner line search iteration.

Definition at line 262 of file lineSearch.H.

References lineSearch::innerIter_.

◆ maxIters()

label maxIters ( ) const
inline

Get max number of iterations.

Definition at line 270 of file lineSearch.H.

References lineSearch::maxIters_.

◆ step()

scalar step ( ) const
inline

Get current step.

Definition at line 278 of file lineSearch.H.

References lineSearch::step_.

◆ converged()

virtual bool converged ( )
pure virtual

Return the correction of the design variables.

Implemented in GCMMA, and ArmijoConditions.

◆ updateStep() [1/2]

virtual void updateStep ( )
pure virtual

Update the line search step based on the specific line search strategy, e.g. bisection, quadratic fit, etc.

Implemented in GCMMA, and ArmijoConditions.

◆ updateStep() [2/2]

void updateStep ( const scalar  newStep)
virtual

Update the step using the supplied value.

Definition at line 186 of file lineSearch.C.

◆ updateCorrection()

void updateCorrection ( scalarField correction)
virtual

Update the correction.

Multiplies with step in base

Reimplemented in GCMMA.

Definition at line 192 of file lineSearch.C.

References Foam::correction().

Here is the call graph for this function:

◆ loop()

bool loop ( )
virtual

Return true if lineSearch should continue and if so increment inner.

iteration

Definition at line 198 of file lineSearch.C.

◆ computeGradient()

bool computeGradient ( ) const
virtual

Does line search need to update the gradient?

Definition at line 211 of file lineSearch.C.

◆ postUpdate()

void postUpdate ( )
virtual

Execute steps at the end of the line search iterations.

Definition at line 217 of file lineSearch.C.

◆ operator++() [1/2]

Foam::lineSearch & operator++ ( )
virtual

Increment iteration number and store merit value corresponding to the previous optimisation cycle.

Definition at line 223 of file lineSearch.C.

References IOstreamOption::ASCII.

◆ operator++() [2/2]

Foam::lineSearch & operator++ ( int  )
virtual

Postfix increment. Necessary for compilation.

Definition at line 242 of file lineSearch.C.

Member Data Documentation

◆ dict_

const dictionary dict_
protected

Subdict within updateMethod.

Definition at line 62 of file lineSearch.H.

Referenced by lineSearch::coeffsDict().

◆ lineSearchDict_

IOdictionary lineSearchDict_
protected

IOdictionary under time/uniform for continuation.

Definition at line 67 of file lineSearch.H.

◆ directionalDeriv_

scalar directionalDeriv_
protected

Directional derivative of the merit function.

Definition at line 72 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ direction_

scalarField direction_
protected

Update direction.

Definition at line 77 of file lineSearch.H.

Referenced by lineSearch::setDirection().

◆ oldMeritValue_

scalar oldMeritValue_
protected

Old merit value from this opt cycle.

Definition at line 82 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ newMeritValue_

scalar newMeritValue_
protected

New merit value from this opt cycle.

Definition at line 87 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ prevMeritDeriv_

scalar prevMeritDeriv_
protected

Merit directional deriv from the previous opt cycle.

Definition at line 92 of file lineSearch.H.

◆ initialStep_

scalar initialStep_
protected

Correction multiplier at the first step of line search.

Definition at line 97 of file lineSearch.H.

◆ minStep_

scalar minStep_
protected

Minimum allowed correction multiplier.

Definition at line 102 of file lineSearch.H.

◆ step_

scalar step_
protected

Correction multiplier.

Definition at line 107 of file lineSearch.H.

Referenced by ArmijoConditions::converged(), and lineSearch::step().

◆ iter_

label iter_
protected

Optimisation cycle.

Definition at line 112 of file lineSearch.H.

Referenced by GCMMA::writeToFiles().

◆ innerIter_

label innerIter_
protected

Inner line search iteration.

Definition at line 117 of file lineSearch.H.

Referenced by lineSearch::innerIter(), and GCMMA::writeToFiles().

◆ maxIters_

label maxIters_
protected

Maximum line search iterations.

Definition at line 122 of file lineSearch.H.

Referenced by lineSearch::maxIters().

◆ extrapolateInitialStep_

bool extrapolateInitialStep_
protected

Whether to extrapolate the correction multiplier for this optimisation cycle based on the previous ones.

Useful for non-quasi Newton methods

Definition at line 130 of file lineSearch.H.

◆ stepUpdate_

autoPtr<stepUpdate> stepUpdate_
protected

Mechanism to update method if line search conditions are not set.

Definition at line 135 of file lineSearch.H.

◆ updateMethod_

updateMethod& updateMethod_
protected

Reference to the update method related to the line search.

Definition at line 140 of file lineSearch.H.


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