Base class for Sequantial Quadratic Programming (SQP) methods. More...
#include <SQPBase.H>


Public Member Functions | |
| TypeName ("SQPBase") | |
| Runtime type information. | |
| SQPBase (const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const updateMethod &UpdateMethod, const word &type) | |
| Construct from components. | |
| virtual | ~SQPBase ()=default |
| Destructor. | |
| virtual bool | addToFile (Ostream &os) const |
| Write continuation info. | |
| virtual bool | writeMeritFunction (const updateMethod &UpdateMethod) |
| Write info about the merit function. | |
| Public Member Functions inherited from constrainedOptimisationMethod | |
| TypeName ("constrainedOptimisationMethod") | |
| Runtime type information. | |
| declareRunTimeSelectionTable (autoPtr, constrainedOptimisationMethod, dictionary,(const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints, const word &type),(mesh, dict, designVars, nConstraints, type)) | |
| constrainedOptimisationMethod (const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints, const word &type) | |
| Construct from components. | |
| virtual | ~constrainedOptimisationMethod ()=default |
| Destructor. | |
Protected Member Functions | |
| virtual scalar | meritFunctionConstraintPart () const =0 |
| Get the part the merit function that depends on the constraints. | |
Protected Attributes | |
| scalarField | LagrangianDerivatives_ |
| Derivatives of the Lagrangian function. | |
| List< scalarField > | constraintDerivativesOld_ |
| The previous constraint derivatives. | |
| scalarField | lamdas_ |
| Lagrange multipliers. | |
| fileName | objFunctionFolder_ |
| Name of the objective folder. | |
| autoPtr< OFstream > | meritFunctionFile_ |
| File including the l1 merit function. | |
| scalar | mu_ |
| Penalty value for the merit function. | |
| scalar | delta_ |
| Safety factor. | |
Base class for Sequantial Quadratic Programming (SQP) methods.
| SQPBase | ( | const fvMesh & | mesh, |
| const dictionary & | dict, | ||
| autoPtr< designVariables > & | designVars, | ||
| const updateMethod & | UpdateMethod, | ||
| const word & | type ) |
Construct from components.
Definition at line 37 of file SQPBase.C.
References constraintDerivativesOld_, delta_, dict, forAll, dictionary::found(), found, LagrangianDerivatives_, lamdas_, UPstream::master(), meritFunctionFile_, mesh, Foam::mkDir(), mu_, Foam::name(), objFunctionFolder_, timeName, and Foam::Zero.

|
virtualdefault |
|
protectedpure virtual |
Get the part the merit function that depends on the constraints.
Implemented in ISQP.
Referenced by writeMeritFunction().

|
virtual |
Write continuation info.
Definition at line 102 of file SQPBase.C.
References constraintDerivativesOld_, forAll, lamdas_, Foam::name(), and os().
Referenced by ISQP::writeData(), and SQP::writeData().


|
virtual |
Write info about the merit function.
Definition at line 115 of file SQPBase.C.
References IOstream::defaultPrecision(), Foam::endl(), forAll, updateMethod::getConstraintValues(), updateMethod::getCycle(), updateMethod::getObjectiveValue(), lamdas_, UPstream::master(), meritFunctionConstraintPart(), meritFunctionFile_, mu_, objFunctionFolder_, and Foam::setw().
Referenced by ISQP::writeAuxiliaryData(), and SQP::writeAuxiliaryData().


|
protected |
Derivatives of the Lagrangian function.
Definition at line 61 of file SQPBase.H.
Referenced by ISQP::computeCorrection(), SQP::computeCorrection(), SQPBase(), and ISQP::updateYS().
|
protected |
The previous constraint derivatives.
Definition at line 66 of file SQPBase.H.
Referenced by addToFile(), SQPBase(), ISQP::storeOldFields(), and ISQP::updateYS().
|
protected |
Lagrange multipliers.
Definition at line 71 of file SQPBase.H.
Referenced by addToFile(), ISQP::computeMeritFunction(), SQP::computeMeritFunction(), ISQP::computeNewtonDirection(), ISQP::computeRHSForDeltaX(), ISQP::diagPreconditioner(), ISQP::initialize(), ISQP::invHFL(), ISQP::lineSearch(), ISQP::matrixVectorProduct(), ISQP::resFExtraVars(), ISQP::resFL(), ISQP::resFlamda(), ISQP::ShermanMorrisonPrecon(), ISQP::solveDeltaPEqn(), ISQP::solveSubproblem(), SQPBase(), ISQP::updateSolution(), ISQP::updateYS(), and writeMeritFunction().
|
protected |
Name of the objective folder.
Definition at line 76 of file SQPBase.H.
Referenced by SQPBase(), and writeMeritFunction().
File including the l1 merit function.
Definition at line 81 of file SQPBase.H.
Referenced by SQPBase(), and writeMeritFunction().
|
protected |
Penalty value for the merit function.
Definition at line 86 of file SQPBase.H.
Referenced by ISQP::computeMeritFunction(), SQP::computeMeritFunction(), ISQP::meritFunctionDirectionalDerivative(), SQP::meritFunctionDirectionalDerivative(), SQPBase(), and writeMeritFunction().
|
protected |
Safety factor.
Definition at line 91 of file SQPBase.H.
Referenced by ISQP::computeMeritFunction(), SQP::computeMeritFunction(), and SQPBase().