Loading...
Searching...
No Matches
LBFGS.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-2020 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
30#include "LBFGS.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39 (
41 LBFGS,
43 );
44}
45
46
47// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48
50{
51 label nVars(activeDesignVars_.size());
52 for (label i = 0; i < nPrevSteps_; ++i)
53 {
54 if (!y_.get(i))
55 {
56 y_.set(i, new scalarField(nVars, Zero));
57 }
58 if (!s_.get(i))
59 {
60 s_.set(i, new scalarField(nVars, Zero));
61 }
62 if (found("y" + Foam::name(i)))
63 {
64 y_[i] = scalarField("y" + Foam::name(i), *this, nVars);
65 }
66 if (found("s" + Foam::name(i)))
67 {
68 s_[i] = scalarField("s" + Foam::name(i), *this, nVars);
69 }
70 }
71}
72
73
75{
76 if (counter_ > nPrevSteps_)
77 {
78 // Reorder list by moving pointers down the line
79 labelList newOrder(nPrevSteps_, -1);
80 newOrder[0] = nPrevSteps_ - 1;
81 for (label i = 1; i < nPrevSteps_; ++i)
82 {
83 newOrder[i] = i - 1;
84 }
85 list.reorder(newOrder);
86
87 // Fill in last element with the provided field
88 list[nPrevSteps_ - 1] = f;
89 }
90 else
91 {
92 list[counter_ - 1] = f;
93 }
94}
95
96
98(
99 const scalarField& derivatives,
100 const scalarField& derivativesOld
101)
102{
103 // Sanity checks
104 if
105 (
106 (derivatives.size() != derivativesOld.size())
107 || (derivatives.size() != designVars_().getVars().size())
108 )
109 {
111 << "Sizes of input derivatives and design variables do not match"
112 << exit(FatalError);
113 }
114
115 // Update list of y. Can only be done here since derivatives
116 // were not known at the end of the previous cycle
117 scalarField yRecent(derivatives - derivativesOld, activeDesignVars_);
118 // Update list of s.
119 // correction_ holds the previous correction
120 scalarField sActive(correctionOld_, activeDesignVars_);
121 applyDamping(yRecent, sActive);
122
123 pivotFields(y_, yRecent);
124 pivotFields(s_, sActive);
125}
126
127
129{
130 const scalar sy(globalSum(s*y));
131 if (useSDamping_)
132 {
133 const scalarField Hy(invHessianVectorProduct(y, counter_ - 1));
134 const scalar yHy(globalSum(y*Hy));
135 scalar theta(1);
136 if (sy < 0.2*yHy)
137 {
139 << "y*s is below threshold. Using damped form" << nl
140 << "sy, yHy " << sy << " " << yHy << endl;
141
142 theta = 0.8*yHy/(yHy - sy);
143 }
144 s = theta*s + (1 - theta)*Hy;
145 }
146 else if (useYDamping_)
147 {
148 const scalarField Bs(HessianVectorProduct(s, counter_ - 1));
149 const scalar sBs(globalSum(s*Bs));
150 scalar theta(1);
151 if (sy < 0.2*sBs)
152 {
154 << "y*s is below threshold. Using damped form" << nl
155 << "sy, sBs " << sy << " " << sBs << endl;
156
157 theta = 0.8*sBs/(sBs - sy);
158 }
159 y = theta*y + (1 - theta)*Bs;
162 << "Curvature index (sy) is " << sy << endl;
163}
164
165
168{
170}
171
172
175(
176 const scalarField& vector,
177 const label counter,
179)
180{
181 // Sanity checks
182 tmp<scalarField> tq(tmp<scalarField>::New(activeDesignVars_.size(), Zero));
183 scalarField& q = tq.ref();
184 label nv = designVars_().getVars().size();
185 label nav = activeDesignVars_.size();
186 if (vector.size() == nv)
187 {
188 q.map(vector, activeDesignVars_);
189 }
190 else if (vector.size() == nav)
191 {
192 q = vector;
193 }
194 else
195 {
197 << "Size of input vector "
198 << "(" << vector.size() << ") "
199 << "is equal to neither the number of design variabes "
200 << "(" << nv << ")"
201 << " nor that of the active design variables "
202 << "(" << nav << ")"
203 << exit(FatalError);
204 }
205
206 if (counter)
207 {
208 // L-BFGS two loop recursion
209 //~~~~~~~~~~~~~~~~~~~~~~~~~~
210 label nSteps(min(counter, nPrevSteps_));
211 label nLast(nSteps - 1);
212 scalarField a(nSteps, 0.);
213 scalarField r(nSteps, 0.);
214 for (label i = nLast; i > -1; --i)
215 {
216 r[i] = 1./globalSum(y_[i]*s_[i]);
217 a[i] = r[i]*globalSum(s_[i]*q);
218 q -= a[i]*y_[i];
219 }
220
221 scalar gamma =
222 globalSum(y_[nLast]*y_[nLast])/globalSum(y_[nLast]*s_[nLast]);
223 if (diag)
224 {
225 q /= (gamma + diag());
226 }
227 else
228 {
229 q /= gamma;
230 }
231
232 scalarField b(activeDesignVars_.size(), Zero);
233 for (label i = 0; i < nSteps; ++i)
234 {
235 b = r[i]*globalSum(y_[i]*q);
236 q += s_[i]*(a[i] - b);
237 }
238 }
239 else if (diag)
240 {
241 q /= (1 + diag());
243
244 return tq;
245}
246
247
250{
252}
253
254
257(
258 const scalarField& vector,
259 const label counter
260)
261{
262 addProfiling(LBFGS, "LBFGS::HessianVectorProduct");
263 // Sanity checks
264 tmp<scalarField> tq(tmp<scalarField>::New(activeDesignVars_.size(), Zero));
265 scalarField& q = tq.ref();
266
267 scalarField source;
268 if (vector.size() == designVars_().getVars().size())
269 {
270 source = scalarField(vector, activeDesignVars_);
271 }
272 else if (vector.size() == activeDesignVars_.size())
273 {
274 source = vector;
275 }
276 else
277 {
279 << "Size of input vector is equal to neither the number of "
280 << " design variabes nor that of the active design variables"
281 << exit(FatalError);
282 }
283
284 if (counter != 0)
285 {
286 const label nSteps(min(counter, nPrevSteps_));
287 const label nLast(nSteps - 1);
288 const scalar delta =
289 globalSum(y_[nLast]*y_[nLast])/globalSum(y_[nLast]*s_[nLast]);
290
291 // Product of the last matrix on the right with the input vector
292 scalarField SKsource(2*nSteps, Zero);
293 for (label i = 0; i < nSteps; ++i)
294 {
295 SKsource[i] = delta*globalSum(s_[i]*source);
296 SKsource[i + nSteps] = globalSum(y_[i]*source);
297 }
298
299 // Form the middle matrix to be inverted
300 SquareMatrix<scalar> M(2*nSteps, 2*nSteps, Zero);
301 for (label i = 0; i < nSteps; ++i)
302 {
303 // Lower diagonal part
304 M[nSteps + i][nSteps + i] = - globalSum(s_[i]*y_[i]);
305 // Upper left part
306 for (label j = 0; j < nSteps; ++j)
307 {
308 M[i][j] = delta*globalSum(s_[i]*s_[j]);
309 }
310 }
311
312 // Upper right and lower left parts
313 for (label j = 0; j < nSteps; ++j)
314 {
315 for (label i = j + 1; i < nSteps; ++i)
316 {
317 scalar value = globalSum(s_[i]*y_[j]);
318 M[i][j + nSteps] = value;
319 M[j + nSteps][i] = value;
320 }
321 }
323
324 // Product of the inverted middle matrix with the right vector
325 scalarField invMSource(rightMult(invM, SKsource));
326
327 // Left vector multiplication with the rest of contributions
328 // vag: parallel comms
329 forAll(q, i)
330 {
331 for (label j = 0; j < nSteps; ++j)
332 {
333 q[i] -=
334 delta*s_[j][i]*invMSource[j]
335 + y_[j][i]*invMSource[j + nSteps];
336 }
337 }
338
339 q += delta*source;
340 }
341 else
342 {
343 q = source;
344 }
345
346 return tq;
347}
348
349
351{
352 // Sanity checks
353 const label n(activeDesignVars_.size());
354 tmp<scalarField> tdiag(tmp<scalarField>::New(n, 1));
355 scalarField& diag = tdiag.ref();
356
357 if (counter_ != 0)
358 {
359 const label nSteps(min(counter_, nPrevSteps_));
360 const label nLast(nSteps - 1);
361 const scalar delta =
362 globalSum(y_[nLast]*y_[nLast])/globalSum(y_[nLast]*s_[nLast]);
363 diag *= delta;
364
365 // Form the middle matrix to be inverted
366 SquareMatrix<scalar> M(2*nSteps, 2*nSteps, Zero);
367 for (label i = 0; i < nSteps; ++i)
368 {
369 // Lower diagonal part
370 M[nSteps + i][nSteps + i] = - globalSum(s_[i]*y_[i]);
371 // Upper left part
372 for (label j = 0; j < nSteps; ++j)
373 {
374 M[i][j] = delta*globalSum(s_[i]*s_[j]);
375 }
376 }
377
378 // Upper right and lower left parts
379 for (label j = 0; j < nSteps; ++j)
380 {
381 for (label i = j + 1; i < nSteps; ++i)
382 {
383 scalar value = globalSum(s_[i]*y_[j]);
384 M[i][j + nSteps] = value;
385 M[j + nSteps][i] = value;
386 }
387 }
388
389 // Invert the matrix
391
392 // Product of the inverse of the middle matrix with the right vector
393 List<scalarField> MR(2*nSteps, scalarField(n, Zero));
394 for (label k = 0; k < n; ++k)
395 {
396 for (label i = 0; i < 2*nSteps; ++i)
397 {
398 for (label j = 0; j < nSteps; ++j)
399 {
400 MR[i][k] +=
401 invM[i][j]*delta*s_[j][k]
402 + invM[i][j + nSteps]*y_[j][k];
403 }
404 }
405 }
406
407 // Part of the Hessian diagonal computed by the multiplication
408 // of the above matrix with the left matrix of the recursive Hessian
409 // reconstruction
410 for (label k = 0; k < n; ++k)
411 {
412 for (label j = 0; j < nSteps; ++j)
413 {
414 diag[k] -=
415 delta*s_[j][k]*MR[j][k] + y_[j][k]*MR[j + nSteps][k];
416 }
417 }
419
420 return tdiag;
421}
422
423
426{
428}
429
430
433(
434 const scalarField& vector,
435 const label counter
436)
437{
438 // Sanity checks
439 tmp<scalarField> tq(tmp<scalarField>::New(activeDesignVars_.size(), Zero));
440 scalarField& q = tq.ref();
441
442 scalarField source;
443 if (vector.size() == designVars_().getVars().size())
444 {
445 source = scalarField(vector, activeDesignVars_);
446 }
447 else if (vector.size() == activeDesignVars_.size())
448 {
449 source = vector;
450 }
451 else
452 {
454 << "Size of input vector is equal to neither the number of "
455 << " design variabes nor that of the active design variables"
456 << exit(FatalError);
457 }
458
459 if (counter != 0)
460 {
461 const label nSteps(min(counter, nPrevSteps_));
462 const label nLast(nSteps - 1);
463 const scalar delta =
464 globalSum(y_[nLast]*y_[nLast])/globalSum(y_[nLast]*s_[nLast]);
465
466 // Product of the last matrix on the right with the input vector
467 scalarField YBSsource(nSteps, Zero);
468 for (label i = 0; i < nSteps; ++i)
469 {
470 YBSsource[i] = globalSum((y_[i] - delta*s_[i])*source);
471 }
472
473 // Form the middle matrix to be inverted
474 SquareMatrix<scalar> M(nSteps, nSteps, Zero);
475 for (label i = 0; i < nSteps; ++i)
476 {
477 // D part
478 M[i][i] += globalSum(s_[i]*y_[i]);
479 // (S^T)BS part
480 for (label j = 0; j < nSteps; ++j)
481 {
482 M[i][j] -= delta*globalSum(s_[i]*s_[j]);
483 }
484 }
485
486 // Upper right and lower left parts
487 for (label j = 0; j < nSteps; ++j)
488 {
489 for (label i = j + 1; i < nSteps; ++i)
490 {
491 scalar value = globalSum(s_[i]*y_[j]);
492 M[i][j] += value;
493 M[j][i] += value;
494 }
495 }
497
498 // Product of the inverted middle matrix with the right vector
499 scalarField invMSource(rightMult(invM, YBSsource));
500
501 // Left vector multiplication with the rest of contributions
502 // vag: parallel comms
503 forAll(q, i)
504 {
505 for (label j = 0; j < nSteps; ++j)
506 {
507 q[i] += (y_[j][i] - delta*s_[j][i])*invMSource[j];
508 }
509 }
510
511 q += delta*source;
512 }
513 else
514 {
515 q = source;
516 }
517
518 return tq;
519}
520
521
523{
524 // Sanity checks
525 const label n(activeDesignVars_.size());
526 tmp<scalarField> tdiag(tmp<scalarField>::New(n, 1));
527 scalarField& diag = tdiag.ref();
528
529 if (counter_ != 0)
530 {
531 const label nSteps(min(counter_, nPrevSteps_));
532 const label nLast(nSteps - 1);
533 const scalar delta =
534 globalSum(y_[nLast]*y_[nLast])/globalSum(y_[nLast]*s_[nLast]);
535 diag *= delta;
536
537 // Form the middle matrix to be inverted
538 SquareMatrix<scalar> M(nSteps, nSteps, Zero);
539 for (label i = 0; i < nSteps; ++i)
540 {
541 // D part
542 M[i][i] += globalSum(s_[i]*y_[i]);
543 // (S^T)BS part
544 for (label j = 0; j < nSteps; ++j)
545 {
546 M[i][j] -= delta*globalSum(s_[i]*s_[j]);
547 }
548 }
549
550 // Upper right and lower left parts
551 for (label j = 0; j < nSteps; ++j)
552 {
553 for (label i = j + 1; i < nSteps; ++i)
554 {
555 scalar value = globalSum(s_[i]*y_[j]);
556 M[i][j] += value;
557 M[j][i] += value;
558 }
559 }
561
562 // Product of the inverse of the middle matrix with the right vector
563 List<scalarField> MR(nSteps, scalarField(n, Zero));
564 for (label k = 0; k < n; ++k)
565 {
566 for (label i = 0; i < nSteps; ++i)
567 {
568 for (label j = 0; j < nSteps; ++j)
569 {
570 MR[i][k] += invM[i][j]*(y_[j][k] - delta*s_[j][k]);
571 }
572 }
573 }
574
575 // Part of the Hessian diagonal computed by the multiplication
576 // of the above matrix with the left matrix of the recursive Hessian
577 // reconstruction
578 for (label k = 0; k < n; ++k)
579 {
580 for (label j = 0; j < nSteps; ++j)
581 {
582 diag[k] += (y_[j][k] - delta*s_[j][k])*MR[j][k];
583 }
585 }
586
587 return tdiag;
588}
589
592{
594}
595
596
598{
599 // In the first few iterations, use steepest descent but update the Hessian
600 // matrix
601 if (counter_ < nSteepestDescent_)
602 {
603 Info<< "Using steepest descent to update design variables" << endl;
604 for (const label varI : activeDesignVars_)
605 {
606 correction_[varI] = -eta_*objectiveDerivatives_[varI];
607 }
608 }
609 // else use LBFGS formula to update the design variables
610 else
611 {
612 scalarField q(invHessianVectorProduct(objectiveDerivatives_));
613 forAll(activeDesignVars_, varI)
614 {
615 correction_[activeDesignVars_[varI]] = -etaHessian_*q[varI];
616 }
617 }
618
619 // Store fields for the next iteration
622}
623
624
625// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
626
627Foam::LBFGS::LBFGS
628(
629 const fvMesh& mesh,
630 const dictionary& dict,
631 autoPtr<designVariables>& designVars,
632 const label nConstraints,
633 const word& type
634)
635:
636 quasiNewton(mesh, dict, designVars, nConstraints, type),
637 nPrevSteps_(coeffsDict(type).getOrDefault<label>("nPrevSteps", 10)),
638 y_(nPrevSteps_),
639 s_(nPrevSteps_),
640 useSDamping_(coeffsDict(type).getOrDefault<bool>("useSDamping", false)),
641 useYDamping_(coeffsDict(type).getOrDefault<bool>("useYDamping", false))
642{
643 // Allocate the correct sizes for y and s
645}
646
647
648// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
649
651{
652 // Write each component of y and s as a separate field so as to allow for
653 // reading them also in binary, since PtrList does not support this
654 forAll(y_, i)
655 {
656 y_[i].writeEntry(word("y" + Foam::name(i)), os);
657 s_[i].writeEntry(word("s" + Foam::name(i)), os);
658 }
659
661}
662
663
664// ************************************************************************* //
scalar y
label k
scalar delta
#define M(I)
bool found
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
The quasi-Newton Limited-memory BFGS formula. Keeps nPrevSteps_ of the y and s vectors and approximat...
Definition LBFGS.H:53
PtrList< scalarField > y_
The previous differences of derivatives. Holds nPrevSteps_ fields.
Definition LBFGS.H:81
void pivotFields(PtrList< scalarField > &list, const scalarField &f)
Move pointers in PtrList to the left and replace last element with given field.
Definition LBFGS.C:67
void applyDamping(scalarField &y, scalarField &s)
Use the damped version of s to ensure positive-definitiveness.
Definition LBFGS.C:121
void updateVectors(const scalarField &derivatives, const scalarField &derivativesOld)
Update y and s vectors.
Definition LBFGS.C:91
void allocateVectors()
Allocate matrices in the first optimisation cycle.
Definition LBFGS.C:42
virtual tmp< scalarField > HessianVectorProduct(const scalarField &vector)
Compute the Hessian - vector product.
Definition LBFGS.C:242
virtual tmp< scalarField > SR1HessianVectorProduct(const scalarField &vector)
Compute the Hessian - vector product.
Definition LBFGS.C:418
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition LBFGS.C:643
PtrList< scalarField > s_
The previous correction. Holds nPrevSteps_ fields.
Definition LBFGS.H:86
tmp< scalarField > HessianDiag()
Return the diagonal of the Hessian.
Definition LBFGS.C:343
bool useYDamping_
Use damping for s to ensure positive-definitiveness.
Definition LBFGS.H:96
label nPrevSteps_
Number of old corrections and grad differences kept.
Definition LBFGS.H:76
bool useSDamping_
Use damping for s to ensure positive-definitiveness.
Definition LBFGS.H:91
virtual tmp< scalarField > invHessianVectorProduct(const scalarField &vector)
Compute the inverse Hessian - vector product.
Definition LBFGS.C:160
virtual void update()
Update design variables.
Definition LBFGS.C:590
tmp< scalarField > SR1HessianDiag()
Return the diagonal of the Hessian.
Definition LBFGS.C:515
virtual void updateHessian()
Update the Hessian matrix by updating the base vectors.
Definition LBFGS.C:584
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
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
void reorder(const labelUList &oldToNew, const bool check=false)
Reorder elements. Reordering must be unique (ie, shuffle).
Definition UPtrList.C:62
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base class for quasi-Newton methods.
Definition quasiNewton.H:52
scalarField correctionOld_
The previous correction.
Definition quasiNewton.H:88
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
label nSteepestDescent_
Number of first steepest descent steps.
Definition quasiNewton.H:65
scalar etaHessian_
Step for the Newton method.
Definition quasiNewton.H:60
scalarField derivativesOld_
The previous derivatives.
Definition quasiNewton.H:83
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Abstract base class for optimisation methods.
scalarField objectiveDerivatives_
Derivatives of the objective functions.
scalarField correction_
Design variables correction.
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.
const labelList & activeDesignVars_
Map to active design variables.
const scalarField rightMult(const SquareMatrix< scalar > &, const scalarField &)
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
scalar eta_
Step multiplying the correction.
autoPtr< designVariables > & designVars_
Design variables.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
const scalar gamma
Definition EEqn.H:9
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
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))
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
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:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
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...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
labelList f(nPoints)
dictionary dict
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299