Loading...
Searching...
No Matches
objectiveManager.C
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019-2021 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\*---------------------------------------------------------------------------*/
29
30#include "objectiveManager.H"
31#include "IOmanip.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40defineTypeNameAndDebug(objectiveManager, 0);
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
44objectiveManager::objectiveManager
45(
46 const fvMesh& mesh,
47 const dictionary& dict,
50)
51:
53 (
55 (
56 "objectiveManager" + adjointSolverName,
57 mesh.time().system(),
58 mesh,
62 )
63 ),
64 mesh_(mesh),
65 dict_(dict),
68 objectives_(0),
70{
71 // Construct objectives
72 //~~~~~~~~~~~~~~~~~~~~~
73 Info << "Constructing objective functions " << nl << endl;
74 const word objectiveType = dict.get<word>("type");
75 const dictionary& objectiveNamesDict(dict.subDict("objectiveNames"));
76 wordList objectiveNames(objectiveNamesDict.toc());
77 objectives_.setSize(objectiveNames.size());
78
79 forAll(objectiveNames, objectivei)
80 {
81 const word& objectiveName = objectiveNames[objectivei];
82
83 objectives_.set
84 (
85 objectivei,
87 (
88 mesh_,
89 objectiveNamesDict.subDict(objectiveName),
90 objectiveType,
93 )
94 );
95 }
96
97 if (objectives_.empty())
98 {
99 FatalIOErrorInFunction(objectiveNamesDict)
100 << "No objectives have been set - cannot perform an optimisation"
101 << exit(FatalIOError);
102 }
103
104 if (Pstream::master())
105 {
106 if (objectives_.size() > 1)
107 {
108 const Time& time = mesh_.time();
110 (
111 new OFstream
112 (
113 time.globalPath()/"optimisation"/"objective"
114 /time.timeName()/"weightedObjective"+adjointSolverName_
115 )
116 );
117
118 unsigned int width = IOstream::defaultPrecision() + 5;
120 << setw(4) << "#" << " "
121 << setw(width) << "weightedObjective" << " ";
122 for (objective& objI : objectives_)
123 {
125 << setw(width) << objI.objectiveName() << " ";
126 }
128 << endl;
130 }
131}
132
133
134// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135
137{
138 dict_ = dict;
139 for (objective& obj : objectives_)
140 {
141 obj.readDict
142 (
143 dict.subDict("objectiveNames").subDict(obj.objectiveName())
144 );
145 }
146
147 return true;
148}
149
150
152{
153 // Update normalization factors for all objectives
154 for (objective& obj : objectives_)
155 {
156 if (obj.normalize())
158 obj.updateNormalizationFactor();
159 }
160 }
161}
162
163
165{
166 // Update all fields related to the objective function
168 {
169 obj.update();
170 }
171}
172
173
175{
176 // Update contributions to adjoint if true, otherwise return nulls
177 for (objective& obj : objectives_)
178 {
179 if (obj.isWithinIntegrationTime())
180 {
181 obj.update();
182 }
183 else
185 obj.nullify();
186 }
187 }
188}
189
190
191void objectiveManager::incrementIntegrationTimes(const scalar timeSpan)
192{
193 // Update start and end integration times by adding the timeSpan
194 // of the optimisation cycle
196 {
197 obj.incrementIntegrationTimes(timeSpan);
198 }
199}
200
201
203{
204 scalar objValue(Zero);
205 Info<< "Adjoint solver " << adjointSolverName_ << endl;
206 for (objective& obj : objectives_)
207 {
208 // This function is used to obtain the value used to figure out if
209 // line search is converged or not. If the objective is not updated
210 // in each iteration of the primal solver, the old objective value
211 // might be returned. Force the update of the objective the next
212 // time objective::J() is called.
213 obj.setComputed(false);
214 const scalar cost = obj.JCycle(negate);
215 objValue += cost;
216
217 Info<< obj.objectiveName() << " : " << cost << endl;
218 }
220 Info<< "Weighted objective : " << objValue << nl << endl;
221
222 return objValue;
223}
224
225
226void objectiveManager::setWrite(const bool shouldWrite)
227{
229 {
230 obj.setWrite(shouldWrite);
231 }
232}
233
234
236{
237 for (const objective& obj : objectives_)
238 {
239 // Write objective function to file
240 if (obj.shouldWrite())
241 {
242 obj.write();
243 obj.writeMeanValue();
245 }
246
247 return true;
248}
249
250
252(
253 const scalar weightedObjective,
254 const bool valid
255)
256{
258 {
259 auto& os = weightedObjectiveFile_();
260
261 unsigned int width = IOstream::defaultPrecision() + 5;
262
263 os
264 << setw(4) << mesh_.time().timeName() << " "
265 << setw(width) << weightedObjective << " ";
266
267 for (objective& objI : objectives_)
268 {
269 os << setw(width) << objI.JCycle() << " ";
270 }
272 }
273
274 return writeObjectives();
275}
276
277
279{
281 update();
282 scalar weightedObjective = print();
283 writeObjectives(weightedObjective);
284}
285
290}
291
296}
297
302}
303
306{
307 return primalSolverName_;
308}
309
310
312{
313 for (const objective& obj : objectives_)
314 {
315 if (!obj.hasIntegrationStartTime() || !obj.hasIntegrationEndTime())
316 {
318 << "Objective function " << obj.objectiveName()
319 << " does not have a defined integration start or end time "
320 << exit(FatalError);
321 }
322 }
323}
324
325
327{
329 {
330 obj.addSource(matrix);
331 }
332}
333
334
336{
337 for (objective& obj : objectives_)
338 {
339 obj.addSource(matrix);
340 }
341}
342
343
344// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
345
346} // End namespace Foam
347
348// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
fileName globalPath() const
The global path for the case = rootPath/globalCaseName.
Definition TimePathsI.H:108
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
wordList toc() const
Return the table of contents.
Definition dictionary.C:587
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
Class for managing objective functions.
PtrList< objective > & getObjectiveFunctions()
Return reference to objective functions.
void updateAndWrite()
Call all functions required prior to the solution of the adjoint equations.
const word & adjointSolverName() const
Return name of the adjointSolver.
virtual void addSource(fvVectorMatrix &matrix)
Add contribution to adjoint momentum PDEs.
virtual bool readDict(const dictionary &dict)
void setWrite(const bool shouldWrite)
Should the objectives be written to file upon calling write()?
virtual bool writeObjectives()
Write objective function history.
scalar print(bool negate=false)
Print to screen.
void checkIntegrationTimes() const
Check integration times for unsteady runs.
PtrList< objective > objectives_
void updateOrNullify()
Update contributions to adjoint if true, otherwise return nulls.
void incrementIntegrationTimes(const scalar timeSpan)
Increment integration times by the optimisation cycle time-span.
void update()
Update objective function related values.
void updateNormalizationFactor()
Update objective function related values.
autoPtr< OFstream > weightedObjectiveFile_
const word & primalSolverName() const
Return name of the primalSolver.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
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:221
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition POSIX.C:1704
fvMatrix< scalar > fvScalarMatrix
messageStream Info
Information stream (stdout output on master, null elsewhere).
Omanip< int > setw(const int i)
Definition IOmanip.H:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void negate(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299