optimisationTypeIncompressible.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28
29Class
30 Foam::incompressible::optimisationType
31
32Description
33 Abstract base class for optimisation methods
34
35SourceFiles
36 optimisationType.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef optimisationTypeIncompressible_H
41#define optimisationTypeIncompressible_H
42
44#include "updateMethod.H"
45#include "lineSearch.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52namespace incompressible
53{
54
55/*---------------------------------------------------------------------------*\
56 Class optimisationType Declaration
57\*---------------------------------------------------------------------------*/
60{
61protected:
62
63 // Protected data
71
72 //- Update the design variables given their correction
74
75 //- Compute eta if not set in the first step
76 virtual void computeEta(scalarField& correction) = 0;
77
78
79private:
80
81 // Private Member Functions
82
83 //- No copy construct
84 optimisationType(const optimisationType&) = delete;
85
86 //- No copy assignment
87 void operator=(const optimisationType&) = delete;
88
89
90public:
91
92 //- Runtime type information
93 TypeName("optimisationType");
94
95
96 // Declare run-time constructor selection table
99 (
100 autoPtr,
103 (
104 fvMesh& mesh,
105 const dictionary& dict,
107 ),
109 );
110
111
112
113 // Constructors
114
115 //- Construct from components
117 (
118 fvMesh& mesh,
119 const dictionary& dict,
121 );
122
123 // Selectors
124
125 //- Return a reference to the selected turbulence model
127 (
128 fvMesh& mesh,
129 const dictionary& dict,
131 );
132
133
134 // Destructor
136 virtual ~optimisationType() = default;
137
138 //- Update design variables
139 virtual void update();
140
141 //- Update design variables based on a given correction
142 virtual void update(scalarField& correction);
143
144 //- Store design variables, as the starting point for line search
145 virtual void storeDesignVariables() = 0;
146
147 //- Reset to starting point of line search
148 virtual void resetDesignVariables() = 0;
149
150 //- Compute update direction
152
153 //- Compute cumulative objective and constraint gradients
154 virtual void updateGradientsAndValues
155 (
156 scalarField& objectiveSens,
157 PtrList<scalarField>& constraintSens,
158 scalar& objectiveValue,
159 scalarField& constraintValues
160 );
161
162 //- Compute the merit function of the optimisation problem.
163 // Could be different than the objective function in case of
164 // constraint optimisation
165 virtual scalar computeMeritFunction();
166
167 //- Derivative of the merit function
168 virtual scalar meritFunctionDirectionalDerivative();
169
170 //- Update old correction. Needed for quasi-Newton Methods
171 virtual void updateOldCorrection(const scalarField&);
172
173 //- Write useful quantities to files
174 virtual void write();
175
176 //- Get source term
178
179 //- Get a reference to the line search object
181};
182
183
184// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186} // End namespace incompressible
187} // End namespace Foam
188
189// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190
191#endif
192
193// ************************************************************************* //
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for optimisation methods.
virtual void storeDesignVariables()=0
Store design variables, as the starting point for line search.
virtual void computeEta(scalarField &correction)=0
Compute eta if not set in the first step.
virtual scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
PtrList< adjointSolverManager > & adjointSolvManagers_
const autoPtr< volScalarField > & sourcePtr()
Get source term.
virtual tmp< scalarField > computeDirection()
Compute update direction.
autoPtr< lineSearch > & getLineSearch()
Get a reference to the line search object.
virtual void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
virtual void updateGradientsAndValues(scalarField &objectiveSens, PtrList< scalarField > &constraintSens, scalar &objectiveValue, scalarField &constraintValues)
Compute cumulative objective and constraint gradients.
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
virtual void write()
Write useful quantities to files.
declareRunTimeSelectionTable(autoPtr, optimisationType, dictionary,(fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers),(mesh, dict, adjointSolverManagers))
virtual void resetDesignVariables()=0
Reset to starting point of line search.
virtual void update()
Update design variables.
TypeName("optimisationType")
Runtime type information.
virtual void updateDesignVariables(scalarField &correction)=0
Update the design variables given their correction.
static autoPtr< optimisationType > New(fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers)
Return a reference to the selected turbulence model.
A class for managing temporary objects.
Definition: tmp.H:65
dynamicFvMesh & mesh
Namespace for OpenFOAM.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes)
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8