46 return dict_.optionalSubDict(
type() +
"Coeffs");
52Foam::lineSearch::lineSearch
73 directionalDeriv_(
Zero),
79 lineSearchDict_.getOrDefault<scalar>(
"prevMeritDeriv",
Zero)
81 initialStep_(
dict.getOrDefault<scalar>(
"initialStep", 1.)),
82 minStep_(
dict.getOrDefault<scalar>(
"minStep", 0.3)),
84 iter_(lineSearchDict_.getOrDefault<label>(
"iter", 0)),
86 maxIters_(
dict.getOrDefault<label>(
"maxIters", 4)),
87 extrapolateInitialStep_
89 dict.getOrDefault<bool>
91 "extrapolateInitialStep",
111 const word modelType(
dict.getOrDefault<
word>(
"type",
"none"));
113 Info<<
"lineSearch type : " << modelType <<
endl;
115 if (modelType !=
"none")
117 auto* ctorPtr = dictionaryConstructorTable(modelType);
126 *dictionaryConstructorTablePtr_
130 lineSrch.reset((ctorPtr(
dict, time, UpdateMethod)).ptr());
134 Info<<
"No line search method specified. "
135 <<
"Proceeding with constant step" <<
endl;
174 if (extrapolateInitialStep_ && iter_ != 0)
180 clamp(step_*prevMeritDeriv_/directionalDeriv_, minStep_, scalar(1));
181 Info<<
"\n------- Computing initial step-------" <<
endl;
182 Info<<
"old dphi(0) " << prevMeritDeriv_ <<
endl;
183 Info<<
"dphi(0) " << directionalDeriv_ <<
endl;
184 Info<<
"Setting initial step value " << step_ <<
endl <<
endl;
207 const bool isRunning = innerIter_ < maxIters_;
233 prevMeritDeriv_ = directionalDeriv_;
234 lineSearchDict_.add<scalar>(
"prevMeritDeriv", prevMeritDeriv_,
true);
235 lineSearchDict_.add<label>(
"iter", iter_,
true);
236 if (lineSearchDict_.time().writeTime())
238 lineSearchDict_.regIOobject::writeObject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A simple container for options an IOstream can normally have.
@ ASCII
"ascii" (normal default)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base class for line search methods.
const dictionary & coeffsDict()
Optional coeffs dict.
virtual void setDeriv(const scalar deriv)
Set directional derivative.
scalar newMeritValue_
New merit value from this opt cycle.
scalar oldMeritValue_
Old merit value from this opt cycle.
bool extrapolateInitialStep_
Whether to extrapolate the correction multiplier for this optimisation cycle based on the previous on...
virtual void updateCorrection(scalarField &correction)
Update the correction.
static autoPtr< lineSearch > New(const dictionary &dict, const Time &time, updateMethod &UpdateMethod)
Return a reference to the selected turbulence model.
virtual void postUpdate()
Execute steps at the end of the line search iterations.
updateMethod & updateMethod_
Reference to the update method related to the line search.
autoPtr< stepUpdate > stepUpdate_
Mechanism to update method if line search conditions are not set.
void setOldMeritValue(const scalar value)
Set old objective value.
void setNewMeritValue(const scalar value)
Set new objective value.
const dictionary dict_
Subdict within updateMethod.
virtual void setNewDeriv(const scalar deriv)
Set new directional derivative.
virtual bool computeGradient() const
Does line search need to update the gradient?
virtual void updateStep()=0
Update the line search step based on the specific line search strategy, e.g. bisection,...
label iter_
Optimisation cycle.
scalarField direction_
Update direction.
scalar step_
Correction multiplier.
scalar directionalDeriv_
Directional derivative of the merit function.
virtual lineSearch & operator++()
Increment iteration number and store merit value corresponding to the previous optimisation cycle.
scalar minStep_
Minimum allowed correction multiplier.
label maxIters_
Maximum line search iterations.
virtual void reset()
Reset step to initial value.
IOdictionary lineSearchDict_
IOdictionary under time/uniform for continuation.
scalar initialStep_
Correction multiplier at the first step of line search.
virtual bool loop()
Return true if lineSearch should continue and if so increment inner.
scalar prevMeritDeriv_
Merit directional deriv from the previous opt cycle.
label innerIter_
Inner line search iteration.
Abstract base class for step update methods used in line search.
Abstract base class for optimisation methods.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.