updateMethod.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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 Class
29  Foam::updateMethod
30 
31 Description
32  Abstract base class for optimisation methods
33 
34 SourceFiles
35  updateMethod.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef updateMethod_H
40 #define updateMethod_H
41 
42 #include "runTimeSelectionTables.H"
43 #include "IOdictionary.H"
44 #include "fvMesh.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class updateMethod Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 class updateMethod
56 {
57 protected:
58 
59  // Protected data
60 
61  const fvMesh& mesh_;
62 
64 
65  //- Used to output values useful for continuation runs
67 
68  //- Derivatives of the objective functions
70 
71  //- Derivatives of the constraints
73 
74  //- Objective value
75  scalar objectiveValue_;
76 
77  //- Constraint values
79 
80  //- Design variables correction
82 
83  //- Cumulative design variables correction throughout the optimisation
84  //- loop
86 
87  //- Step multiplying the correction
88  scalar eta_;
89 
90  //- Is initially set?
91  bool initialEtaSet_;
92 
93  //- Folder storing the corrections to file
94  // For some optimisation methods with a very high number of
95  // design variables (e.g. topology), it doesn't make much sense
96  // to write all updates in the updateMethodDict. Hence, a
97  // separate file is used to write the corrections, in case they are
98  // needed for post-processing
100 
101  //- Whether to use gSum or sum in the inner products
102  bool globalSum_;
103 
104  // Scalar -- matrix multiplications
105  const scalarField leftMult
106  (
107  const scalarField&,
108  const SquareMatrix<scalar>&
109  );
110 
111  const scalarField rightMult
112  (
113  const SquareMatrix<scalar>&,
114  const scalarField&
115  );
116 
118  (
119  const scalarField&,
120  const scalarField&
121  );
122 
124 
125  //- Compute either global or local sum, based on globalSum flag
126  scalar globalSum(const scalarField& field);
127 
128  //- Compute either global or local sum, based on globalSum flag
129  scalar globalSum(tmp<scalarField>& tfield);
130 
131  //- Return optional dictionary with parameters specific to each method
133 
134 
135 private:
136 
137  // Private Member Functions
138 
139  //- No copy construct
140  updateMethod(const updateMethod&) = delete;
141 
142  //- No copy assignment
143  void operator=(const updateMethod&) = delete;
144 
145 
146 public:
147 
148  //- Runtime type information
149  TypeName("updateMethod");
150 
151 
152  // Declare run-time constructor selection table
153 
155  (
156  autoPtr,
157  updateMethod,
158  dictionary,
159  (
160  const fvMesh& mesh,
161  const dictionary& dict
162  ),
163  (mesh, dict)
164  );
165 
166 
167  // Constructors
168 
169  //- Construct from components
171  (
172  const fvMesh& mesh,
173  const dictionary& dict
174  );
175 
176 
177  // Selectors
178 
179  //- Return a reference to the selected turbulence model
181  (
182  const fvMesh& mesh,
183  const dictionary& dict
184  );
185 
186 
187  //- Destructor
188  virtual ~updateMethod() = default;
189 
190 
191  // Member Functions
192 
193  //- Set objective derivative
194  void setObjectiveDeriv(const scalarField& derivs);
195 
196  //- Set constraints derivative
197  void setConstraintDeriv(const PtrList<scalarField>& derivs);
198 
199  //- Set constraints derivative
200  void setObjectiveValue(const scalar value);
201 
202  //- Set constraints derivative
204 
205  //- Set step for optimisation methods
206  void setStep(const scalar eta);
207 
208  //- Set globalSum variable.
209  // Should be set by the optimisationType owining the updateMethod
210  void setGlobalSum(const bool useGlobalSum);
211 
212  //- Return the correction of the design variables
213  virtual void computeCorrection()=0;
214 
215  //- Return the correction of the design variables
216  //const scalarField& returnCorrection() const;
217 
218  //- Return the correction of the design variables
220 
221  void writeCorrection();
222 
223  //- Compute merit function. Could be different than the objective
224  //- in the presence of constraints
225  virtual scalar computeMeritFunction();
226 
227  //- Directional derivative of the merit function, in the direction of
228  //- the correction. Could be different than the objective directional
229  //- derivative in the presence of constraints
230  virtual scalar meritFunctionDirectionalDerivative();
231 
232  //- Return whether initial eta was set
233  bool& initialEtaSet();
234 
235  //- Update old correction. useful for quasi-newton methods coupled with
236  //- line search
237  virtual void updateOldCorrection(const scalarField& oldCorrection);
238 
239  //- Write useful quantities to files
240  virtual void write();
241 };
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 } // End namespace Foam
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 #endif
251 
252 // ************************************************************************* //
void setConstraintDeriv(const PtrList< scalarField > &derivs)
Set constraints derivative.
Definition: updateMethod.C:291
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition: updateMethod.C:364
dictionary dict
rDeltaTY field()
virtual scalar meritFunctionDirectionalDerivative()
Directional derivative of the merit function, in the direction of the correction. Could be different ...
Definition: updateMethod.C:370
scalarField cValues_
Constraint values.
Definition: updateMethod.H:83
TypeName("updateMethod")
Runtime type information.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
declareRunTimeSelectionTable(autoPtr, updateMethod, dictionary,(const fvMesh &mesh, const dictionary &dict),(mesh, dict))
Abstract base class for optimisation methods.
Definition: updateMethod.H:50
const dictionary dict_
Definition: updateMethod.H:58
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
Definition: updateMethod.C:65
bool & initialEtaSet()
Return whether initial eta was set.
Definition: updateMethod.C:376
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
Definition: updateMethod.C:38
void setConstraintValues(const scalarField &values)
Set constraints derivative.
Definition: updateMethod.C:305
virtual void computeCorrection()=0
Return the correction of the design variables.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
Definition: updateMethod.H:73
scalar eta_
Step multiplying the correction.
Definition: updateMethod.H:99
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
static autoPtr< updateMethod > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition: updateMethod.C:256
scalar objectiveValue_
Objective value.
Definition: updateMethod.H:78
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. useful for quasi-newton methods coupled with line search.
Definition: updateMethod.C:383
dynamicFvMesh & mesh
void setStep(const scalar eta)
Set step for optimisation methods.
Definition: updateMethod.C:311
A class for handling words, derived from Foam::string.
Definition: word.H:63
SquareMatrix< scalar > inv(SquareMatrix< scalar > A)
Definition: updateMethod.C:118
void setObjectiveValue(const scalar value)
Set constraints derivative.
Definition: updateMethod.C:299
void setObjectiveDeriv(const scalarField &derivs)
Set objective derivative.
Definition: updateMethod.C:284
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
Definition: updateMethod.C:168
void setGlobalSum(const bool useGlobalSum)
Set globalSum variable.
Definition: updateMethod.C:317
const fvMesh & mesh_
Definition: updateMethod.H:56
bool initialEtaSet_
Is initially set?
Definition: updateMethod.H:104
bool globalSum_
Whether to use gSum or sum in the inner products.
Definition: updateMethod.H:120
scalarField objectiveDerivatives_
Derivatives of the objective functions.
Definition: updateMethod.H:68
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
Definition: updateMethod.C:92
scalarField correction_
Design variables correction.
Definition: updateMethod.H:88
virtual void write()
Write useful quantities to files.
Definition: updateMethod.C:391
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
scalarField cumulativeCorrection_
Cumulative design variables correction throughout the optimisation loop.
Definition: updateMethod.H:94
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
word correctionFolder_
Folder storing the corrections to file.
Definition: updateMethod.H:115
IOdictionary optMethodIODict_
Used to output values useful for continuation runs.
Definition: updateMethod.H:63
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
scalarField & returnCorrection()
Return the correction of the design variables.
Definition: updateMethod.C:323
virtual ~updateMethod()=default
Destructor.
Namespace for OpenFOAM.
dictionary coeffsDict()
Return optional dictionary with parameters specific to each method.
Definition: updateMethod.C:247