objective.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-2020 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
28Class
29 Foam::objective
30
31Description
32 Abstract base class for objective functions. No point in making this
33 runTime selectable since its children will have different constructors.
34
35SourceFiles
36 objective.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef objective_H
41#define objective_H
42
43#include "localIOdictionary.H"
44#include "autoPtr.H"
46#include "OFstream.H"
47#include "boundaryFieldsFwd.H"
48#include "solverControl.H"
49#include "objectiveFwd.H"
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace Foam
54{
55
56/*---------------------------------------------------------------------------*\
57 Class objective Declaration
58\*---------------------------------------------------------------------------*/
60class objective
61:
63{
64protected:
65
66 // Protected data
68 const fvMesh& mesh_;
75 bool normalize_;
76
77 //- Objective function value and weight
78 scalar J_;
79
80 //- Average objective value
81 scalar JMean_;
82
83 //- Objective weight
84 scalar weight_;
85
86 //- Normalization factor
88
89 //- Target value, in case the objective is used as a constraint
90 // Should be used in caution and with updateMethods than get affected
91 // by the target value, without requiring a sqr (e.g. SQP, MMA)
93
94 //- Objective integration start and end times (for unsteady flows)
97
98 //- Contribution to field sensitivity derivatives
99 // Topology optimisation or other variants with
100 // as many design variables as the mesh cells
102
103 // Contribution to surface sensitivity derivatives
104 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
105
106 //- Term from material derivative
108
109 //- Term multiplying delta(n dS)/delta b
111
112 //- Term multiplying delta(n)/delta b
114
115 //- Term multiplying delta(x)/delta b at the boundary
117
118 //- Term multiplying delta(x)/delta b at the boundary
119 //- for objectives that directly depend on x, e.g. moment
120 //- Needed in both FI and SI computations
122
123 //- Contribution located in specific parts of a patch.
124 //- Mainly intended for patch boundary edges contributions, e.g.
125 //- NURBS surfaces G1 continuity
127
128 // Contribution to volume-based sensitivities from volume-based
129 // objective functions
130 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
131
132 //- Multiplier of d(Volume)/db
134
135 //- Emerging from volume objectives that include spatial derivatives
137
138 //- Output file variables
140
141 //- File to keep the objective values after the end of the primal solver
143
144 //- File to keep the objective values at each iteration of the primal
145 //- solver
147
148 //- File to keep the average objective values after the end of the
149 //- primal solver
151
152 //- Default width of entries when writing in the objective files
153 unsigned int width_;
154
155
156 // Protected Member Functions
157
158 //- Return objective dictionary
159 const dictionary& dict() const;
160
161 //- Set the output file ptr
162 void setObjectiveFilePtr() const;
163
164 //- Set the output file ptr for the instantaneous value
165 void setInstantValueFilePtr() const;
166
167 //- Set the output file ptr for the mean value
168 void setMeanValueFilePtr() const;
169
170
171private:
172
173 // Private Member Functions
174
175 //- No copy construct
176 objective(const objective&) = delete;
177
178 //- No copy assignment
179 void operator=(const objective&) = delete;
180
181 //- Make objective Function Folder
182 void makeFolder();
183
184
185public:
186
187 //- Runtime type information
188 TypeName("objective");
189
190
191 // Declare run-time constructor selection table
194 (
195 autoPtr,
196 objective,
197 objective,
198 (
199 const fvMesh& mesh,
200 const dictionary& dict,
201 const word& adjointSolverName,
202 const word& primalSolverName
203 ),
204 (mesh, dict, adjointSolverName, primalSolverName)
205 );
206
207
208 // Constructors
209
210 //- Construct from components
212 (
213 const fvMesh& mesh,
214 const dictionary& dict,
215 const word& adjointSolverName,
216 const word& primalSolverName
217 );
218
219
220 // Selectors
221
222 //- Return a reference to the selected turbulence model
224 (
225 const fvMesh& mesh,
226 const dictionary& dict,
227 const word& objectiveType,
228 const word& adjointSolverName,
229 const word& primalSolverName
230 );
231
232
233 //- Destructor
234 virtual ~objective() = default;
235
236
237 // Member Functions
238
239 virtual bool readDict(const dictionary& dict);
240
241 //- Return the instantaneous objective function value
242 virtual scalar J() = 0;
243
244 //- Return the mean objective function value, if it exists,
245 //- otherwise the mean one
246 scalar JCycle() const;
247
248 //- Accumulate contribution for the mean objective value
249 // For steady-state runs
251
252 //- Accumulate contribution for the mean objective value
253 // For unsteady runs
254 void accumulateJMean();
255
256 //- Return the objective function weight
257 scalar weight() const;
258
259 //- Is the objective normalized
260 bool normalize() const;
261
262 //- Normalize all fields allocated by the objective
263 virtual void doNormalization();
264
265 //- Check whether this is an objective integration time
266 bool isWithinIntegrationTime() const;
267
268 //- Increment integration times
269 void incrementIntegrationTimes(const scalar timeSpan);
270
271 //- Contribution to field sensitivities
272 const volScalarField& dJdb();
273
274 //- Contribution to surface sensitivities for a specific patch
275 const fvPatchVectorField& boundarydJdb(const label);
276
277 //- Multiplier of delta(n dS)/delta b
278 const fvPatchVectorField& dSdbMultiplier(const label);
279
280 //- Multiplier of delta(n dS)/delta b
281 const fvPatchVectorField& dndbMultiplier(const label);
282
283 //- Multiplier of delta(x)/delta b
284 const fvPatchVectorField& dxdbMultiplier(const label);
285
286 //- Multiplier of delta(x)/delta b
287 const fvPatchVectorField& dxdbDirectMultiplier(const label);
288
289 //- Multiplier located at patch boundary edges
291 (
292 const label patchI,
293 const label edgeI
294 );
295
296 //- Contribution to surface sensitivities for all patches
298
299 //- Multiplier of delta(n dS)/delta b for all patches
301
302 //- Multiplier of delta(n dS)/delta b for all patches
304
305 //- Multiplier of delta(x)/delta b for all patches
307
308 //- Multiplier of delta(x)/delta b for all patches
310
311 //- Multiplier located at patch boundary edges
313
314 //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
316
317 //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
319
320 //- Update objective function derivatives
321 virtual void update() = 0;
322
323 //- Nullify adjoint contributions
324 virtual void nullify();
325
326 //- Update normalization factors, for objectives in
327 //- which the factor is not known a priori
328 virtual void updateNormalizationFactor();
329
330 //- Update objective function derivative term
331 virtual void update_boundarydJdb()
332 {}
333
334 //- Update d (normal dS) / db multiplier. Surface-based sensitivity term
335 virtual void update_dSdbMultiplier()
336 {}
337
338 //- Update d (normal) / db multiplier. Surface-based sensitivity term
339 virtual void update_dndbMultiplier()
340 {}
341
342 //- Update d (x) / db multiplier. Surface-based sensitivity term
343 virtual void update_dxdbMultiplier()
344 {}
345
346 //- Update d (x) / db multiplier. Surface and volume-based sensitivity
347 //- term
348 virtual void update_dxdbDirectMultiplier()
349 {}
350
351 //- Update boundary edge contributions
353 {}
354
355 //- Update div( dx/db multiplier). Volume-based sensitivity term
356 virtual void update_divDxDbMultiplier()
357 {}
358
359 //- Update grad( dx/db multiplier). Volume-based sensitivity term
360 virtual void update_gradDxDbMultiplier()
361 {}
362
363 //- Write objective function history
364 virtual bool write(const bool valid = true) const;
365
366 //- Write objective function history at each primal solver iteration
367 virtual void writeInstantaneousValue() const;
368
369 //- Append a blank line after the end of optimisation cycle to the
370 //- file holding the instantaneous objective values to facilitate
371 //- plotting
372 virtual void writeInstantaneousSeparator() const;
373
374 //- Write mean objective function history
375 virtual void writeMeanValue() const;
376
377 //- Write averaged objective for continuation
378 virtual bool writeData(Ostream& os) const;
379
380 // Helper write functions
381
382 //- Write any information that needs to go the header of the file
383 // (e.g. targets, directions, etc)
384 virtual void addHeaderInfo() const;
385
386 //- Write headers for additional columns
387 virtual void addHeaderColumns() const;
388
389 //- Write information to additional columns
390 virtual void addColumnValues() const;
391
392 //- Return the objective name
393 inline const word& objectiveName() const;
394
395 // Inline functions for checking whether pointers are set or not
396 inline bool hasdJdb() const;
397 inline bool hasBoundarydJdb() const;
398 inline bool hasdSdbMult() const;
399 inline bool hasdndbMult() const;
400 inline bool hasdxdbMult() const;
401 inline bool hasdxdbDirectMult() const;
402 inline bool hasBoundaryEdgeContribution() const;
403 inline bool hasDivDxDbMult() const;
404 inline bool hasGradDxDbMult() const;
405
406 // Inline functions for checking whether integration times are set
407 inline bool hasIntegrationStartTime() const;
408 inline bool hasIntegrationEndTime() const;
409};
410
411
412// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
413
414} // End namespace Foam
415
416// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
417
418#include "objectiveI.H"
419
420// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
421
422#endif
423
424// ************************************************************************* //
Useful typenames for fields defined only at the boundaries.
Generic GeometricBoundaryField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
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
A class for handling file names.
Definition: fileName.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:62
virtual void writeInstantaneousSeparator() const
Definition: objective.C:721
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:704
scalar JCycle() const
Definition: objective.C:248
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:342
autoPtr< OFstream > instantValueFilePtr_
Definition: objective.H:145
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:359
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:106
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:515
virtual void update_dxdbDirectMultiplier()
Definition: objective.H:347
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:589
void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:302
const fvMesh & mesh_
Definition: objective.H:67
scalar JMean_
Average objective value.
Definition: objective.H:80
bool hasdSdbMult() const
Definition: objectiveI.H:51
bool hasdndbMult() const
Definition: objectiveI.H:57
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:109
const word objectiveName_
Definition: objective.H:71
bool hasBoundaryEdgeContribution() const
Definition: objectiveI.H:75
bool hasGradDxDbMult() const
Definition: objectiveI.H:87
static autoPtr< objective > New(const fvMesh &mesh, const dictionary &dict, const word &objectiveType, const word &adjointSolverName, const word &primalSolverName)
Return a reference to the selected turbulence model.
Definition: objective.C:205
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:80
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:733
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:141
autoPtr< OFstream > meanValueFilePtr_
Definition: objective.H:149
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:355
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition: objective.C:777
scalar weight() const
Return the objective function weight.
Definition: objective.C:324
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
Definition: objective.C:760
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Definition: objective.H:120
virtual void update_boundarydJdb()
Update objective function derivative term.
Definition: objective.H:330
bool hasIntegrationStartTime() const
Definition: objectiveI.H:93
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:609
dictionary dict_
Definition: objective.H:68
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
Definition: objective.H:94
const word & objectiveName() const
Return the objective name.
Definition: objectiveI.H:33
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
Definition: objective.C:505
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:59
autoPtr< volScalarField > dJdbPtr_
Contribution to field sensitivity derivatives.
Definition: objective.H:100
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:86
scalar weight_
Objective weight.
Definition: objective.H:83
bool computeMeanFields_
Definition: objective.H:72
virtual bool readDict(const dictionary &dict)
Definition: objective.C:241
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:568
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:545
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:152
virtual scalar J()=0
Return the instantaneous objective function value.
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:95
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:334
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
Definition: objective.C:771
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:338
bool hasDivDxDbMult() const
Definition: objectiveI.H:81
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:132
virtual void update()=0
Update objective function derivatives.
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:525
scalar J_
Objective function value and weight.
Definition: objective.H:77
virtual ~objective()=default
Destructor.
TypeName("objective")
Runtime type information.
bool hasdxdbDirectMult() const
Definition: objectiveI.H:69
bool hasdxdbMult() const
Definition: objectiveI.H:63
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:68
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:403
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:336
bool normalize() const
Is the objective normalized.
Definition: objective.C:330
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:94
const volScalarField & dJdb()
Contribution to field sensitivities.
Definition: objective.C:419
const word primalSolverName_
Definition: objective.H:70
bool hasBoundarydJdb() const
Definition: objectiveI.H:45
bool hasIntegrationEndTime() const
Definition: objectiveI.H:99
virtual void updateNormalizationFactor()
Definition: objective.C:276
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:138
bool hasdJdb() const
Definition: objectiveI.H:39
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
Definition: objective.H:91
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:112
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:351
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:115
const word adjointSolverName_
Definition: objective.H:69
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:535
virtual void addColumnValues() const
Write information to additional columns.
Definition: objective.C:783
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
Definition: objective.C:555
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:382
autoPtr< vectorField3 > bEdgeContribution_
Definition: objective.H:125
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:135
Base class for solver control classes.
Definition: solverControl.H:52
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
runTime write()
Macros to ease declaration of run-time selection tables.
#define declareRunTimeNewSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection for derived classes.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73