lineSearch Class Referenceabstract

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

Inheritance diagram for lineSearch:
[legend]
Collaboration diagram for lineSearch:
[legend]

Public Member Functions

 TypeName ("lineSearch")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, lineSearch, dictionary,(const dictionary &dict, const Time &time),(dict, time))
 
 lineSearch (const dictionary &dict, const Time &time)
 Construct from components. More...
 
virtual ~lineSearch ()=default
 Destructor. More...
 
virtual void setDeriv (const scalar deriv)
 Set objective 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 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
 
virtual void updateStep (const scalar newStep)
 Update the step using the supplied value. More...
 
virtual lineSearchoperator++ ()
 
virtual lineSearchoperator++ (int)
 Postfix increment. Necessary for compilation. More...
 

Static Public Member Functions

static autoPtr< lineSearchNew (const dictionary &dict, const Time &time)
 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_
 Inner line search iteration. More...
 
label maxIters_
 Maximum line search iterations. More...
 
bool extrapolateInitialStep_
 
autoPtr< stepUpdatestepUpdate_
 Mechanism to update method if line search conditions are not set. More...
 

Detailed Description

Abstract base class for line search methods.

Source files

Definition at line 56 of file lineSearch.H.

Constructor & Destructor Documentation

◆ lineSearch()

lineSearch ( const dictionary dict,
const Time time 
)

Construct from components.

Definition at line 52 of file lineSearch.C.

◆ ~lineSearch()

virtual ~lineSearch ( )
virtualdefault

Destructor.

Member Function Documentation

◆ coeffsDict()

const Foam::dictionary & coeffsDict ( )
protected

Optional coeffs dict.

Definition at line 44 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)  ,
(dict, time)   
)

◆ New()

Foam::autoPtr< Foam::lineSearch > New ( const dictionary dict,
const Time time 
)
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 objective derivative.

Definition at line 136 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ setDirection()

void setDirection ( const scalarField direction)

Set direction.

Definition at line 143 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ setNewMeritValue()

void setNewMeritValue ( const scalar  value)

Set new objective value.

Definition at line 149 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ setOldMeritValue()

void setOldMeritValue ( const scalar  value)

Set old objective value.

Definition at line 156 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ reset()

void reset ( )
virtual

Reset step to initial value.

Definition at line 163 of file lineSearch.C.

References Foam::endl(), Foam::Info, Foam::max(), and Foam::min().

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ maxIters()

Foam::label maxIters ( ) const

Get max number of iterations.

Definition at line 183 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ step()

Foam::scalar step ( ) const

Get current step.

Definition at line 189 of file lineSearch.C.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ converged()

virtual bool converged ( )
pure virtual

Return the correction of the design variables.

Implemented in ArmijoConditions.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ 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 ArmijoConditions.

Referenced by steadyOptimisation::lineSearchUpdate().

Here is the caller graph for this function:

◆ updateStep() [2/2]

void updateStep ( const scalar  newStep)
virtual

Update the step using the supplied value.

Definition at line 195 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 201 of file lineSearch.C.

References IOstreamOption::ASCII.

◆ operator++() [2/2]

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

Postfix increment. Necessary for compilation.

Definition at line 217 of file lineSearch.C.

Member Data Documentation

◆ dict_

const dictionary dict_
protected

Subdict within updateMethod.

Definition at line 63 of file lineSearch.H.

Referenced by lineSearch::coeffsDict().

◆ lineSearchDict_

IOdictionary lineSearchDict_
protected

IOdictionary under time/uniform for continuation.

Definition at line 66 of file lineSearch.H.

◆ directionalDeriv_

scalar directionalDeriv_
protected

Directional derivative of the merit function.

Definition at line 69 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ direction_

scalarField direction_
protected

Update direction.

Definition at line 72 of file lineSearch.H.

◆ oldMeritValue_

scalar oldMeritValue_
protected

Old merit value from this opt cycle.

Definition at line 75 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ newMeritValue_

scalar newMeritValue_
protected

New merit value from this opt cycle.

Definition at line 78 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ prevMeritDeriv_

scalar prevMeritDeriv_
protected

Merit directional deriv from the previous opt cycle.

Definition at line 81 of file lineSearch.H.

◆ initialStep_

scalar initialStep_
protected

Correction multiplier at the first step of line search.

Definition at line 84 of file lineSearch.H.

◆ minStep_

scalar minStep_
protected

Minimum allowed correction multiplier.

Definition at line 87 of file lineSearch.H.

◆ step_

scalar step_
protected

Correction multiplier.

Definition at line 90 of file lineSearch.H.

Referenced by ArmijoConditions::converged().

◆ iter_

label iter_
protected

Inner line search iteration.

Definition at line 93 of file lineSearch.H.

◆ maxIters_

label maxIters_
protected

Maximum line search iterations.

Definition at line 96 of file lineSearch.H.

◆ 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 101 of file lineSearch.H.

◆ stepUpdate_

autoPtr<stepUpdate> stepUpdate_
protected

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

Definition at line 104 of file lineSearch.H.


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