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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 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 
32 #include "runTimeSelectionTables.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 namespace incompressible
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(optimisationType, 0);
45 defineRunTimeSelectionTable(optimisationType, dictionary);
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 optimisationType::optimisationType
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  // Has constraints but is not a constraint optimisation method
90  auto cstTablePtr
91  (
92  constrainedOptimisationMethod::dictionaryConstructorTablePtr_
93  );
95  << "Found " << nConstraints << " adjoint solvers corresponding to "
96  << "constraints but the optimisation method used "
97  << "(" << updateMethod_().type() << ") "
98  << "is not a constrainedOptimisationMethod. " << nl
99  << "Available constrainedOptimisationMethods are :" << nl
100  << cstTablePtr->sortedToc()
101  << exit(FatalError);
102  }
103  else if
104  (
105  !nConstraints
106  && isA<constrainedOptimisationMethod>(updateMethod_())
107  )
108  {
109  // Does not have constraints but is a constrained optimisation method
111  << "Did not find any adjoint solvers corresponding to "
112  << "constraints but the optimisation method used "
113  << "(" << updateMethod_().type() << ") "
114  << "is a constrainedOptimisationMethod. " << nl << nl
115  << "This can cause some constraintOptimisationMethods to misbehave."
116  << nl << nl
117  << "Either the isConstraint bool is not set in one of the adjoint "
118  << "solvers or you should consider using an updateMethod "
119  << "that is not a constrainedOptimisationMethod"
120  << nl << endl;
121  }
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
126 
128 (
129  fvMesh& mesh,
130  const dictionary& dict,
132 )
133 {
134  const word modelType(dict.subDict("optimisationType").get<word>("type"));
135 
136  Info<< "optimisationType type : " << modelType << endl;
137 
138  auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
139 
140  if (!cstrIter.found())
141  {
143  (
144  dict,
145  "optimisationType",
146  modelType,
147  *dictionaryConstructorTablePtr_
148  ) << exit(FatalIOError);
149  }
150 
152  (
153  cstrIter()(mesh, dict, adjointSolverManagers)
154  );
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
159 
161 {
162  // Sum contributions
163  scalarField objectiveSens;
164  PtrList<scalarField> constraintSens;
165  scalar objectiveValue(Zero);
166  scalarField constraintValues;
167 
168  for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
169  {
170  const scalar opWeight = adjSolvManager.operatingPointWeight();
171 
172  // Allocate objective sens size if necessary
173  tmp<scalarField> tadjointSolverManagerSens =
174  adjSolvManager.aggregateSensitivities();
175 
176  if (objectiveSens.empty())
177  {
178  objectiveSens.setSize(tadjointSolverManagerSens().size(), Zero);
179  }
180 
181  objectiveSens += opWeight*tadjointSolverManagerSens();
182  objectiveValue += opWeight*adjSolvManager.objectiveValue();
183 
184  // Allocate constraint sens size if necessary
185  PtrList<scalarField> adjointSolverManagerConstSens =
186  adjSolvManager.constraintSensitivities();
187 
188  tmp<scalarField> cValues = adjSolvManager.constraintValues();
189 
190  if (constraintSens.empty())
191  {
192  constraintSens.setSize(adjointSolverManagerConstSens.size());
193  forAll(constraintSens, cI)
194  {
195  constraintSens.set
196  (
197  cI,
198  new scalarField
199  (
200  adjointSolverManagerConstSens[cI].size(),
201  Zero
202  )
203  );
204  constraintValues.setSize(cValues().size());
205  constraintValues = Zero;
206  }
207  }
208 
209  forAll(constraintSens, cI)
210  {
211  constraintSens[cI] += opWeight*adjointSolverManagerConstSens[cI];
212  }
213  constraintValues += opWeight*cValues();
214  }
215 
216  // Based on the sensitivities, return design variables correction
217  updateMethod_->setObjectiveDeriv(objectiveSens);
218  updateMethod_->setConstraintDeriv(constraintSens);
219  updateMethod_->setObjectiveValue(objectiveValue);
220  updateMethod_->setConstraintValues(constraintValues);
221  tmp<scalarField> tcorrection
222  (
223  new scalarField(objectiveSens.size(), Zero)
224  );
225  scalarField& correction = tcorrection.ref();
226  correction = updateMethod_->returnCorrection();
227 
228  return tcorrection;
229 }
230 
231 
233 {
234  // Compute new objective and constraint values and update the ones
235  // in updateMethod
236  scalar objectiveValue(Zero);
237  scalarField constraintValues;
238 
239  for (adjointSolverManager& adjSolvManager : adjointSolvManagers_)
240  {
241  const scalar opWeight = adjSolvManager.operatingPointWeight();
242 
243  objectiveValue += opWeight*adjSolvManager.objectiveValue();
244  tmp<scalarField> cValues = adjSolvManager.constraintValues();
245 
246  if (constraintValues.empty())
247  {
248  constraintValues.setSize(cValues().size(), Zero);
249  }
250  constraintValues += opWeight*cValues();
251  }
252  updateMethod_->setObjectiveValue(objectiveValue);
253  updateMethod_->setConstraintValues(constraintValues);
254 
255  return updateMethod_->computeMeritFunction();
256 }
257 
258 
260 {
261  return updateMethod_->meritFunctionDirectionalDerivative();
262 }
263 
264 
266 {
267  updateMethod_->updateOldCorrection(oldCorrection);
268 }
269 
270 
272 {
273  updateMethod_->write();
274 }
275 
276 
278 {
279  return sourcePtr_;
280 }
281 
282 
284 {
285  return lineSearch_;
286 }
287 
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace incompressible
292 } // End namespace Foam
293 
294 // ************************************************************************* //
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::incompressible::optimisationType::getLineSearch
autoPtr< lineSearch > & getLineSearch()
Get a reference to the line search object.
Definition: optimisationTypeIncompressible.C:283
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
Foam::incompressible::optimisationType::meritFunctionDirectionalDerivative
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function.
Definition: optimisationTypeIncompressible.C:259
Foam::incompressible::optimisationType::write
virtual void write()
Write useful quantities to files.
Definition: optimisationTypeIncompressible.C:271
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::incompressible::optimisationType::sourcePtr
const autoPtr< volScalarField > & sourcePtr()
Get source term.
Definition: optimisationTypeIncompressible.C:277
Foam::incompressible::optimisationType::sourcePtr_
autoPtr< volScalarField > sourcePtr_
Definition: optimisationTypeIncompressible.H:68
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::incompressible::optimisationType::computeDirection
virtual tmp< scalarField > computeDirection()
Compute update direction.
Definition: optimisationTypeIncompressible.C:160
Foam::incompressible::defineRunTimeSelectionTable
defineRunTimeSelectionTable(adjointSensitivity, dictionary)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:380
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::lineSearch::New
static autoPtr< lineSearch > New(const dictionary &dict, const Time &time)
Return a reference to the selected turbulence model.
Definition: lineSearch.C:96
Foam::updateMethod::New
static autoPtr< updateMethod > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected turbulence model.
Definition: updateMethod.C:263
Foam::PtrList::setSize
void setSize(const label newLen)
Same as resize()
Definition: PtrListI.H:108
adjointSolverManagers
PtrList< adjointSolverManager > & adjointSolverManagers
Definition: createFields.H:8
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
Foam::adjointSolverManager
Class for managing adjoint solvers, which may be more than one per operating point.
Definition: adjointSolverManager.H:54
Foam::incompressible::optimisationType::adjointSolvManagers_
PtrList< adjointSolverManager > & adjointSolvManagers_
Definition: optimisationTypeIncompressible.H:66
dict
dictionary dict
Definition: searchingEngine.H:14
constrainedOptimisationMethod.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (if set) or nullptr.
Definition: PtrListI.H:143
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::incompressible::optimisationType::updateMethod_
autoPtr< updateMethod > updateMethod_
Definition: optimisationTypeIncompressible.H:67
Foam::nl
constexpr char nl
Definition: Ostream.H:372
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
optimisationTypeIncompressible.H
Foam::incompressible::optimisationType::New
static autoPtr< optimisationType > New(fvMesh &mesh, const dictionary &dict, PtrList< adjointSolverManager > &adjointSolverManagers)
Return a reference to the selected turbulence model.
Definition: optimisationTypeIncompressible.C:128
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::incompressible::optimisationType::updateOldCorrection
virtual void updateOldCorrection(const scalarField &)
Update old correction. Needed for quasi-Newton Methods.
Definition: optimisationTypeIncompressible.C:265
Foam::incompressible::optimisationType::lineSearch_
autoPtr< lineSearch > lineSearch_
Definition: optimisationTypeIncompressible.H:69
Foam::incompressible::optimisationType::computeMeritFunction
virtual scalar computeMeritFunction()
Compute the merit function of the optimisation problem.
Definition: optimisationTypeIncompressible.C:232