optimisationType Class Referenceabstract

Abstract base class for optimisation methods. More...

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

Public Member Functions

 TypeName ("optimisationType")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, optimisationType, dictionary,(fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers),(mesh, dict, adjointSolverManagers))
 
 optimisationType (fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers)
 Construct from components. More...
 
virtual ~optimisationType ()=default
 
virtual void update ()
 Update design variables. More...
 
virtual void update (scalarField &correction)
 Update design variables based on a given correction. More...
 
virtual void storeDesignVariables ()=0
 Store design variables, as the starting point for line search. More...
 
virtual void resetDesignVariables ()=0
 Reset to starting point of line search. More...
 
virtual tmp< scalarFieldcomputeDirection ()
 Compute update direction. More...
 
virtual void updateGradientsAndValues (scalarField &objectiveSens, PtrList< scalarField > &constraintSens, scalar &objectiveValue, scalarField &constraintValues)
 Compute cumulative objective and constraint gradients. More...
 
virtual scalar computeMeritFunction ()
 Compute the merit function of the optimisation problem. More...
 
virtual scalar meritFunctionDirectionalDerivative ()
 Derivative of the merit function. More...
 
virtual void updateOldCorrection (const scalarField &)
 Update old correction. Needed for quasi-Newton Methods. More...
 
virtual void write ()
 Write useful quantities to files. More...
 
const autoPtr< volScalarField > & sourcePtr ()
 Get source term. More...
 
autoPtr< lineSearch > & getLineSearch ()
 Get a reference to the line search object. More...
 

Static Public Member Functions

static autoPtr< optimisationTypeNew (fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers)
 Return a reference to the selected turbulence model. More...
 

Protected Member Functions

virtual void updateDesignVariables (scalarField &correction)=0
 Update the design variables given their correction. More...
 
virtual void computeEta (scalarField &correction)=0
 Compute eta if not set in the first step. More...
 

Protected Attributes

fvMeshmesh_
 
const dictionary dict_
 
PtrList< adjointSolverManager > & adjointSolvManagers_
 
autoPtr< updateMethodupdateMethod_
 
autoPtr< volScalarFieldsourcePtr_
 
autoPtr< lineSearchlineSearch_
 

Detailed Description

Abstract base class for optimisation methods.

Source files

Definition at line 58 of file optimisationTypeIncompressible.H.

Constructor & Destructor Documentation

◆ optimisationType()

optimisationType ( fvMesh mesh,
const dictionary dict,
PtrList< adjointSolverManager > &  adjointSolverManagers 
)

Construct from components.

Definition at line 49 of file optimisationTypeIncompressible.C.

References optimisationType::adjointSolvManagers_, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::nl, optimisationType::updateMethod_, and WarningInFunction.

Here is the call graph for this function:

◆ ~optimisationType()

virtual ~optimisationType ( )
virtualdefault

Member Function Documentation

◆ updateDesignVariables()

virtual void updateDesignVariables ( scalarField correction)
protectedpure virtual

Update the design variables given their correction.

Implemented in shapeOptimisation.

Referenced by optimisationType::update().

Here is the caller graph for this function:

◆ computeEta()

virtual void computeEta ( scalarField correction)
protectedpure virtual

Compute eta if not set in the first step.

Implemented in shapeOptimisation.

Referenced by optimisationType::computeDirection().

Here is the caller graph for this function:

◆ TypeName()

TypeName ( "optimisationType"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

◆ New()

autoPtr< optimisationType > New ( fvMesh mesh,
const dictionary dict,
PtrList< adjointSolverManager > &  adjointSolverManagers 
)
static

Return a reference to the selected turbulence model.

Definition at line 126 of file optimisationTypeIncompressible.C.

References adjointSolverManagers, dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, dictionary::get(), Foam::Info, mesh, and dictionary::subDict().

Here is the call graph for this function:

◆ update() [1/2]

void update ( )
virtual

Update design variables.

Definition at line 159 of file optimisationTypeIncompressible.C.

References optimisationType::computeDirection(), Foam::correction(), tmp< T >::ref(), optimisationType::update(), optimisationType::updateOldCorrection(), and optimisationType::write().

Referenced by optimisationType::update().

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

◆ update() [2/2]

void update ( scalarField correction)
virtual

Update design variables based on a given correction.

Definition at line 175 of file optimisationTypeIncompressible.C.

References Foam::correction(), optimisationType::lineSearch_, and optimisationType::updateDesignVariables().

Here is the call graph for this function:

◆ storeDesignVariables()

virtual void storeDesignVariables ( )
pure virtual

Store design variables, as the starting point for line search.

Implemented in shapeOptimisation.

◆ resetDesignVariables()

virtual void resetDesignVariables ( )
pure virtual

Reset to starting point of line search.

Implemented in shapeOptimisation.

◆ computeDirection()

tmp< scalarField > computeDirection ( )
virtual

Compute update direction.

Definition at line 189 of file optimisationTypeIncompressible.C.

References optimisationType::computeEta(), Foam::correction(), tmp< T >::ref(), UList< T >::size(), optimisationType::updateGradientsAndValues(), optimisationType::updateMethod_, and Foam::Zero.

Referenced by optimisationType::update().

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

◆ updateGradientsAndValues()

void updateGradientsAndValues ( scalarField objectiveSens,
PtrList< scalarField > &  constraintSens,
scalar &  objectiveValue,
scalarField constraintValues 
)
virtual

Compute cumulative objective and constraint gradients.

Definition at line 224 of file optimisationTypeIncompressible.C.

References optimisationType::adjointSolvManagers_, UList< T >::empty(), UPtrList< T >::empty(), forAll, PtrList< T >::set(), List< T >::setSize(), PtrList< T >::setSize(), UPtrList< T >::size(), and Foam::Zero.

Referenced by optimisationType::computeDirection().

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

◆ computeMeritFunction()

scalar computeMeritFunction ( )
virtual

Compute the merit function of the optimisation problem.

Could be different than the objective function in case of constraint optimisation

Definition at line 282 of file optimisationTypeIncompressible.C.

References optimisationType::adjointSolvManagers_, UList< T >::empty(), List< T >::setSize(), optimisationType::updateMethod_, and Foam::Zero.

Here is the call graph for this function:

◆ meritFunctionDirectionalDerivative()

scalar meritFunctionDirectionalDerivative ( )
virtual

Derivative of the merit function.

Definition at line 309 of file optimisationTypeIncompressible.C.

References optimisationType::updateMethod_.

◆ updateOldCorrection()

void updateOldCorrection ( const scalarField oldCorrection)
virtual

Update old correction. Needed for quasi-Newton Methods.

Definition at line 315 of file optimisationTypeIncompressible.C.

References optimisationType::updateMethod_.

Referenced by optimisationType::update().

Here is the caller graph for this function:

◆ write()

void write ( )
virtual

Write useful quantities to files.

Reimplemented in shapeOptimisation.

Definition at line 321 of file optimisationTypeIncompressible.C.

References optimisationType::updateMethod_.

Referenced by optimisationType::update(), and shapeOptimisation::write().

Here is the caller graph for this function:

◆ sourcePtr()

const autoPtr< volScalarField > & sourcePtr ( )

Get source term.

Definition at line 327 of file optimisationTypeIncompressible.C.

References optimisationType::sourcePtr_.

◆ getLineSearch()

autoPtr< lineSearch > & getLineSearch ( )

Get a reference to the line search object.

Definition at line 333 of file optimisationTypeIncompressible.C.

References optimisationType::lineSearch_.

Member Data Documentation

◆ mesh_

◆ dict_

const dictionary dict_
protected

◆ adjointSolvManagers_

◆ updateMethod_

◆ sourcePtr_

autoPtr<volScalarField> sourcePtr_
protected

Definition at line 68 of file optimisationTypeIncompressible.H.

Referenced by optimisationType::sourcePtr().

◆ lineSearch_

autoPtr<lineSearch> lineSearch_
protected

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