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-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::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 "autoPtr.H"
44 #include "runTimeSelectionTables.H"
45 #include "OFstream.H"
46 #include "boundaryFieldsFwd.H"
47 #include "solverControl.H"
48 #include "objectiveFwd.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class objective Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class objective
60 {
61 protected:
62 
63  // Protected data
64 
65  const fvMesh& mesh_;
71  bool nullified_;
72 
73  // Objective function value and weight
74  scalar J_;
75  scalar JMean_; //average value
76  scalar weight_;
77 
78  // Objective integration start and end times (for unsteady flows)
81 
82  // Contribution to field sensitivity derivatives
83  // Topology optimisation
84  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
85 
87 
88  // Contribution to surface sensitivity derivatives
89  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
90 
91  //- Term from material derivative
93 
94  //- Term multiplying delta(n dS)/delta b
96 
97  //- Term multiplying delta(n)/delta b
99 
100  //- Term multiplying delta(x)/delta b at the boundary
102 
103  //- Term multiplying delta(x)/delta b at the boundary
104  //- for objectives that directly depend on x, e.g. moment
105  //- Needed in both FI and SI computations
107 
108  //- Contribution located in specific parts of a patch.
109  //- Mainly intended for patch boundary edges contributions, e.g.
110  //- NURBS surfaces G1 continuity
112 
113  //- For use with discrete-like sensitivities
115 
116  // Contribution to volume-based sensitivities from volume-based
117  // objective functions
118  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
119 
120  //- Multiplier of d(Volume)/db
122 
123  //- Emerging from volume objectives that include spatial derivatives
125 
126  //- Output file variables
128 
129  //- File to keep the objective values after the end of the primal solver
131 
132  //- File to keep the objective values at each iteration of the primal
133  //- solver
135 
136  //- File to keep the average objective values after the end of the
137  //- primal solver
139 
140 
141  // Protected Member Functions
142 
143  //- Return objective dictionary
144  const dictionary& dict() const;
145 
146  //- Set the output file ptr
147  void setObjectiveFilePtr() const;
148 
149  //- Set the output file ptr for the instantaneous value
150  void setInstantValueFilePtr() const;
151 
152  //- Set the output file ptr for the mean value
153  void setMeanValueFilePtr() const;
154 
155 
156 private:
157 
158  // Private Member Functions
159 
160  //- No copy construct
161  objective(const objective&) = delete;
162 
163  //- No copy assignment
164  void operator=(const objective&) = delete;
165 
166  //- Make objective Function Folder
167  void makeFolder();
168 
169 
170 public:
171 
172  //- Runtime type information
173  TypeName("objective");
174 
175 
176  // Declare run-time constructor selection table
177 
179  (
180  autoPtr,
181  objective,
182  objective,
183  (
184  const fvMesh& mesh,
185  const dictionary& dict,
186  const word& adjointSolverName,
187  const word& primalSolverName
188  ),
189  (mesh, dict, adjointSolverName, primalSolverName)
190  );
191 
192 
193  // Constructors
194 
195  //- Construct from components
196  objective
197  (
198  const fvMesh& mesh,
199  const dictionary& dict,
200  const word& adjointSolverName,
201  const word& primalSolverName
202  );
203 
204 
205  // Selectors
206 
207  //- Return a reference to the selected turbulence model
208  static autoPtr<objective> New
209  (
210  const fvMesh& mesh,
211  const dictionary& dict,
212  const word& objectiveType,
213  const word& adjointSolverName,
214  const word& primalSolverName
215  );
216 
217 
218  //- Destructor
219  virtual ~objective() = default;
220 
221 
222  // Member Functions
223 
224  virtual bool readDict(const dictionary& dict);
225 
226  //- Return the instantaneous objective function value
227  virtual scalar J() = 0;
228 
229  //- Return the mean objective function value, if it exists,
230  //- otherwise the mean one
231  scalar JCycle() const;
232 
233  //- Accumulate contribution for the mean objective value
234  // For steady-state runs
236 
237  //- Accumulate contribution for the mean objective value
238  // For unsteady runs
239  void accumulateJMean();
240 
241  //- Return the objective function weight
242  scalar weight() const;
243 
244  //- Check whether this is an objective integration time
245  bool isWithinIntegrationTime() const;
246 
247  //- Increment integration times
248  void incrementIntegrationTimes(const scalar timeSpan);
249 
250  //- Contribution to field sensitivities
251  const volScalarField& dJdb();
252 
253  //- Contribution to surface sensitivities for a specific patch
254  const fvPatchVectorField& boundarydJdb(const label);
255 
256  //- Multiplier of delta(n dS)/delta b
258 
259  //- Multiplier of delta(n dS)/delta b
261 
262  //- Multiplier of delta(x)/delta b
264 
265  //- Multiplier of delta(x)/delta b
267 
268  //- Multiplier located at patch boundary edges
270  (
271  const label patchI,
272  const label edgeI
273  );
274 
275  //- Objective partial deriv wrt stress tensor
277 
278  //- Contribution to surface sensitivities for all patches
280 
281  //- Multiplier of delta(n dS)/delta b for all patches
283 
284  //- Multiplier of delta(n dS)/delta b for all patches
286 
287  //- Multiplier of delta(x)/delta b for all patches
289 
290  //- Multiplier of delta(x)/delta b for all patches
292 
293  //- Multiplier located at patch boundary edges
295 
296  //- Objective partial deriv wrt stress tensor
298 
299  //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
301 
302  //- Multiplier of grad( delta(x)/delta b) for volume-based sensitivities
304 
305  //- Update objective function derivatives
306  virtual void update() = 0;
307 
308  //- Nullify adjoint contributions
309  virtual void nullify();
310 
311  //- Update normalization factors, for objectives in
312  //- which the factor is not known a priori
313  virtual void updateNormalizationFactor();
314 
315  //- Update objective function derivative term
316  virtual void update_boundarydJdb()
317  {}
318 
319  //- Update d (normal dS) / db multiplier. Surface-based sensitivity term
320  virtual void update_dSdbMultiplier()
321  {}
322 
323  //- Update d (normal) / db multiplier. Surface-based sensitivity term
324  virtual void update_dndbMultiplier()
325  {}
326 
327  //- Update d (x) / db multiplier. Surface-based sensitivity term
328  virtual void update_dxdbMultiplier()
329  {}
330 
331  //- Update d (x) / db multiplier. Surface and volume-based sensitivity
332  //- term
333  virtual void update_dxdbDirectMultiplier()
334  {}
335 
336  //- Update boundary edge contributions
337  virtual void update_boundaryEdgeContribution()
338  {}
339 
340  //- Update dJ/dStress field
341  virtual void update_dJdStressMultiplier()
342  {}
343 
344  //- Update div( dx/db multiplier). Volume-based sensitivity term
345  virtual void update_divDxDbMultiplier()
346  {}
347 
348  //- Update grad( dx/db multiplier). Volume-based sensitivity term
349  virtual void update_gradDxDbMultiplier()
350  {}
351 
352  //- Write objective function history
353  virtual void write() const;
354 
355  //- Write objective function history at each primal solver iteration
356  virtual void writeInstantaneousValue() const;
357 
358  //- Write mean objective function history
359  virtual void writeMeanValue() const;
360 
361  // Inline functions for checking whether pointers are set or not
362  inline const word& objectiveName() const;
363  inline bool hasdJdb() const;
364  inline bool hasBoundarydJdb() const;
365  inline bool hasdSdbMult() const;
366  inline bool hasdndbMult() const;
367  inline bool hasdxdbMult() const;
368  inline bool hasdxdbDirectMult() const;
369  inline bool hasBoundaryEdgeContribution() const;
370  inline bool hasBoundarydJdStress() const;
371  inline bool hasDivDxDbMult() const;
372  inline bool hasGradDxDbMult() const;
373 
374  // Inline functions for checking whether integration times are set
375  inline bool hasIntegrationStartTime() const;
376  inline bool hasIntegrationEndTime() const;
377 };
378 
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 } // End namespace Foam
383 
384 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
385 
386 #include "objectiveI.H"
387 
388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389 
390 #endif
391 
392 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::objective::nullify
virtual void nullify()
Nullify adjoint contributions.
Definition: objective.C:538
Foam::objective::bdxdbDirectMultPtr_
autoPtr< boundaryVectorField > bdxdbDirectMultPtr_
Definition: objective.H:105
Foam::objective::bdJdbPtr_
autoPtr< boundaryVectorField > bdJdbPtr_
Term from material derivative.
Definition: objective.H:91
Foam::objective::boundaryEdgeMultiplier
const vectorField3 & boundaryEdgeMultiplier()
Multiplier located at patch boundary edges.
Definition: objective.C:474
Foam::objective::hasBoundaryEdgeContribution
bool hasBoundaryEdgeContribution() const
Definition: objectiveI.H:75
Foam::objective::integrationStartTimePtr_
autoPtr< scalar > integrationStartTimePtr_
Definition: objective.H:78
Foam::objective::gradDxDbMultPtr_
autoPtr< volTensorField > gradDxDbMultPtr_
Emerging from volume objectives that include spatial derivatives.
Definition: objective.H:123
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
objectiveI.H
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::objective::declareRunTimeNewSelectionTable
declareRunTimeNewSelectionTable(autoPtr, objective, objective,(const fvMesh &mesh, const dictionary &dict, const word &adjointSolverName, const word &primalSolverName),(mesh, dict, adjointSolverName, primalSolverName))
Foam::objective::bdJdStressPtr_
autoPtr< boundaryTensorField > bdJdStressPtr_
For use with discrete-like sensitivities.
Definition: objective.H:113
Foam::objective::objectiveName_
const word objectiveName_
Definition: objective.H:68
Foam::objective::dndbMultiplier
const boundaryVectorField & dndbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:444
Foam::objective::hasdSdbMult
bool hasdSdbMult() const
Definition: objectiveI.H:51
Foam::objective::dxdbDirectMultiplier
const boundaryVectorField & dxdbDirectMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:464
Foam::objective::New
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:182
Foam::objective::setInstantValueFilePtr
void setInstantValueFilePtr() const
Set the output file ptr for the instantaneous value.
Definition: objective.C:67
Foam::objective::bdxdbMultPtr_
autoPtr< boundaryVectorField > bdxdbMultPtr_
Term multiplying delta(x)/delta b at the boundary.
Definition: objective.H:100
Foam::objective::setObjectiveFilePtr
void setObjectiveFilePtr() const
Set the output file ptr.
Definition: objective.C:58
Foam::objective::update_divDxDbMultiplier
virtual void update_divDxDbMultiplier()
Update div( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:344
Foam::solverControl
Base class for solver control classes.
Definition: solverControl.H:51
Foam::objective::JCycle
scalar JCycle() const
Definition: objective.C:225
Foam::objective::update_boundarydJdb
virtual void update_boundarydJdb()
Update objective function derivative term.
Definition: objective.H:315
objectiveFwd.H
Foam::objective::updateNormalizationFactor
virtual void updateNormalizationFactor()
Definition: objective.C:240
Foam::objective::update_dSdbMultiplier
virtual void update_dSdbMultiplier()
Update d (normal dS) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:319
Foam::objective::update_dxdbMultiplier
virtual void update_dxdbMultiplier()
Update d (x) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:327
Foam::objective::hasdxdbMult
bool hasdxdbMult() const
Definition: objectiveI.H:63
Foam::objective::objFunctionFilePtr_
autoPtr< OFstream > objFunctionFilePtr_
File to keep the objective values after the end of the primal solver.
Definition: objective.H:129
Foam::objective::hasdJdb
bool hasdJdb() const
Definition: objectiveI.H:39
Foam::objective::setMeanValueFilePtr
void setMeanValueFilePtr() const
Set the output file ptr for the mean value.
Definition: objective.C:79
Foam::objective::adjointSolverName_
const word adjointSolverName_
Definition: objective.H:66
OFstream.H
Foam::objective::writeMeanValue
virtual void writeMeanValue() const
Write mean objective function history.
Definition: objective.C:627
Foam::objective::~objective
virtual ~objective()=default
Destructor.
Foam::objective::writeInstantaneousValue
virtual void writeInstantaneousValue() const
Write objective function history at each primal solver iteration.
Definition: objective.C:610
Foam::objective::dJdbPtr_
autoPtr< volScalarField > dJdbPtr_
Definition: objective.H:85
Foam::objective::hasIntegrationStartTime
bool hasIntegrationStartTime() const
Definition: objectiveI.H:99
Foam::objective::hasIntegrationEndTime
bool hasIntegrationEndTime() const
Definition: objectiveI.H:105
Foam::objective::hasDivDxDbMult
bool hasDivDxDbMult() const
Definition: objectiveI.H:81
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::objective::bEdgeContribution_
autoPtr< vectorField3 > bEdgeContribution_
Definition: objective.H:110
Foam::Field< vector >
Foam::objective::accumulateJMean
void accumulateJMean()
Accumulate contribution for the mean objective value.
Definition: objective.C:263
Foam::objective::J_
scalar J_
Definition: objective.H:73
Foam::objective::hasdndbMult
bool hasdndbMult() const
Definition: objectiveI.H:57
Foam::objective::divDxDbMultiplier
const volScalarField & divDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:497
Foam::objective::gradDxDbMultiplier
const volTensorField & gradDxDbMultiplier()
Multiplier of grad( delta(x)/delta b) for volume-based sensitivities.
Definition: objective.C:518
Foam::objective::bdndbMultPtr_
autoPtr< boundaryVectorField > bdndbMultPtr_
Term multiplying delta(n)/delta b.
Definition: objective.H:97
Foam::objective::update_boundaryEdgeContribution
virtual void update_boundaryEdgeContribution()
Update boundary edge contributions.
Definition: objective.H:336
Foam::objective::dxdbMultiplier
const boundaryVectorField & dxdbMultiplier()
Multiplier of delta(x)/delta b for all patches.
Definition: objective.C:454
Foam::objective::primalSolverName_
const word primalSolverName_
Definition: objective.H:67
Foam::objective::weight_
scalar weight_
Definition: objective.H:75
Foam::objective::update_gradDxDbMultiplier
virtual void update_gradDxDbMultiplier()
Update grad( dx/db multiplier). Volume-based sensitivity term.
Definition: objective.H:348
Foam::objective::divDxDbMultPtr_
autoPtr< volScalarField > divDxDbMultPtr_
Multiplier of d(Volume)/db.
Definition: objective.H:120
Foam::objective::meanValueFilePtr_
autoPtr< OFstream > meanValueFilePtr_
Definition: objective.H:137
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::objective::integrationEndTimePtr_
autoPtr< scalar > integrationEndTimePtr_
Definition: objective.H:79
Foam::objective::boundarydJdStress
const boundaryTensorField & boundarydJdStress()
Objective partial deriv wrt stress tensor.
Definition: objective.C:487
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::objective::instantValueFilePtr_
autoPtr< OFstream > instantValueFilePtr_
Definition: objective.H:133
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::objective::hasBoundarydJdb
bool hasBoundarydJdb() const
Definition: objectiveI.H:45
Foam::objective::hasBoundarydJdStress
bool hasBoundarydJdStress() const
Definition: objectiveI.H:93
Foam::objective::bdSdbMultPtr_
autoPtr< boundaryVectorField > bdSdbMultPtr_
Term multiplying delta(n dS)/delta b.
Definition: objective.H:94
Foam::objective::update_dndbMultiplier
virtual void update_dndbMultiplier()
Update d (normal) / db multiplier. Surface-based sensitivity term.
Definition: objective.H:323
Foam::objective::objectiveName
const word & objectiveName() const
Definition: objectiveI.H:33
Foam::autoPtr< scalar >
Foam::objective::nullified_
bool nullified_
Definition: objective.H:70
Foam::objective::update_dxdbDirectMultiplier
virtual void update_dxdbDirectMultiplier()
Definition: objective.H:332
Foam::objective::computeMeanFields_
bool computeMeanFields_
Definition: objective.H:69
Foam::GeometricField::Boundary
Definition: GeometricField.H:114
boundaryFieldsFwd.H
Useful typenames for fields defined only at the boundaries.
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::objective::dSdbMultiplier
const boundaryVectorField & dSdbMultiplier()
Multiplier of delta(n dS)/delta b for all patches.
Definition: objective.C:434
Foam::objective::objFunctionFolder_
fileName objFunctionFolder_
Output file variables.
Definition: objective.H:126
Foam::objective::hasGradDxDbMult
bool hasGradDxDbMult() const
Definition: objectiveI.H:87
Foam::objective::mesh_
const fvMesh & mesh_
Definition: objective.H:64
Foam::objective::boundarydJdb
const boundaryVectorField & boundarydJdb()
Contribution to surface sensitivities for all patches.
Definition: objective.C:424
Foam::objective::dict
const dictionary & dict() const
Return objective dictionary.
Definition: objective.C:93
Foam::objective::dict_
dictionary dict_
Definition: objective.H:65
Foam::objective::JMean_
scalar JMean_
Definition: objective.H:74
Foam::objective::update_dJdStressMultiplier
virtual void update_dJdStressMultiplier()
Update dJ/dStress field.
Definition: objective.H:340
Foam::objective::TypeName
TypeName("objective")
Runtime type information.
Foam::objective::weight
scalar weight() const
Return the objective function weight.
Definition: objective.C:285
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::objective::hasdxdbDirectMult
bool hasdxdbDirectMult() const
Definition: objectiveI.H:69
Foam::objective
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:58
Foam::objective::incrementIntegrationTimes
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times.
Definition: objective.C:312
Foam::objective::readDict
virtual bool readDict(const dictionary &dict)
Definition: objective.C:218
Foam::objective::isWithinIntegrationTime
bool isWithinIntegrationTime() const
Check whether this is an objective integration time.
Definition: objective.C:291
Foam::objective::J
virtual scalar J()=0
Return the instantaneous objective function value.
Foam::objective::dJdb
const volScalarField & dJdb()
Contribution to field sensitivities.
Definition: objective.C:328
solverControl.H
Foam::objective::write
virtual void write() const
Write objective function history.
Definition: objective.C:593
Foam::objective::update
virtual void update()=0
Update objective function derivatives.
autoPtr.H