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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 \*---------------------------------------------------------------------------*/
29 
30 #include "adjointSolverManager.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 Foam::adjointSolverManager::adjointSolverManager
43 (
44  fvMesh& mesh,
45  const word& managerType,
46  const dictionary& dict,
47  bool overrideUseSolverName
48 )
49 :
51  (
52  IOobject
53  (
54  "adjointSolverManager" + dict.dictName(),
55  mesh.time().system(),
56  mesh,
57  IOobject::NO_READ,
58  IOobject::NO_WRITE,
59  IOobject::REGISTER
60  )
61  ),
62  mesh_(mesh),
63  dict_(dict),
64  managerName_(dict.dictName()),
65  primalSolverName_(dict.get<word>("primalSolver")),
66  adjointSolvers_(0),
67  objectiveSolverIDs_(0),
68  constraintSolverIDs_(0),
69  operatingPointWeight_
70  (
71  dict.getOrDefault<scalar>("operatingPointWeight", 1)
72  )
73 {
74  dictionary& adjointSolversDict =
75  const_cast<dictionary&>(dict.subDict("adjointSolvers"));
76 
77  const wordList adjSolverNames = adjointSolversDict.toc();
78  adjointSolvers_.setSize(adjSolverNames.size());
79  objectiveSolverIDs_.setSize(adjSolverNames.size());
80  constraintSolverIDs_.setSize(adjSolverNames.size());
81  label nObjectives(0);
82  label nConstraints(0);
83  forAll(adjSolverNames, namei)
84  {
85  dictionary& solverDict =
86  adjointSolversDict.subDict(adjSolverNames[namei]);
87  if (overrideUseSolverName || adjointSolvers_.size() > 1)
88  {
89  solverDict.add<bool>("useSolverNameForFields", true);
90  }
91  adjointSolvers_.set
92  (
93  namei,
95  (
96  mesh_,
97  managerType,
98  solverDict,
100  )
101  );
102 
103  if (adjointSolvers_[namei].isConstraint())
104  {
106  }
107  else
108  {
109  objectiveSolverIDs_[nObjectives++] = namei;
110  }
111  }
114 
115  Info<< "Found " << nConstraints
116  << " adjoint solvers acting as constraints" << endl;
117 
118  // Having more than one non-aggregated objectives per operating point
119  // is needlessly expensive. Issue a warning
120  if (objectiveSolverIDs_.size() > 1)
121  {
123  << "Number of adjoint solvers corresponding to objectives "
124  << "is greater than 1 (" << objectiveSolverIDs_.size() << ")" << nl
125  << "Consider aggregating your objectives to one" << endl;
126  }
127 }
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
132 bool Foam::adjointSolverManager::readDict(const dictionary& dict)
133 {
134  dict_ = dict;
135 
136  const dictionary& adjointSolversDict = dict.subDict("adjointSolvers");
137 
138  // Note: only updating existing solvers
139  for (adjointSolver& solver : adjointSolvers_)
140  {
141  solver.readDict(adjointSolversDict.subDict(solver.name()));
142  }
143 
144  return true;
145 }
146 
149 {
150  return managerName_;
151 }
152 
155 {
156  return primalSolverName_;
157 }
158 
159 
161 {
162  return dict_;
163 }
164 
165 
168 {
169  return adjointSolvers_;
170 }
171 
172 
175 {
176  return adjointSolvers_;
177 }
178 
181 {
182  return operatingPointWeight_;
183 }
184 
187 {
188  return constraintSolverIDs_.size();
189 }
190 
192 Foam::label Foam::adjointSolverManager::nObjectives() const
193 {
194  return objectiveSolverIDs_.size();
195 }
196 
199 {
200  return nConstraints() + nObjectives();
201 }
202 
203 
205 {
206  for (adjointSolver& solver : adjointSolvers_)
207  {
208  // Solve the adjoint equations taking into consideration the weighted
209  // contribution of possibly multiple objectives
210  solver.solve();
211  }
212 }
213 
214 
217 {
218  tmp<scalarField> tsens(new scalarField(0));
219  scalarField& sens = tsens.ref();
220 
221  // Sum sensitivities from all objectives expect the constraints
222  for (const label solveri : objectiveSolverIDs_)
223  {
224  // Sum contributions
225  const scalarField& solverSens =
226  adjointSolvers_[solveri].getObjectiveSensitivities();
227 
228  if (sens.empty())
229  {
230  sens = scalarField(solverSens.size(), Zero);
231  }
232  sens += solverSens;
233  }
234 
235  return tsens;
236 }
237 
238 
241 {
242  PtrList<scalarField> constraintSens(constraintSolverIDs_.size());
243  forAll(constraintSens, cI)
244  {
245  label consI = constraintSolverIDs_[cI];
246  constraintSens.set
247  (
248  cI,
249  new scalarField(adjointSolvers_[consI].getObjectiveSensitivities())
250  );
251  }
252 
253  return constraintSens;
254 }
255 
256 
258 {
259  for (adjointSolver& adjSolver : adjointSolvers_)
260  {
261  adjSolver.computeObjectiveSensitivities();
262  }
263 }
264 
265 
267 {
268  for (adjointSolver& adjSolver : adjointSolvers_)
269  {
270  adjSolver.clearSensitivities();
271  }
272 }
273 
274 
276 {
277  scalar objValue(Zero);
278  for (const label solveri : objectiveSolverIDs_)
279  {
280  objectiveManager& objManager =
281  adjointSolvers_[solveri].getObjectiveManager();
282  objValue += objManager.print();
283  }
284 
285  return objValue;
286 }
287 
288 
290 {
291  tmp<scalarField> tconstraintValues
292  (
293  new scalarField(constraintSolverIDs_.size(), Zero)
294  );
295  scalarField& constraintValues = tconstraintValues.ref();
296  forAll(constraintValues, cI)
297  {
298  objectiveManager& objManager =
299  adjointSolvers_[constraintSolverIDs_[cI]].getObjectiveManager();
300  constraintValues[cI] = objManager.print();
301  }
302 
303  return tconstraintValues;
304 }
305 
306 
308 {
309  if (primalSolverName_ == name)
310  {
311  for (adjointSolver& solver : adjointSolvers_)
312  {
313  solver.updatePrimalBasedQuantities();
314  }
315  }
316 }
317 
318 
319 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
void updatePrimalBasedQuantities(const word &name)
Update fields related to primal solution.
scalar objectiveValue()
Get objective value.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual PtrList< scalarField > constraintSensitivities()
Get constraint sensitivities. One scalarField per constraint.
virtual tmp< scalarField > constraintValues()
Get constraint values.
label nObjectives() const
Number of adjoint solvers corresponding to objectives.
Base class for solution control classes.
Definition: solver.H:45
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
label nAdjointSolvers() const
Total number of adjoint solvers.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
scalar operatingPointWeight() const
Const access to adjoint solvers.
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void clearSensitivities()
Clear sensitivity fields from all adjoint solvers.
virtual tmp< scalarField > aggregateSensitivities()
Aggregate sensitivities from various adjoint solvers.
Base class for adjoint solvers.
Definition: adjointSolver.H:51
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:637
const word & managerName() const
Const access to the manager name.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
wordList toc() const
Return the table of contents.
Definition: dictionary.C:599
PtrList< adjointSolver > adjointSolvers_
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
label nConstraints() const
Number of adjoint solvers corresponding to constraints.
void computeAllSensitivities()
Compute sensitivities for all adjoint solvers (both objective- and constraint-related ones) ...
void setSize(const label n)
Alias for resize()
Definition: List.H:289
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static autoPtr< adjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName)
Return a reference to the selected turbulence model.
Definition: adjointSolver.C:72
const PtrList< adjointSolver > & adjointSolvers() const
Const access to adjoint solvers.
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
const dictionary & dict() const
Const access to the construction dictionary.
virtual bool readDict(const dictionary &dict)
Definition: solver.C:75
virtual void solve()=0
Main control loop.
defineTypeNameAndDebug(combustionModel, 0)
virtual void solveAdjointEquations()
Update objective function-related values and solve adjoint equations.
Class for managing adjoint solvers, which may be more than one per operating point.
#define WarningInFunction
Report a warning using Foam::Warning.
virtual bool readDict(const dictionary &dict)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:65
messageStream Info
Information stream (stdout output on master, null elsewhere)
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition: POSIX.C:1702
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const word & primalSolverName() const
Const access to the primal solver name.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133