adjointSolver.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 -------------------------------------------------------------------------------
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 "adjointSolver.H"
31 #include "adjointSensitivity.H"
32 #include "designVariables.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
40 }
41 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
48  {
50  (
52  );
53  }
54 }
55 
56 
58 {
59  // Re-read optimisationDict here to cover multi-region cases
60  return
61  IOdictionary
62  (
63  IOobject
64  (
65  "optimisationDict",
66  mesh_.time().globalPath()/"system",
67  mesh_,
70  false
71  )
72  ).subDict("optimisation").subDict("designVariables");
73 }
74 
75 
76 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
77 
78 Foam::adjointSolver::adjointSolver
79 (
80  fvMesh& mesh,
81  const word& managerType,
82  const dictionary& dict,
83  const word& primalSolverName,
84  const word& solverName
85 )
86 :
87  solver(mesh, managerType, dict, solverName),
88  primalSolverName_(primalSolverName),
89  objectiveManager_
90  (
91  mesh,
92  dict.subDict("objectives"),
93  solverName_,
94  primalSolverName
95  ),
96  sensitivities_(nullptr),
97  computeSensitivities_
98  (
99  dict.getOrDefault<bool>("computeSensitivities", true)
100  ),
101  isConstraint_(dict.getOrDefault<bool>("isConstraint", false)),
102  isDoubleSidedConstraint_
103  (dict.getOrDefault<bool>("isDoubleSidedConstraint", false)),
104  adjointSensitivity_(nullptr)
105 {
106  // Force solver to not be a (single-sided) contraint if flagged as
107  // double-sided
109  {
110  isConstraint_ = false;
111  }
112  // Update objective-related quantities to get correct derivatives
113  // in case of continuation
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
119 
121 (
122  fvMesh& mesh,
123  const word& managerType,
124  const dictionary& dict,
125  const word& primalSolverName,
126  const word& solverName
127 )
128 {
129  const word solverType(dict.get<word>("type"));
130 
131  auto* ctorPtr = adjointSolverConstructorTable(solverType);
132 
133  if (!ctorPtr)
134  {
136  (
137  dict,
138  "adjointSolver",
139  solverType,
140  *adjointSolverConstructorTablePtr_
141  ) << exit(FatalIOError);
142  }
143 
144  return autoPtr<adjointSolver>
145  (
146  ctorPtr(mesh, managerType, dict, primalSolverName, solverName)
147  );
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 {
155  if (solver::readDict(dict))
156  {
157  computeSensitivities_ =
158  dict.getOrDefault<bool>("computeSensitivities", true);
159 
160  objectiveManager_.readDict(dict.subDict("objectives"));
161 
162  if (adjointSensitivity_)
163  {
164  adjointSensitivity_().readDict(designVarsDict());
165  }
166 
167  return true;
168  }
169 
170  return false;
171 }
172 
175 {
176  return false;
177 }
178 
179 
181 {
183  return dimless;
184 }
185 
186 
188 {
190  return dimless;
191 }
192 
195 {
196  return nullptr;
197 }
198 
199 
201 {
202  return nullptr;
203 }
204 
205 
207 (
208  autoPtr<designVariables>& designVars
209 )
210 {
211  if (computeSensitivities_)
212  {
213  preCalculateSensitivities();
214  const scalarField& sens =
215  adjointSensitivity_->calculateSensitivities(designVars);
216  if (!sensitivities_)
217  {
218  sensitivities_.reset(new scalarField(sens.size(), Zero));
219  }
220  sensitivities_.ref() = sens;
221  }
222  else
223  {
224  sensitivities_.reset(new scalarField());
225  }
226 }
227 
228 
230 (
231  autoPtr<designVariables>& designVars
232 )
233 {
234  if (!sensitivities_)
235  {
236  // Read sensitivities from file in case of continuation
237  // Done here and not in allocateSensitivities since the size of the
238  // design variables and, hence, the sensitivities is not known there
239  if (dictionary::found("sensitivities"))
240  {
241  sensitivities_ =
243  ("sensitivities", *this, designVars().size());
244  }
245  else
246  {
247  computeObjectiveSensitivities(designVars);
248  }
249  }
250 
251  return sensitivities_();
252 }
253 
254 
256 {
257  if (computeSensitivities_)
258  {
259  adjointSensitivity_->clearSensitivities();
260  sensitivities_.clear();
261  }
262 }
263 
266 {
267  // Does nothing in base
268 }
269 
270 
272 {
273  if (sensitivities_.valid())
274  {
275  sensitivities_().writeEntry("sensitivities", os);
276  }
277  return true;
278 }
279 
280 
281 // ************************************************************************* //
bool computeSensitivities_
Are sensitivities computed.
Definition: adjointSolver.H:94
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Base solver class.
Definition: solver.H:45
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void allocateSensitivities()
Allocate the sensitivity derivatives.
Definition: adjointSolver.C:38
dictionary designVarsDict() const
Return the dictionary corresponding to the design variables.
Definition: adjointSolver.C:50
Base class for adjoint solvers.
Definition: adjointSolver.H:53
autoPtr< adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
virtual void clearSensitivities()
Clears the sensitivity field known by the adjoint solver.
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
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.
virtual bool writeData(Ostream &os) const
Write the sensitivity derivatives.
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool isConstraint_
Is the adjoint solver used to tackle a constraint.
Definition: adjointSolver.H:99
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual bool readDict(const dictionary &dict)
Definition: solver.C:70
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
void update()
Update objective function related values.
virtual void computeObjectiveSensitivities(autoPtr< designVariables > &designVars)
Compute sensitivities of the underlaying objectives.
bool isDoubleSidedConstraint_
Is the adjoint solver used to tackle a double-sided constraint.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual bool readDict(const dictionary &dict)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE. ...
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
objectiveManager objectiveManager_
Object to manage objective functions.
Definition: adjointSolver.H:84
virtual const scalarField & getObjectiveSensitivities(autoPtr< designVariables > &designVars)
Grab a reference to the computed sensitivities.
virtual bool includeDistance() const
Does the adjoint to an equation computing distances need to taken into consideration.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127