updateMethod.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-2022 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 "updateMethod.H"
31 #include "OFstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(updateMethod, 0);
39 }
40 
41 
42 // * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
43 
45 (
46  const scalarField& s,
47  const SquareMatrix<scalar>& m
48 )
49 {
50  if (s.size() != m.n())
51  {
53  << "scalar derivative and HessianInv matrix do not have the "
54  << "same dimension"
55  << abort(FatalError);
56  }
57 
58  scalarField res(s.size(), Zero);
59  forAll(s, i)
60  {
61  forAll(s, j)
62  {
63  res[i] += s[j]*m[j][i];
64  }
65  }
66 
67  return (res);
68 }
69 
70 
72 (
73  const SquareMatrix<scalar>& m,
74  const scalarField& s
75 )
76 {
77  if (s.size() != m.n())
78  {
80  << "scalar derivative and HessianInv matrix do not have the "
81  << "same dimension"
82  << abort(FatalError);
83  }
84 
85  scalarField res(s.size(), Zero);
86  forAll(s, i)
87  {
88  forAll(s, j)
89  {
90  res[i] += m[i][j]*s[j];
91  }
92  }
93 
94  return (res);
95 }
96 
97 
99 (
100  const scalarField& a,
101  const scalarField& b
102 )
103 {
104  if (a.size() != b.size())
105  {
107  << "operands of outerProduct do not have the same dimension"
108  << abort(FatalError);
109  }
110 
111  SquareMatrix<scalar> res(a.size(), Zero);
112  forAll(a, i)
113  {
114  forAll(a, j)
115  {
116  res[i][j] = a[i]*b[j];
117  }
118  }
119 
120  return (res);
121 }
122 
123 
125 Foam::updateMethod::inv(SquareMatrix<scalar> A)
126 {
127  label n(A.n());
128  SquareMatrix<scalar> invA(n, Zero);
129 
130  // LU decomposition of A
131  labelList pivotIndices(n, Zero);
132  LUDecompose(A, pivotIndices);
133  DebugInfo
134  << "LU decomposed A " << A << endl;
135 
136  // Compute inverse of A by successive back-substitutions.
137  for (label j = 0; j < n; j++)
138  {
139  scalarField rhs(n, Zero);
140  rhs[j] = scalar(1);
141  LUBacksubstitute(A, pivotIndices, rhs);
142  // After LUBacksubstitute, rhs contains the j-th column of the inverse
143  for (label i = 0; i < n; i++)
144  {
145  invA[i][j] = rhs[i];
146  }
147  }
148 
149 
150  /*
151  // Alternative using SVD. Slower and less accurate
152  tempscalarRectangularMatrix Atemp(n, n, 0);
153  for (label i = 0; i < n; i++)
154  {
155  for (label j = 0; j < n; j++)
156  {
157  Atemp[i][j] = A[i][j];
158  }
159  }
160  scalarRectangularMatrix invTemp = SVDinv(Atemp);
161  scalarSquareMatrix invA(n, n, 0);
162  for (label i = 0; i < n; i++)
163  {
164  for (label j = 0; j < n; j++)
165  {
166  invA[i][j] = invTemp[i][j];
167  }
168  }
169  */
170 
171  return invA;
172 }
173 
174 
176 {
177  if (globalSum_)
178  {
179  return gSum(field);
180  }
181  return sum(field);
182 }
183 
184 
185 Foam::scalar Foam::updateMethod::globalSum(tmp<scalarField>& tfield)
186 {
187  scalar value = globalSum(tfield());
188  tfield.clear();
189  return value;
190 }
191 
192 
193 Foam::label Foam::updateMethod::globalSum(const label size)
194 {
195  label res(0);
196  if (globalSum_)
197  {
198  res = returnReduce(size, sumOp<label>());
199  }
200  else
201  {
202  res = size;
203  }
204  return res;
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209 
210 Foam::updateMethod::updateMethod
211 (
212  const fvMesh& mesh,
213  const dictionary& dict,
214  autoPtr<designVariables>& designVars,
215  const label nConstraints,
216  const word& type
217 )
218 :
220  (
221  IOobject
222  (
223  "updateMethodDict",
224  mesh.time().timeName(),
225  "uniform",
226  mesh,
227  IOobject::READ_IF_PRESENT,
228  IOobject::AUTO_WRITE
229  ),
230  word::null // avoid type checking
231  ),
232  mesh_(mesh),
233  dict_(dict),
234  designVars_(designVars),
235  nConstraints_(nConstraints),
236  activeDesignVars_(designVars().activeDesignVariables()),
237  objectiveDerivatives_(designVars().getVars().size(), Zero),
238  constraintDerivatives_(0),
239  objectiveValue_(0),
240  objectiveValueOld_(0),
241  cValues_(0),
242  correction_(readOrZeroField("correction", designVars().getVars().size())),
243  cumulativeCorrection_(0),
244  eta_(1),
245  counter_(getOrDefault<label>("counter", Zero)),
246  initialEtaSet_(false),
247  correctionFolder_(mesh_.time().globalPath()/"optimisation"/"correction"),
248  globalSum_(designVars_->globalSum())
249 {
250  // Set initial eta, if present. It might be set either in the
251  // optimisationDict or in the specific dictionary dedicated to the
252  // updateMethod
253  if (dict.readIfPresent("eta", eta_))
254  {
255  initialEtaSet_ = true;
256  }
257  else if (this->readIfPresent("eta", eta_))
258  {
259  initialEtaSet_ = true;
260  }
261 }
262 
263 
265 (
266  const word& name,
267  const label size
268 )
269 {
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * //
275 
277 (
278  const fvMesh& mesh,
279  const dictionary& dict,
280  autoPtr<designVariables>& designVars,
281  const label nConstraints
282 )
283 {
284  const word modelType(dict.get<word>("method"));
285 
286  Info<< "updateMethod type : " << modelType << endl;
287 
288  auto* ctorPtr = dictionaryConstructorTable(modelType);
289 
290  if (!ctorPtr)
291  {
293  (
294  dict,
295  "updateMethod",
296  modelType,
297  *dictionaryConstructorTablePtr_
298  ) << exit(FatalIOError);
299  }
300 
302  (ctorPtr(mesh, dict, designVars, nConstraints, modelType));
303 }
304 
305 
306 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
309 {
310  return dict_.optionalSubDict(type);
311 }
312 
313 
315 {
316  objectiveDerivatives_ = derivs;
317 }
318 
319 
321 (
322  const PtrList<scalarField>& derivs
323 )
324 {
325  constraintDerivatives_ = derivs;
326 }
327 
329 void Foam::updateMethod::setObjectiveValue(const scalar value)
330 {
331  objectiveValue_ = value;
332 }
333 
335 void Foam::updateMethod::setObjectiveValueOld(const scalar value)
336 {
337  objectiveValueOld_ = value;
338 }
339 
342 {
343  cValues_ = values;
344 }
345 
347 Foam::scalar Foam::updateMethod::getObjectiveValue() const
348 {
349  return objectiveValue_;
350 }
351 
353 Foam::scalar Foam::updateMethod::getObjectiveValueOld() const
354 {
355  return objectiveValueOld_;
356 }
357 
360 {
361  return cValues_;
362 }
363 
365 Foam::label Foam::updateMethod::getCycle() const
366 {
367  return counter_;
368 }
369 
371 void Foam::updateMethod::setStep(const scalar eta)
372 {
373  eta_ = eta;
374 }
375 
377 void Foam::updateMethod::modifyStep(const scalar multiplier)
378 {
379  eta_ *= multiplier;
380 }
381 
383 void Foam::updateMethod::setGlobalSum(const bool useGlobalSum)
384 {
385  globalSum_ = useGlobalSum;
386 }
387 
389 void Foam::updateMethod::setConstaintsNumber(const label nConstraints)
390 {
391  nConstraints_ = nConstraints;
392 }
393 
395 Foam::label Foam::updateMethod::nConstraints() const
396 {
397  return nConstraints_;
398 }
399 
402 {
403  return correction_;
404 }
405 
406 
408 {
409  if (Pstream::master())
410  {
411  // Allocate cumulativeCorrection if necessary
412  if (cumulativeCorrection_.empty())
413  {
414  cumulativeCorrection_.setSize(correction_.size(), Zero);
415  }
416  // Accumulate correction
417  cumulativeCorrection_ += correction_;
418 
419  fileName correctionFile
420  (
421  correctionFolder_/"correction" + mesh_.time().timeName()
422  );
423  fileName cumulativeCorrectionFile
424  (
425  correctionFolder_/"cumulativeCorrection" + mesh_.time().timeName()
426  );
427 
428  OFstream corFile(correctionFile);
429  OFstream cumulCorFile(cumulativeCorrectionFile);
430  forAll(correction_, cI)
431  {
432  corFile
433  << cI << " " << correction_[cI] << endl;
434  cumulCorFile
435  << cI << " " << cumulativeCorrection_[cI] << endl;
436  }
437  }
438 }
439 
442 {
443  return objectiveValue_;
444 }
445 
448 {
449  return globalSum(objectiveDerivatives_*correction_);
450 }
451 
452 
454 {
455  return initialEtaSet_;
456 }
457 
458 
460 (
461  const scalarField& oldCorrection
462 )
463 {
464  correction_ = oldCorrection;
465 }
466 
467 
469 {
470  // Insert eta if set
471  if (initialEtaSet_)
472  {
473  os.writeEntry("eta", eta_);
474  }
475 
476  os.writeEntry("counter", counter_);
477  correction_.writeEntry("correction", os);
478 
479  return true;
480 }
481 
482 
484 {
485  // Does nothing in base
486  return true;
487 }
488 
489 
490 // ************************************************************************* //
void setConstraintDeriv(const PtrList< scalarField > &derivs)
Set constraints derivative.
Definition: updateMethod.C:314
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition: updateMethod.C:434
dictionary dict
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
void setConstaintsNumber(const label nConstraints)
Set the number of constraints.
Definition: updateMethod.C:382
rDeltaTY field()
virtual scalar meritFunctionDirectionalDerivative()
Directional derivative of the merit function, in the direction of the correction. Could be different ...
Definition: updateMethod.C:440
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
label nConstraints() const
Get the number of constraints.
Definition: updateMethod.C:388
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
static autoPtr< updateMethod > New(const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints)
Return a reference to the selected turbulence model.
Definition: updateMethod.C:270
Abstract base class for optimisation methods.
Definition: updateMethod.H:50
volSymmTensorField invA(inv(I *UEqn.A()+drag->Dcu()))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
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.
scalar getObjectiveValueOld() const
Get old objective value.
Definition: updateMethod.C:346
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
Definition: updateMethod.C:65
bool & initialEtaSet()
Return whether initial eta was set.
Definition: updateMethod.C:446
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
scalar getObjectiveValue() const
Get objective value.
Definition: updateMethod.C:340
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
Definition: updateMethod.C:38
void setConstraintValues(const scalarField &values)
Set values of constraints.
Definition: updateMethod.C:334
scalar eta_
Step multiplying the correction.
Definition: updateMethod.H:118
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
localIOdictionary is derived from IOdictionary but excludes parallel master reading.
label n() const noexcept
The number of columns.
Definition: Matrix.H:262
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
word timeName
Definition: getTimeIndex.H:3
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:560
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-newton methods coupled with line search.
Definition: updateMethod.C:453
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void setStep(const scalar eta)
Set step for optimisation methods.
Definition: updateMethod.C:364
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
void modifyStep(const scalar multiplier)
Multiply step.
Definition: updateMethod.C:370
virtual bool writeAuxiliaryData()
Write non-continuation data, usually under the optimisation folder.
Definition: updateMethod.C:476
void setObjectiveValueOld(const scalar value)
Set old objective value.
Definition: updateMethod.C:328
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
SquareMatrix< scalar > inv(SquareMatrix< scalar > A)
Definition: updateMethod.C:118
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void setObjectiveValue(const scalar value)
Set objective value.
Definition: updateMethod.C:322
void setObjectiveDeriv(const scalarField &derivs)
Set objective derivative.
Definition: updateMethod.C:307
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
Definition: updateMethod.C:168
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
Definition: updateMethod.C:301
void setGlobalSum(const bool useGlobalSum)
Set globalSum variable.
Definition: updateMethod.C:376
bool initialEtaSet_
Is initially set?
Definition: updateMethod.H:128
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
const scalarField & getConstraintValues() const
Get values of constraints.
Definition: updateMethod.C:352
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
Definition: updateMethod.C:92
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
label getCycle() const
Get optimisation cycle.
Definition: updateMethod.C:358
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Reading is optional [identical to READ_IF_PRESENT].
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
scalarField & returnCorrection()
Return the correction of the design variables.
Definition: updateMethod.C:394
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:645
tmp< scalarField > readOrZeroField(const word &name, const label size)
Helper function to either read a scalarField of certain size from a dictionary, or construct a zero f...
Definition: updateMethod.C:258
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
Definition: updateMethod.C:461
Namespace for OpenFOAM.
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
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