Abstract base class for optimisation methods. More...
Public Member Functions | |
TypeName ("updateMethod") | |
declareRunTimeSelectionTable (autoPtr, updateMethod, dictionary,(const fvMesh &mesh, const dictionary &dict),(mesh, dict)) | |
updateMethod (const fvMesh &mesh, const dictionary &dict) | |
virtual | ~updateMethod ()=default |
void | setObjectiveDeriv (const scalarField &derivs) |
void | setConstraintDeriv (const PtrList< scalarField > &derivs) |
void | setObjectiveValue (const scalar value) |
void | setConstraintValues (const scalarField &values) |
void | setStep (const scalar eta) |
void | setGlobalSum (const bool useGlobalSum) |
virtual void | computeCorrection ()=0 |
scalarField & | returnCorrection () |
void | writeCorrection () |
virtual scalar | computeMeritFunction () |
virtual scalar | meritFunctionDirectionalDerivative () |
bool & | initialEtaSet () |
virtual void | updateOldCorrection (const scalarField &oldCorrection) |
virtual void | write () |
Static Public Member Functions | |
static autoPtr< updateMethod > | New (const fvMesh &mesh, const dictionary &dict) |
Protected Member Functions | |
const scalarField | leftMult (const scalarField &, const SquareMatrix< scalar > &) |
const scalarField | rightMult (const SquareMatrix< scalar > &, const scalarField &) |
SquareMatrix< scalar > | outerProd (const scalarField &, const scalarField &) |
SquareMatrix< scalar > | inv (SquareMatrix< scalar > A) |
scalar | globalSum (const scalarField &field) |
scalar | globalSum (tmp< scalarField > &tfield) |
dictionary | coeffsDict () |
Protected Attributes | |
const fvMesh & | mesh_ |
const dictionary | dict_ |
IOdictionary | optMethodIODict_ |
scalarField | objectiveDerivatives_ |
PtrList< scalarField > | constraintDerivatives_ |
scalar | objectiveValue_ |
scalarField | cValues_ |
scalarField | correction_ |
scalarField | cumulativeCorrection_ |
scalar | eta_ |
bool | initialEtaSet_ |
word | correctionFolder_ |
bool | globalSum_ |
Abstract base class for optimisation methods.
Definition at line 50 of file updateMethod.H.
updateMethod | ( | const fvMesh & | mesh, |
const dictionary & | dict | ||
) |
Definition at line 194 of file updateMethod.C.
References dict, UPstream::master(), Foam::mkDir(), and dictionary::readIfPresent().
|
virtualdefault |
|
protected |
Definition at line 38 of file updateMethod.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, Matrix::n(), s, and Foam::Zero.
|
protected |
Definition at line 65 of file updateMethod.C.
References Foam::abort(), Foam::FatalError, FatalErrorInFunction, forAll, Matrix::n(), s, and Foam::Zero.
|
protected |
Definition at line 92 of file updateMethod.C.
References Foam::abort(), Foam::constant::physicoChemical::b, Foam::FatalError, FatalErrorInFunction, forAll, and Foam::Zero.
|
protected |
Definition at line 118 of file updateMethod.C.
References A, DebugInfo, Foam::endl(), invA(), Foam::LUBacksubstitute(), Foam::LUDecompose(), n, and Foam::Zero.
Referenced by constraintProjection::computeCorrection().
|
protected |
Definition at line 168 of file updateMethod.C.
References field(), Foam::gSum(), and Foam::sum().
Referenced by constraintProjection::computeCorrection().
|
protected |
Definition at line 183 of file updateMethod.C.
References tmp::clear().
|
protected |
Definition at line 247 of file updateMethod.C.
References dictionary::subOrEmptyDict(), and Foam::type().
Referenced by SQP::SQP(), and SR1::SR1().
TypeName | ( | "updateMethod" | ) |
declareRunTimeSelectionTable | ( | autoPtr | , |
updateMethod | , | ||
dictionary | , | ||
(const fvMesh &mesh, const dictionary &dict) | , | ||
(mesh, dict) | |||
) |
|
static |
Definition at line 256 of file updateMethod.C.
References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::get(), Foam::Info, and mesh.
void setObjectiveDeriv | ( | const scalarField & | derivs | ) |
Definition at line 284 of file updateMethod.C.
void setConstraintDeriv | ( | const PtrList< scalarField > & | derivs | ) |
Definition at line 291 of file updateMethod.C.
void setObjectiveValue | ( | const scalar | value | ) |
Definition at line 299 of file updateMethod.C.
void setConstraintValues | ( | const scalarField & | values | ) |
Definition at line 305 of file updateMethod.C.
References Foam::HashTableOps::values().
void setStep | ( | const scalar | eta | ) |
Definition at line 311 of file updateMethod.C.
void setGlobalSum | ( | const bool | useGlobalSum | ) |
Definition at line 317 of file updateMethod.C.
|
pure virtual |
Implemented in constrainedOptimisationMethod, SQP, LBFGS, DBFGS, BFGS, SR1, conjugateGradient, constraintProjection, and steepestDescent.
Foam::scalarField & returnCorrection | ( | ) |
Definition at line 323 of file updateMethod.C.
void writeCorrection | ( | ) |
Definition at line 330 of file updateMethod.C.
References Foam::endl(), forAll, UPstream::master(), and Foam::Zero.
|
virtual |
Reimplemented in SQP.
Definition at line 364 of file updateMethod.C.
|
virtual |
Reimplemented in SQP.
Definition at line 370 of file updateMethod.C.
bool & initialEtaSet | ( | ) |
Definition at line 376 of file updateMethod.C.
|
virtual |
Reimplemented in SQP, LBFGS, DBFGS, BFGS, SR1, and conjugateGradient.
Definition at line 383 of file updateMethod.C.
Referenced by BFGS::updateOldCorrection(), SR1::updateOldCorrection(), DBFGS::updateOldCorrection(), LBFGS::updateOldCorrection(), and SQP::updateOldCorrection().
|
virtual |
Reimplemented in SQP, LBFGS, DBFGS, BFGS, SR1, and conjugateGradient.
Definition at line 391 of file updateMethod.C.
References IOstreamOption::ASCII.
Referenced by conjugateGradient::write(), BFGS::write(), SR1::write(), DBFGS::write(), LBFGS::write(), and SQP::write().
|
protected |
Definition at line 56 of file updateMethod.H.
|
protected |
Definition at line 58 of file updateMethod.H.
|
protected |
Definition at line 61 of file updateMethod.H.
|
protected |
Definition at line 64 of file updateMethod.H.
Referenced by conjugateGradient::allocateFields(), SR1::allocateMatrices(), DBFGS::allocateMatrices(), LBFGS::allocateMatrices(), BFGS::allocateMatrices(), steepestDescent::computeCorrection(), and constraintProjection::computeCorrection().
|
protected |
Definition at line 67 of file updateMethod.H.
Referenced by constraintProjection::computeCorrection().
|
protected |
Definition at line 70 of file updateMethod.H.
|
protected |
Definition at line 73 of file updateMethod.H.
Referenced by constraintProjection::computeCorrection().
|
protected |
Definition at line 76 of file updateMethod.H.
Referenced by steepestDescent::computeCorrection(), and constraintProjection::computeCorrection().
|
protected |
Definition at line 80 of file updateMethod.H.
|
protected |
Definition at line 83 of file updateMethod.H.
Referenced by steepestDescent::computeCorrection(), and constraintProjection::computeCorrection().
|
protected |
Definition at line 86 of file updateMethod.H.
|
protected |
Definition at line 94 of file updateMethod.H.
|
protected |
Definition at line 97 of file updateMethod.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.