Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
31#include "OFstream.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37 defineTypeNameAndDebug(updateMethod, 0);
39}
40
41
42// * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
43
45(
46 const scalarField& s,
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
98Foam::SquareMatrix<Foam::scalar> Foam::updateMethod::outerProd
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 }
119
120 return (res);
121}
122
123
125Foam::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);
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 {
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 }
169 */
170
171 return invA;
172}
173
174
176{
177 if (globalSum_)
179 return gSum(field);
180 }
181 return sum(field);
182}
183
184
185Foam::scalar Foam::updateMethod::globalSum(tmp<scalarField>& tfield)
187 scalar value = globalSum(tfield());
188 tfield.clear();
189 return value;
190}
191
192
193Foam::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;
204 return res;
205}
206
207
208// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
209
210Foam::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:
219 localIOdictionary
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_(nullptr),
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)
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
315{
316 objectiveDerivatives_ = derivs;
317}
318
319
321(
323)
324{
325 constraintDerivatives_ = derivs;
326}
327
329void Foam::updateMethod::setObjectiveValue(const scalar value)
330{
331 objectiveValue_ = value;
332}
333
334
335void Foam::updateMethod::setObjectiveValueOld(const scalar value)
336{
337 if (!objectiveValueOld_)
339 objectiveValueOld_.reset(new scalar(Zero));
340 }
341 objectiveValueOld_.ref() = value;
342}
343
346{
347 cValues_ = values;
348}
349
350
352{
353 return objectiveValue_;
354}
355
356
361}
362
367}
368
370Foam::label Foam::updateMethod::getCycle() const
371{
372 return counter_;
373}
374
376void Foam::updateMethod::setStep(const scalar eta)
377{
378 eta_ = eta;
379}
380
382void Foam::updateMethod::modifyStep(const scalar multiplier)
383{
384 eta_ *= multiplier;
385}
386
388void Foam::updateMethod::setGlobalSum(const bool useGlobalSum)
389{
390 globalSum_ = useGlobalSum;
391}
392
397}
398
400Foam::label Foam::updateMethod::nConstraints() const
401{
402 return nConstraints_;
403}
404
407{
408 return correction_;
409}
410
411
413{
414 if (Pstream::master())
415 {
416 // Allocate cumulativeCorrection if necessary
417 if (cumulativeCorrection_.empty())
418 {
419 cumulativeCorrection_.setSize(correction_.size(), Zero);
420 }
421 // Accumulate correction
422 cumulativeCorrection_ += correction_;
423
424 fileName correctionFile
425 (
426 correctionFolder_/"correction" + mesh_.time().timeName()
427 );
428 fileName cumulativeCorrectionFile
429 (
430 correctionFolder_/"cumulativeCorrection" + mesh_.time().timeName()
431 );
432
433 OFstream corFile(correctionFile);
434 OFstream cumulCorFile(cumulativeCorrectionFile);
435 forAll(correction_, cI)
436 {
437 corFile
438 << cI << " " << correction_[cI] << endl;
439 cumulCorFile
440 << cI << " " << cumulativeCorrection_[cI] << endl;
441 }
442 }
443}
444
449}
450
455}
456
459{
460 return initialEtaSet_;
461}
462
463
465(
466 const scalarField& oldCorrection
467)
468{
469 correction_ = oldCorrection;
470}
471
472
474{
475 // Insert eta if set
476 if (initialEtaSet_)
477 {
478 os.writeEntry("eta", eta_);
479 }
480
481 os.writeEntry("counter", counter_);
482 correction_.writeEntry("correction", os);
483
484 return true;
485}
486
487
489{
490 // Does nothing in base
491 return true;
492}
493
494
495// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
label n
if(patchID !=-1)
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ LAZY_READ
Reading is optional [identical to READ_IF_PRESENT].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
fileName globalPath() const
The complete global path for the object (with instance, local,...).
Definition IOobject.C:512
label n() const noexcept
The number of columns.
Definition Matrix.H:271
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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...
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...
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
A class for handling file names.
Definition fileName.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
localIOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Abstract base class for optimisation methods.
void setConstraintValues(const scalarField &values)
Set values of constraints.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
label nConstraints_
Number of constraints.
virtual scalar meritFunctionDirectionalDerivative()
Directional derivative of the merit function, in the direction of the correction. Could be different ...
scalarField correction_
Design variables correction.
bool globalSum_
Whether to use gSum or sum in the inner products.
const fvMesh & mesh_
dictionary coeffsDict(const word &type) const
Return optional dictionary with parameters specific to each method.
scalar globalSum(const scalarField &field)
Compute either global or local sum, based on globalSum flag.
scalarField cValues_
Constraint values.
scalar objectiveValue_
Objective value.
virtual bool writeAuxiliaryData()
Write non-continuation data, usually under the optimisation folder.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
void setStep(const scalar eta)
Set step for optimisation methods.
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
const scalarField leftMult(const scalarField &, const SquareMatrix< scalar > &)
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...
scalarField & returnCorrection()
Return the correction of the design variables.
SquareMatrix< scalar > inv(SquareMatrix< scalar > A)
const labelList & activeDesignVars_
Map to active design variables.
scalar getObjectiveValue() const
Get objective value.
label getCycle() const
Get optimisation cycle.
bool initialEtaSet_
Is initially set?
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
void setObjectiveValue(const scalar value)
Set objective value.
autoPtr< scalar > objectiveValueOld_
Old objective value.
static autoPtr< updateMethod > New(const fvMesh &mesh, const dictionary &dict, autoPtr< designVariables > &designVars, const label nConstraints)
Return a reference to the selected turbulence model.
const autoPtr< scalar > & getObjectiveValueOld() const
Get old objective value.
void setObjectiveValueOld(const scalar value)
Set old objective value.
label counter_
Optimisation cycle count.
void setConstaintsNumber(const label nConstraints)
Set the number of constraints.
void setConstraintDeriv(const PtrList< scalarField > &derivs)
Set constraints derivative.
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
const dictionary dict_
const scalarField & getConstraintValues() const
Get values of constraints.
label nConstraints() const
Get the number of constraints.
void modifyStep(const scalar multiplier)
Multiply step.
scalarField cumulativeCorrection_
Cumulative design variables correction throughout the optimisation loop.
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-newton methods coupled with line search.
scalar eta_
Step multiplying the correction.
bool & initialEtaSet()
Return whether initial eta was set.
autoPtr< designVariables > & designVars_
Design variables.
word correctionFolder_
Folder storing the corrections to file.
SquareMatrix< scalar > outerProd(const scalarField &, const scalarField &)
void setObjectiveDeriv(const scalarField &derivs)
Set objective derivative.
void setGlobalSum(const bool useGlobalSum)
Set globalSum variable.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volSymmTensorField invA(inv(I *UEqn.A()+drag->Dcu()))
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
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))
word timeName
Definition getTimeIndex.H:3
#define DebugInfo
Report an information message using Foam::Info.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition List.H:62
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution in the source.
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
errorManip< error > abort(error &err)
Definition errorManip.H:139
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299