Loading...
Searching...
No Matches
adjointSolverManager.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-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
28\*---------------------------------------------------------------------------*/
29
31#include "primalSolver.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
37 defineTypeNameAndDebug(adjointSolverManager, 0);
38}
39
40
41// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42
43Foam::adjointSolverManager::adjointSolverManager
44(
45 fvMesh& mesh,
46 autoPtr<designVariables>& designVars,
47 const word& managerType,
48 const dictionary& dict,
49 bool overrideUseSolverName
50)
51:
53 (
55 (
56 "adjointSolverManager" + dict.dictName(),
57 mesh.time().system(),
58 mesh,
59 IOobject::NO_READ,
60 IOobject::NO_WRITE,
61 IOobject::REGISTER
62 )
63 ),
64 mesh_(mesh),
65 dict_(dict),
66 managerName_(dict.dictName()),
67 managerType_(managerType),
68 primalSolverName_(dict.get<word>("primalSolver")),
69 adjointSolvers_(0),
70 objectiveSolverIDs_(0),
71 oneSidedConstraintSolverIDs_(0),
72 doubleSidedConstraintSolverIDs_(0),
73 operatingPointWeight_
74 (
75 dict.getOrDefault<scalar>("operatingPointWeight", 1)
76 ),
77 nActiveAdjointSolvers_(0),
78 designVars_(designVars)
79{
80 dictionary& adjointSolversDict =
81 const_cast<dictionary&>(dict.subDict("adjointSolvers"));
82
83 const wordList adjSolverNames = adjointSolversDict.toc();
84 adjointSolvers_.setSize(adjSolverNames.size());
85 objectiveSolverIDs_.setSize(adjSolverNames.size());
88 label nObjectives(0);
89 label nOneSidedConstraints(0);
91 forAll(adjSolverNames, namei)
92 {
93 dictionary& solverDict =
94 adjointSolversDict.subDict(adjSolverNames[namei]);
95 if (overrideUseSolverName)
96 {
97 solverDict.add<bool>("useSolverNameForFields", true);
98 }
100 (
101 namei,
103 (
104 mesh_,
105 managerType,
106 solverDict,
108 adjSolverNames[namei]
109 )
110 );
111 if (adjointSolvers_[namei].active())
112 {
114 }
115 if (adjointSolvers_[namei].isDoubleSidedConstraint())
116 {
118 }
119 else if (adjointSolvers_[namei].isConstraint())
120 {
122 }
123 else
124 {
126 }
127 }
131
132 Info<< "Found " << nOneSidedConstraints
133 << " adjoint solvers acting as single-sided constraints" << endl;
134
135 Info<< "Found " << nDoubleSidedConstraints
136 << " adjoint solvers acting as double-sided constraints" << endl;
137
138 Info<< "Found " << nActiveAdjointSolvers_
139 << " active adjoint solvers" << endl;
140
141 // Having more than one non-aggregated objectives per operating point
142 // is needlessly expensive. Issue a warning
143 if (objectiveSolverIDs_.size() > 1)
144 {
146 << "Number of adjoint solvers corresponding to objectives "
147 << "is greater than 1 (" << objectiveSolverIDs_.size() << ")" << nl
148 << "Consider aggregating your objectives to one\n" << endl;
149 }
150}
151
152
153// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154
156{
157 dict_ = dict;
158
159 const dictionary& adjointSolversDict = dict.subDict("adjointSolvers");
160
161 // Note: only updating existing solvers
162 for (adjointSolver& solver : adjointSolvers_)
163 {
164 solver.readDict(adjointSolversDict.subDict(solver.name()));
165 }
166
167 return true;
168}
169
174}
175
180}
181
182
184{
185 return dict_;
186}
187
188
194
195
205 wordList names(adjointSolvers_.size());
206 forAll(adjointSolvers_, sI)
208 names[sI] = adjointSolvers_[sI].name();
209 }
210 return names;
211}
212
217}
218
221{
223}
224
225
227(
228 const dictionary& dict
229)
230{
231 const dictionary& adjointSolversDict = dict.subDict("adjointSolvers");
232 const wordList adjSolverNames = adjointSolversDict.toc();
233 label n(0);
234 Switch active(true);
235 forAll(adjSolverNames, namei)
236 {
237 active = adjointSolversDict.subDict(adjSolverNames[namei]).
238 getOrDefault<bool>("active", true);
239 if (active)
240 {
241 n++;
242 }
243 }
244 return n;
245}
246
251}
252
257}
258
263}
264
267{
268 return objectiveSolverIDs_.size();
269}
270
273{
275}
276
277
279{
280 // Solve all adjoint equations of this adjointSolverManager
281 for (adjointSolver& solver : adjointSolvers_)
282 {
283 // Update all primal-based quantities needed by the adjoint equations
284 solver.updatePrimalBasedQuantities();
285
286 // Solve the adjoint equations taking into consideration the weighted
287 // contribution of possibly multiple objectives
288 solver.solve();
289
290 // Compute sensitivities and force writing to the adjoint dictionary
291 // if this an output time
292 solver.computeObjectiveSensitivities(designVars_);
293 if (mesh_.time().writeTime())
294 {
295 solver.regIOobject::write(true);
296 }
297 }
298}
299
300
301Foam::tmp<Foam::scalarField>
303{
304 auto tsens = tmp<scalarField>::New();
305 auto& sens = tsens.ref();
306
307 // Sum sensitivities from all objectives expect the constraints
308 for (const label solveri : objectiveSolverIDs_)
309 {
310 // Sum contributions
311 const scalarField& solverSens =
312 adjointSolvers_[solveri].getObjectiveSensitivities(designVars_);
313
314 if (sens.empty())
315 {
316 sens = scalarField(solverSens.size(), Zero);
317 }
318 sens += solverSens;
320
321 return tsens;
322}
323
324
327{
328 PtrList<scalarField> constraintSens(nConstraints());
329 // Only one-sided constraints
330 label cI(0);
331 for (const label consI : oneSidedConstraintSolverIDs_)
332 {
333 constraintSens.set
334 (
335 cI++,
336 new scalarField
337 (adjointSolvers_[consI].getObjectiveSensitivities(designVars_))
338 );
339 }
340
341 // Two-sided constraints. Negated left-most side sensitivities
342 for (const label consI : doubleSidedConstraintSolverIDs_)
343 {
344 scalarField sens
345 (adjointSolvers_[consI].getObjectiveSensitivities(designVars_));
346 constraintSens.set(cI++, new scalarField( sens));
347 constraintSens.set(cI++, new scalarField(- sens));
348 }
349
350 return constraintSens;
351}
352
353
355{
357 {
358 adjSolver.computeObjectiveSensitivities(designVars_);
359 }
360}
361
362
364{
366 {
367 adjSolver.clearSensitivities();
368 }
369}
370
371
373{
374 scalar objValue(Zero);
375 for (const label solveri : objectiveSolverIDs_)
376 {
377 objectiveManager& objManager =
378 adjointSolvers_[solveri].getObjectiveManager();
379 objValue += objManager.print();
380 }
381
382 return objValue;
383}
384
385
387{
388 auto tconstraintValues(tmp<scalarField>::New(nConstraints(), Zero));
389 scalarField& constraintValues = tconstraintValues.ref();
390 label cI(0);
391 // One-sided constraints only
392 for (const label consI : oneSidedConstraintSolverIDs_)
393 {
394 objectiveManager& objManager =
395 adjointSolvers_[consI].getObjectiveManager();
396 constraintValues[cI++] = objManager.print();
397 }
398 // Double-sided constraints
399 // Objective value of the left-most side is negated
400 for (const label consI : doubleSidedConstraintSolverIDs_)
401 {
402 objectiveManager& objManager =
403 adjointSolvers_[consI].getObjectiveManager();
404 constraintValues[cI++] = objManager.print(false);
405 constraintValues[cI++] = objManager.print(true);
406 }
407
408 return tconstraintValues;
409}
410
411
413{
414 if (primalSolverName_ == name)
415 {
416 for (adjointSolver& solver : adjointSolvers_)
418 solver.updatePrimalBasedQuantities();
419 }
420 }
421}
422
423
425{
426 return mesh_.lookupObject<primalSolver>(primalSolverName_).isMaster();
427}
428
429
430// ************************************************************************* //
label n
@ 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
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
void setSize(label n)
Alias for resize().
Definition List.H:536
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Class for managing adjoint solvers, which may be more than one per operating point.
void clearSensitivities()
Clear sensitivity fields from all adjoint solvers.
label nOneSidedConstraints() const
Number of adjoint solvers corresponding to one-sided constraints.
scalar objectiveValue()
Get objective value.
PtrList< scalarField > constraintSensitivities()
Get constraint sensitivities. One scalarField per constraint.
label nAdjointSolvers() const
Total number of adjoint solvers.
tmp< scalarField > constraintValues()
Get constraint values.
bool isMaster() const
Whether the primal solver corresponding to the adjointSolverManager is the master one,...
bool readDict(const dictionary &dict)
tmp< scalarField > aggregateSensitivities()
Aggregate sensitivities from various adjoint solvers.
void computeAllSensitivities()
Compute sensitivities for all adjoint solvers (both objective- and constraint-related ones).
const word & primalSolverName() const
Const access to the primal solver name.
void updatePrimalBasedQuantities(const word &name)
Update fields related to primal solution.
wordList adjointSolversNames() const
Return the names of all adjointSolvers.
PtrList< adjointSolver > adjointSolvers_
label nConstraints() const
Number of constraints.
void solveAdjointEquations()
Update objective function-related values and solve adjoint equations.
const dictionary & dict() const
Const access to the construction dictionary.
label nDoubleSidedConstraints() const
Number of adjoint solvers corresponding to double-sided constraints.
scalar operatingPointWeight() const
Const access to adjoint solvers.
autoPtr< designVariables > & designVars_
const PtrList< adjointSolver > & adjointSolvers() const
Const access to adjoint solvers.
label nActiveAdjointSolvers() const
Return number of active adjoint solvers, either corresponding.
label nObjectives() const
Number of adjoint solvers corresponding to objectives.
const word & managerName() const
Const access to the manager name.
Base class for adjoint solvers.
static autoPtr< adjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected turbulence model.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
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
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
wordList toc() const
Return the table of contents.
Definition dictionary.C:587
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class for managing objective functions.
scalar print(bool negate=false)
Print to screen.
Base class for primal solvers.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition regIOobject.H:71
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
Base solver class.
Definition solver.H:48
virtual bool readDict(const dictionary &dict)
Definition solver.C:70
virtual void solve()=0
Main control loop.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
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
const word dictName("faMeshDefinition")
auto & name
auto & names
#define WarningInFunction
Report a warning using Foam::Warning.
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
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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