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 -------------------------------------------------------------------------------
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::objective
30 
31 Description
32  Abstract base class for objective functions. No point in making this
33  runTime selectable since its children will have different constructors.
34 
35 SourceFiles
36  objective.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef objective_H
41 #define objective_H
42 
43 #include "localIOdictionary.H"
44 #include "autoPtr.H"
45 #include "runTimeSelectionTables.H"
46 #include "OFstream.H"
47 #include "boundaryFieldsFwd.H"
48 #include "solverControl.H"
49 #include "objectiveFwd.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class objective Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 class objective
61 :
62  public localIOdictionary
63 {
64 protected:
65 
66  // Protected data
67 
68  const fvMesh& mesh_;
71  const word primalSolverName_;
72  const word objectiveName_;
73  bool computeMeanFields_;
74  bool nullified_;
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;
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 
171 private:
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 
185 public:
186 
187  //- Runtime type information
188  TypeName("objective");
189 
191  // Declare run-time constructor selection table
192 
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
211  objective
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
223  static autoPtr<objective> New
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
352  virtual void update_boundaryEdgeContribution()
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 // ************************************************************************* //
bool hasGradDxDbMult() const
Definition: objectiveI.H:80
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:61
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:454
virtual void updateNormalizationFactor()
Update normalization factors, for objectives in which the factor is not known a priori.
Definition: objective.C:269
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:135
autoPtr< scalar > target_
Target value, in case the objective is used as a constraint.
Definition: objective.H:98
A class for handling file names.
Definition: fileName.H:71
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:396
const word adjointSolverName_
Definition: objective.H:65
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:87
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:726
dictionary dict_
Definition: objective.H:64
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
virtual void update()=0
Update objective function derivatives.
bool hasdxdbMult() const
Definition: objectiveI.H:56
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:173
autoPtr< scalar > normFactor_
Normalization factor.
Definition: objective.H:90
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:466
bool hasdndbMult() const
Definition: objectiveI.H:50
bool hasdJdb() const
Definition: objectiveI.H:32
virtual void doNormalization()
Normalize all fields allocated by the objective.
Definition: objective.C:329
virtual scalar J()=0
Return the instantaneous objective function value.
bool hasDivDxDbMult() const
Definition: objectiveI.H:74
void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:295
autoPtr< vectorField3 > bEdgeContribution_
Contribution located in specific parts of a patch. Mainly intended for patch boundary edges contribut...
Definition: objective.H:149
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:120
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:518
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:582
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:538
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:528
const word objectiveName_
Definition: objective.H:67
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:491
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
bool hasBoundarydJdb() const
Definition: objectiveI.H:38
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:485
virtual bool write(const bool valid=true) const
Write objective function history.
Definition: objective.C:653
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:198
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:104
dynamicFvMesh & mesh
bool computeMeanFields_
Definition: objective.H:68
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:73
virtual void addHeaderColumns() const
Write headers for additional columns.
Definition: objective.C:770
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:460
virtual ~objective()=default
Destructor.
scalar JCycle() const
Return the mean objective function value, if it exists, otherwise the mean one.
Definition: objective.C:241
Base class for solver control classes.
Definition: solverControl.H:45
virtual void update_boundarydJdb()
Update objective function derivative term.
Definition: objective.H:448
const fvMesh & mesh_
Definition: objective.H:63
scalar JMean_
Average objective value.
Definition: objective.H:80
scalar J_
Objective function value and weight.
Definition: objective.H:75
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:508
TypeName("objective")
Runtime type information.
scalar weight() const
Return the objective function weight.
Definition: objective.C:317
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:130
virtual bool writeData(Ostream &os) const
Write averaged objective for continuation.
Definition: objective.C:753
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:697
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
Definition: objective.C:498
autoPtr< volScalarField > dJdbPtr_
Contribution to field sensitivity derivatives.
Definition: objective.H:112
OBJstream os(runTime.globalPath()/outputName)
bool hasdxdbDirectMult() const
Definition: objectiveI.H:62
bool hasIntegrationEndTime() const
Definition: objectiveI.H:92
bool hasIntegrationStartTime() const
Definition: objectiveI.H:86
const word primalSolverName_
Definition: objective.H:66
Useful typenames for fields defined only at the boundaries.
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:375
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:158
virtual void addColumnValues() const
Write information to additional columns.
Definition: objective.C:776
autoPtr< OFstream > meanValueFilePtr_
File to keep the average objective values after the end of the primal solver.
Definition: objective.H:185
const volScalarField & dJdb()
Contribution to field sensitivities.
Definition: objective.C:412
unsigned int width_
Default width of entries when writing in the objective files.
Definition: objective.H:190
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:561
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:479
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:602
virtual void writeInstantaneousSeparator() const
Append a blank line after the end of optimisation cycle to the file holding the instantaneous objecti...
Definition: objective.C:714
bool normalize() const
Is the objective normalized.
Definition: objective.C:323
scalar weight_
Objective weight.
Definition: objective.H:85
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
Definition: objective.C:548
virtual void addHeaderInfo() const
Write any information that needs to go the header of the file.
Definition: objective.C:764
virtual bool readDict(const dictionary &dict)
Definition: objective.C:234
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:163
autoPtr< scalar > integrationStartTimePtr_
Objective integration start and end times (for unsteady flows)
Definition: objective.H:103
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
autoPtr< OFstream > instantValueFilePtr_
File to keep the objective values at each iteration of the primal solver.
Definition: objective.H:179
Macros to ease declaration of run-time selection tables.
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Term multiplying delta(x)/delta b at the boundary for objectives that directly depend on x...
Definition: objective.H:142
virtual void update_dxdbDirectMultiplier()
Update d (x) / db multiplier. Surface and volume-based sensitivity term.
Definition: objective.H:473
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:125
bool hasdSdbMult() const
Definition: objectiveI.H:44
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:168
bool hasBoundaryEdgeContribution() const
Definition: objectiveI.H:68
const word & objectiveName() const
Return the objective name.
Definition: objectiveI.H:26
Namespace for OpenFOAM.