optimisationTypeIncompressible.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-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 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
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39namespace incompressible
40{
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
50(
51 fvMesh& mesh,
52 const dictionary& dict,
54)
55:
56 mesh_(mesh),
57 dict_(dict),
58 adjointSolvManagers_(adjointSolverManagers),
59 updateMethod_
60 (
61 updateMethod::New(mesh_, dict_.subDict("updateMethod"))
62 ),
63 sourcePtr_(nullptr),
64 lineSearch_
65 (
67 (
68 dict_.subDict("updateMethod").subOrEmptyDict("lineSearch"),
69 mesh.time()
70 )
71 )
72{
73 // Figure out number of adjoint solvers corresponding to constraints.
74 // Looks in all operating poitns
75 label nConstraints(0);
76 for (const adjointSolverManager& adjManagerI : adjointSolvManagers_)
77 {
78 nConstraints += adjManagerI.nConstraints();
79 }
80
81 // Sanity checks for combinations of number of constraints and
82 // optimisation methods
83 if
84 (
85 nConstraints
86 && !isA<constrainedOptimisationMethod>(updateMethod_())
87 )
88 {
89 const auto& cnstrTable =
90 *(constrainedOptimisationMethod::dictionaryConstructorTablePtr_);
91
92 // Has constraints but is not a constraint optimisation method
94 << "Found " << nConstraints << " adjoint solvers corresponding to "
95 << "constraints but the optimisation method ("
96 << updateMethod_().type()
97 << ") is not a constrainedOptimisationMethod." << nl
98 << "Available constrainedOptimisationMethods:" << nl
99 << cnstrTable.sortedToc()
100 << exit(FatalError);
101 }
102 else if
103 (
104 !nConstraints
105 && isA<constrainedOptimisationMethod>(updateMethod_())
106 )
107 {
108 // Does not have constraints but is a constrained optimisation method
110 << "Did not find any adjoint solvers corresponding to "
111 << "constraints but the optimisation method ("
112 << updateMethod_().type()
113 << ") is a constrainedOptimisationMethod." << nl << nl
114 << "This can cause some constraintOptimisationMethods to misbehave."
115 << nl << nl
116 << "Either the isConstraint bool is not set in one of the adjoint "
117 << "solvers or you should consider using an updateMethod "
118 << "that is not a constrainedOptimisationMethod"
119 << nl << endl;
120 }
121}
122
123
124// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
125
127(
128 fvMesh& mesh,
129 const dictionary& dict,
131)
132{
133 const word modelType(dict.subDict("optimisationType").get<word>("type"));
134
135 Info<< "optimisationType type : " << modelType << endl;
136
137 auto* ctorPtr = dictionaryConstructorTable(modelType);
138
139 if (!ctorPtr)
140 {
142 (
143 dict,
144 "optimisationType",
145 modelType,
146 *dictionaryConstructorTablePtr_
147 ) << exit(FatalIOError);
148 }
149
151 (
153 );
154}
155
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
158
160{
161 // Compute update of the design variables
162 tmp<scalarField> tcorrection(computeDirection());
163 scalarField& correction = tcorrection.ref();
164
165 // Update design variables given the correction
167
168 // If direction has been scaled (say by setting the initial eta), the
169 // old correction has to be updated
171 write();
172}
173
174
176{
177 // Multiply with line search step, if necessary
179 if (lineSearch_)
180 {
181 correction *= lineSearch_->step();
182 }
183
184 // Update the design variables
186}
187
188
190{
191 // Sum contributions for sensitivities and objective/constraint values
192 scalarField objectiveSens;
193 PtrList<scalarField> constraintSens;
194 scalar objectiveValue(Zero);
195 scalarField constraintValues;
196
198 (
199 objectiveSens,
200 constraintSens,
201 objectiveValue,
202 constraintValues
203 );
204
205 // Based on the sensitivities, return design variables correction
206 updateMethod_->setObjectiveDeriv(objectiveSens);
207 updateMethod_->setConstraintDeriv(constraintSens);
208 updateMethod_->setObjectiveValue(objectiveValue);
209 updateMethod_->setConstraintValues(constraintValues);
210 tmp<scalarField> tcorrection
211 (
212 new scalarField(objectiveSens.size(), Zero)
213 );
214 scalarField& correction = tcorrection.ref();
215 correction = updateMethod_->returnCorrection();
216
217 // Compute eta if needed
219
220 return tcorrection;
221}
222
223
225(
226 scalarField& objectiveSens,
227 PtrList<scalarField>& constraintSens,
228 scalar& objectiveValue,
229 scalarField& constraintValues
230)
231{
232 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
233 {
234 const scalar opWeight = adjSolvManager.operatingPointWeight();
235
236 // Allocate objective sens size if necessary
237 tmp<scalarField> tadjointSolverManagerSens =
238 adjSolvManager.aggregateSensitivities();
239
240 if (objectiveSens.empty())
241 {
242 objectiveSens.setSize(tadjointSolverManagerSens().size(), Zero);
243 }
244
245 objectiveSens += opWeight*tadjointSolverManagerSens();
246 objectiveValue += opWeight*adjSolvManager.objectiveValue();
247
248 // Allocate constraint sens size if necessary
249 PtrList<scalarField> adjointSolverManagerConstSens =
250 adjSolvManager.constraintSensitivities();
251
252 tmp<scalarField> cValues = adjSolvManager.constraintValues();
253
254 if (constraintSens.empty())
255 {
256 constraintSens.setSize(adjointSolverManagerConstSens.size());
257 forAll(constraintSens, cI)
258 {
259 constraintSens.set
260 (
261 cI,
262 new scalarField
263 (
264 adjointSolverManagerConstSens[cI].size(),
265 Zero
266 )
267 );
268 constraintValues.setSize(cValues().size());
269 constraintValues = Zero;
270 }
271 }
272
273 forAll(constraintSens, cI)
274 {
275 constraintSens[cI] += opWeight*adjointSolverManagerConstSens[cI];
276 }
277 constraintValues += opWeight*cValues();
278 }
279}
280
281
283{
284 // Compute new objective and constraint values and update the ones
285 // in updateMethod
286 scalar objectiveValue(Zero);
287 scalarField constraintValues;
288
289 for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
290 {
291 const scalar opWeight = adjSolvManager.operatingPointWeight();
292
293 objectiveValue += opWeight*adjSolvManager.objectiveValue();
294 tmp<scalarField> cValues = adjSolvManager.constraintValues();
295
296 if (constraintValues.empty())
297 {
298 constraintValues.setSize(cValues().size(), Zero);
299 }
300 constraintValues += opWeight*cValues();
301 }
302 updateMethod_->setObjectiveValue(objectiveValue);
303 updateMethod_->setConstraintValues(constraintValues);
304
305 return updateMethod_->computeMeritFunction();
306}
307
308
310{
311 return updateMethod_->meritFunctionDirectionalDerivative();
312}
313
314
316{
317 updateMethod_->updateOldCorrection(oldCorrection);
318}
319
320
322{
323 updateMethod_->write();
324}
325
326
328{
329 return sourcePtr_;
330}
331
332
334{
335 return lineSearch_;
336}
337
338
339// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340
341} // End namespace incompressible
342} // End namespace Foam
343
344// ************************************************************************* //
void setSize(const label n)
Alias for resize()
Definition: List.H:218
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
bool empty() const noexcept
True if the list is empty (ie, size() is zero)
Definition: UPtrListI.H:113
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Class for managing adjoint solvers, which may be more than one per operating point.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for optimisation methods.
virtual void computeEta(scalarField &correction)=0
Compute eta if not set in the first step.
virtual scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
PtrList< adjointSolverManager > & adjointSolvManagers_
const autoPtr< volScalarField > & sourcePtr()
Get source term.
virtual tmp< scalarField > computeDirection()
Compute update direction.
autoPtr< lineSearch > & getLineSearch()
Get a reference to the line search object.
virtual void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
virtual void updateGradientsAndValues(scalarField &objectiveSens, PtrList< scalarField > &constraintSens, scalar &objectiveValue, scalarField &constraintValues)
Compute cumulative objective and constraint gradients.
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
virtual void write()
Write useful quantities to files.
virtual void update()
Update design variables.
virtual void updateDesignVariables(scalarField &correction)=0
Update the design variables given their correction.
Abstract base class for line search methods.
Definition: lineSearch.H:57
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Abstract base class for optimisation methods.
Definition: updateMethod.H:55
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Namespace for OpenFOAM.
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:372
uint8_t direction
Definition: direction.H:56
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Macros to ease declaration of run-time selection tables.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8