Abstract base class for line search methods. More...
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 |
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 lineSearch & | operator++ () |
Increment iteration number and store merit value corresponding to the previous optimisation cycle. More... | |
virtual lineSearch & | operator++ (int) |
Postfix increment. Necessary for compilation. More... | |
Static Public Member Functions | |
static autoPtr< lineSearch > | New (const dictionary &dict, const Time &time) |
Return a reference to the selected turbulence model. More... | |
Protected Member Functions | |
const dictionary & | coeffsDict () |
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_ |
Whether to extrapolate the correction multiplier for this optimisation cycle based on the previous ones. More... | |
autoPtr< stepUpdate > | stepUpdate_ |
Mechanism to update method if line search conditions are not set. More... | |
Abstract base class for line search methods.
Definition at line 52 of file lineSearch.H.
lineSearch | ( | const dictionary & | dict, |
const Time & | time | ||
) |
Construct from components.
Definition at line 45 of file lineSearch.C.
|
virtualdefault |
Destructor.
|
protected |
Optional coeffs dict.
Definition at line 37 of file lineSearch.C.
References lineSearch::dict_, dictionary::optionalSubDict(), and Foam::type().
TypeName | ( | "lineSearch" | ) |
Runtime type information.
declareRunTimeSelectionTable | ( | autoPtr | , |
lineSearch | , | ||
dictionary | , | ||
(const dictionary &dict, const Time &time) | , | ||
(dict, time) | |||
) |
|
static |
Return a reference to the selected turbulence model.
Definition at line 89 of file lineSearch.C.
References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::getOrDefault(), Foam::Info, and autoPtr< T >::reset().
|
virtual |
Set objective derivative.
Definition at line 129 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
void setDirection | ( | const scalarField & | direction | ) |
Set direction.
Definition at line 136 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
void setNewMeritValue | ( | const scalar | value | ) |
Set new objective value.
Definition at line 142 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
void setOldMeritValue | ( | const scalar | value | ) |
Set old objective value.
Definition at line 149 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
|
virtual |
Reset step to initial value.
Definition at line 156 of file lineSearch.C.
References Foam::endl(), Foam::Info, Foam::max(), and Foam::min().
Referenced by steadyOptimisation::lineSearchUpdate().
Foam::label maxIters | ( | ) | const |
Get max number of iterations.
Definition at line 176 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
Foam::scalar step | ( | ) | const |
Get current step.
Definition at line 182 of file lineSearch.C.
Referenced by steadyOptimisation::lineSearchUpdate().
|
pure virtual |
Return the correction of the design variables.
Implemented in ArmijoConditions.
Referenced by steadyOptimisation::lineSearchUpdate().
|
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().
|
virtual |
Update the step using the supplied value.
Definition at line 188 of file lineSearch.C.
|
virtual |
Increment iteration number and store merit value corresponding to the previous optimisation cycle.
Definition at line 194 of file lineSearch.C.
References IOstreamOption::ASCII.
|
virtual |
Postfix increment. Necessary for compilation.
Definition at line 210 of file lineSearch.C.
|
protected |
Subdict within updateMethod.
Definition at line 61 of file lineSearch.H.
Referenced by lineSearch::coeffsDict().
|
protected |
IOdictionary under time/uniform for continuation.
Definition at line 66 of file lineSearch.H.
|
protected |
Directional derivative of the merit function.
Definition at line 71 of file lineSearch.H.
Referenced by ArmijoConditions::converged().
|
protected |
Update direction.
Definition at line 76 of file lineSearch.H.
|
protected |
Old merit value from this opt cycle.
Definition at line 81 of file lineSearch.H.
Referenced by ArmijoConditions::converged().
|
protected |
New merit value from this opt cycle.
Definition at line 86 of file lineSearch.H.
Referenced by ArmijoConditions::converged().
|
protected |
Merit directional deriv from the previous opt cycle.
Definition at line 91 of file lineSearch.H.
|
protected |
Correction multiplier at the first step of line search.
Definition at line 96 of file lineSearch.H.
|
protected |
Minimum allowed correction multiplier.
Definition at line 101 of file lineSearch.H.
|
protected |
Correction multiplier.
Definition at line 106 of file lineSearch.H.
Referenced by ArmijoConditions::converged().
|
protected |
Inner line search iteration.
Definition at line 111 of file lineSearch.H.
|
protected |
Maximum line search iterations.
Definition at line 116 of file lineSearch.H.
|
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 124 of file lineSearch.H.
|
protected |
Mechanism to update method if line search conditions are not set.
Definition at line 129 of file lineSearch.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.