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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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 "updateMethod.H"
31 #include "OFstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(updateMethod, 0);
38  defineRunTimeSelectionTable(updateMethod, dictionary);
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 
126 {
127  label n(A.n());
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  scalar value(0);
178  if (globalSum_)
179  {
180  value = gSum(field);
181  }
182  else
183  {
184  value = sum(field);
185  }
186  return value;
187 }
188 
189 
191 {
192  scalar value = globalSum(tfield());
193  tfield.clear();
194  return value;
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199 
200 Foam::updateMethod::updateMethod
201 (
202  const fvMesh& mesh,
203  const dictionary& dict
204 )
205 :
206  mesh_(mesh),
207  dict_(dict),
208  optMethodIODict_
209  (
210  IOobject
211  (
212  "updateMethodDict",
213  mesh_.time().timeName(),
214  "uniform",
215  mesh_,
218  )
219  ),
220  objectiveDerivatives_(0),
221  constraintDerivatives_(0),
222  objectiveValue_(0),
223  cValues_(0),
224  correction_(0),
225  cumulativeCorrection_(0),
226  eta_(1),
227  initialEtaSet_(false),
228  correctionFolder_(mesh_.time().globalPath()/"optimisation"/"correction"),
229  globalSum_
230  (
231  dict.getOrDefault<bool>("globalSum", false)
232  )
233 {
234  // Create folder to store corrections
235  if (Pstream::master())
236  {
237  mkDir(correctionFolder_);
238  }
239 
240  // Set initial eta, if present. It might be set either in the
241  // optimisationDict or in the specific dictionary dedicated to the
242  // updateMethod
243  if (dict.readIfPresent("eta", eta_))
244  {
245  initialEtaSet_ = true;
246  }
247  else if (optMethodIODict_.readIfPresent("eta", eta_))
248  {
249  initialEtaSet_ = true;
250  }
251 }
252 
253 
255 {
256  return dict_.subOrEmptyDict(type());
257 }
258 
259 
260 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * //
261 
263 (
264  const fvMesh& mesh,
265  const dictionary& dict
266 )
267 {
268  const word modelType(dict.get<word>("method"));
269 
270  Info<< "updateMethod type : " << modelType << endl;
271 
272  auto* ctorPtr = dictionaryConstructorTable(modelType);
273 
274  if (!ctorPtr)
275  {
277  (
278  dict,
279  "updateMethod",
280  modelType,
281  *dictionaryConstructorTablePtr_
282  ) << exit(FatalIOError);
283  }
284 
285  return autoPtr<updateMethod>(ctorPtr(mesh, dict));
286 }
287 
288 
289 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * //
290 
292 {
293  objectiveDerivatives_ = derivs;
294 }
295 
296 
298 (
299  const PtrList<scalarField>& derivs
300 )
301 {
302  constraintDerivatives_ = derivs;
303 }
304 
305 
307 {
308  objectiveValue_ = value;
309 }
310 
311 
313 {
314  cValues_ = values;
315 }
316 
317 
318 void Foam::updateMethod::setStep(const scalar eta)
319 {
320  eta_ = eta;
321 }
322 
323 
324 void Foam::updateMethod::setGlobalSum(const bool useGlobalSum)
325 {
326  globalSum_ = useGlobalSum;
327 }
328 
329 
331 {
332  computeCorrection();
333  return correction_;
334 }
335 
336 
338 {
339  if (Pstream::master())
340  {
341  // Allocate cumulativeCorrection if necessary
342  if (cumulativeCorrection_.empty())
343  {
344  cumulativeCorrection_.setSize(correction_.size(), Zero);
345  }
346  // Accumulate correction
347  cumulativeCorrection_ += correction_;
348 
349  fileName correctionFile
350  (
351  correctionFolder_/"correction" + mesh_.time().timeName()
352  );
353  fileName cumulativeCorrectionFile
354  (
355  correctionFolder_/"cumulativeCorrection" + mesh_.time().timeName()
356  );
357 
358  OFstream corFile(correctionFile);
359  OFstream cumulCorFile(cumulativeCorrectionFile);
360  forAll(correction_, cI)
361  {
362  corFile
363  << cI << " " << correction_[cI] << endl;
364  cumulCorFile
365  << cI << " " << cumulativeCorrection_[cI] << endl;
366  }
367  }
368 }
369 
370 
372 {
373  return objectiveValue_;
374 }
375 
376 
378 {
379  return globalSum(objectiveDerivatives_*correction_);
380 }
381 
382 
384 {
385  return initialEtaSet_;
386 }
387 
388 
390 (
391  const scalarField& oldCorrection
392 )
393 {
394  correction_ = oldCorrection;
395 }
396 
397 
399 {
400  // Insert eta if set
401  if (initialEtaSet_)
402  {
403  optMethodIODict_.add<scalar>("eta", eta_, true);
404  }
405 
406  optMethodIODict_.add<scalarField>("correction", correction_, true);
407 
408  // Write IOdictionary
409  // Always write in ASCII format.
410  // Even when choosing to write in binary through controlDict,
411  // the content is written in ASCII format but with a binary header.
412  // This creates problems when the content is read back in
413  // (e.g. continuation)
414  optMethodIODict_.regIOobject::writeObject
415  (
416  IOstreamOption(IOstream::ASCII, mesh_.time().writeCompression()),
417  true
418  );
419 }
420 
421 
422 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::updateMethod::returnCorrection
scalarField & returnCorrection()
Return the correction of the design variables.
Definition: updateMethod.C:330
updateMethod.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::updateMethod::rightMult
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
Definition: updateMethod.C:72
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
s
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))
Definition: gmvOutputSpray.H:25
Foam::updateMethod::setGlobalSum
void setGlobalSum(const bool useGlobalSum)
Set globalSum variable.
Definition: updateMethod.C:324
Foam::updateMethod::setConstraintValues
void setConstraintValues(const scalarField &values)
Set constraints derivative.
Definition: updateMethod.C:312
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Matrix::n
label n() const noexcept
The number of columns.
Definition: MatrixI.H:103
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::updateMethod::meritFunctionDirectionalDerivative
virtual scalar meritFunctionDirectionalDerivative()
Definition: updateMethod.C:377
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::updateMethod::outerProd
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
Definition: updateMethod.C:99
invA
volSymmTensorField invA(inv(I *UEqn.A()+drag->Dcu()))
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::updateMethod::computeMeritFunction
virtual scalar computeMeritFunction()
Definition: updateMethod.C:371
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::updateMethod::setConstraintDeriv
void setConstraintDeriv(const PtrList< scalarField > &derivs)
Set constraints derivative.
Definition: updateMethod.C:298
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::updateMethod::writeCorrection
void writeCorrection()
Definition: updateMethod.C:337
n
label n
Definition: TABSMDCalcMethod2.H:31
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::updateMethod::updateOldCorrection
virtual void updateOldCorrection(const scalarField &oldCorrection)
Definition: updateMethod.C:390
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:63
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
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
dict
dictionary dict
Definition: searchingEngine.H:14
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:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::updateMethod::initialEtaSet
bool & initialEtaSet()
Return whether initial eta was set.
Definition: updateMethod.C:383
Foam::LUBacksubstitute
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
Definition: scalarMatricesTemplates.C:120
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:540
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::LUDecompose
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
Definition: scalarMatrices.C:34
Foam::updateMethod::setStep
void setStep(const scalar eta)
Set step for optimisation methods.
Definition: updateMethod.C:318
Foam::updateMethod::coeffsDict
dictionary coeffsDict()
Return optional dictionary with parameters specific to each method.
Definition: updateMethod.C:254
Foam::autoPtr< Foam::updateMethod >
Foam::SquareMatrix< scalar >
Foam::updateMethod::write
virtual void write()
Write useful quantities to files.
Definition: updateMethod.C:398
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::updateMethod::setObjectiveValue
void setObjectiveValue(const scalar value)
Set constraints derivative.
Definition: updateMethod.C:306
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::updateMethod::globalSum
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
Definition: updateMethod.C:175
Foam::List< label >
Foam::updateMethod::leftMult
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
Definition: updateMethod.C:45
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::updateMethod::inv
SquareMatrix< scalar > inv(SquareMatrix< scalar > A)
Definition: updateMethod.C:125
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::updateMethod::setObjectiveDeriv
void setObjectiveDeriv(const scalarField &derivs)
Set objective derivative.
Definition: updateMethod.C:291
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405