Loading...
Searching...
No Matches
ISQP.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) 2020-2023 PCOpt/NTUA
9 Copyright (C) 2020-2023 FOSS GP
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ISQP.H"
30#include "IOmanip.H"
31#include "Constant.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
42 ISQP,
44 );
46 (
48 ISQP,
50 );
51}
52
53
56({
57 { preconditioners::diag, "diag" },
58 { preconditioners::invHessian, "invHessian" },
59 { preconditioners::ShermanMorrison, "ShermanMorrison" }
60});
61
62
63// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64
66{
67 const label n = activeDesignVars_.size();
68 if (n != deltaP_.size())
69 {
70 // Correction fields
71 p_.setSize(n, Zero);
72 deltaP_.setSize(n, Zero);
73
74 // Lagrange multipliers and slack variables for bound constraints
76 {
77 lTilda_().setSize(n, Zero);
78 ls_().setSize(n, Zero);
79 uTilda_().setSize(n, Zero);
80 us_().setSize(n, Zero);
81
82 deltaLTilda_().setSize(n, Zero);
83 deltaLs_().setSize(n, Zero);
84 deltaUTilda_().setSize(n, Zero);
85 deltaUs_().setSize(n, Zero);
86 }
87
88 // Fields used to compute the Hessian
89 for (label i = 0; i < nPrevSteps_; ++i)
90 {
91 y_[i].setSize(n, Zero);
92 s_[i].setSize(n, Zero);
93 }
94 }
95}
96
97
99{
100 if (includeBoundConstraints_)
101 {
102 // Number of constraints
103 const label n(activeDesignVars_.size());
104
105 if (!lTilda_)
106 {
107 lTilda_.reset(autoPtr<scalarField>::New(n, Zero));
108 }
109 ls_.reset(autoPtr<scalarField>::New(n, Zero));
110 if (!uTilda_)
111 {
112 uTilda_.reset(autoPtr<scalarField>::New(n, Zero));
113 }
114 us_.reset(autoPtr<scalarField>::New(n, Zero));
115
116 deltaLTilda_.reset(autoPtr<scalarField>::New(n, Zero));
120 }
121}
122
123
125{
126 // Number of constraints
127 const label m(nConstraints_);
128 // Allocate the extra variables ensuring the feasibility
129 if (includeExtraVars_)
130 {
131 extraVars_.reset(autoPtr<scalarField>::New(m, 1));
132 const scalar t = mesh_.time().timeOutputValue();
133 z_.reset(autoPtr<scalarField>::New(m, max(1, 0.5*c_->value(t))));
134
135 deltaExtraVars_.reset(autoPtr<scalarField>::New(m, Zero));
137 }
138
139 doAllocateLagrangeMultipliers_ = false;
140}
141
142
144{
145 // Compute the Lagranian and old Lagrangian derivatives
146 scalarField LagrangianDerivativesOld(derivativesOld_);
147 forAll(constraintDerivatives_, cI)
148 {
149 LagrangianDerivatives_ += lamdas_[cI]*constraintDerivatives_[cI];
150 LagrangianDerivativesOld += lamdas_[cI]*constraintDerivativesOld_[cI];
151 }
152
153 if (includeBoundConstraints_)
154 {
155 forAll(activeDesignVars_, aI)
156 {
157 const label varI(activeDesignVars_[aI]);
158 const scalar contr(uTilda_()[aI] - lTilda_()[aI]);
159 LagrangianDerivatives_[varI] += contr;
160 LagrangianDerivativesOld[varI] += contr;
161 }
163
164 // Update vectors associated to the inverse Hessian matrix
165 updateVectors(LagrangianDerivatives_, LagrangianDerivativesOld);
166}
167
168
170{
171 const scalarField x(designVars_().getVars(), activeDesignVars_);
172
173 // Quantities related to design variables
174 p_ = Zero;
175 if (includeBoundConstraints_)
176 {
177 lTilda_() = scalar(1);
178 uTilda_() = scalar(1);
179 ls_() = scalar(1);
180 us_() = scalar(1);
181 }
182
183 // Quantities related to constraints
184 lamdas_ = scalar(1);
185 gs_ = scalar(1);
186
187 if (includeExtraVars_)
188 {
189 extraVars_() = scalar(1);
190 const scalar c = c_->value(mesh_.time().timeOutputValue());
191 z_() = max(1, 0.5*c);
192 Info<< "Penalty multiplier (c):: " << c << endl;
194
195 // Reset eps
196 eps_ = 1;
197}
198
199
201{
202 deltaP_ = Zero;
203 deltaLamda_ = Zero;
204 deltaGs_ = Zero;
205
206 if (includeBoundConstraints_)
207 {
208 deltaLTilda_() = Zero;
209 deltaLs_() = Zero;
210 deltaUTilda_() = Zero;
211 deltaUs_() = Zero;
212 }
213
214 if (includeExtraVars_)
217 deltaZ_() = Zero;
218 }
219}
220
221
223{
224 addProfiling(ISQP, "ISQP::solveDeltaPEqn");
225 // Explicit part of the right hand side of the deltaX equation
226 scalarField FDx(-resFL());
227 if (includeBoundConstraints_)
228 {
229 FDx +=
230 (uTilda_()*resFus() + resFuTilda())/us_()
231 - (lTilda_()*resFls() + resFlTilda())/ls_();
232 }
233 scalarField AMult(resFlamda()/lamdas_ - resFGs());
234 scalarField mult(gs_/lamdas_);
235 if (includeExtraVars_)
236 {
237 mult += extraVars_()/z_();
238 AMult -= (extraVars_()*resFExtraVars() + resFz())/z_();
239 }
240 AMult /= mult;
241 forAll(FDx, aI)
242 {
243 const label varI(activeDesignVars_[aI]);
244 forAll(constraintDerivatives_, cI)
245 {
246 FDx[aI] += constraintDerivatives_[cI][varI]*AMult[cI];
247 }
248 }
249 CGforDeltaP(FDx);
250}
251
252
253Foam::tmp<Foam::scalarField> Foam::ISQP::computeRHSForDeltaX
254(
255 const scalarField& FDx
256)
257{
259 scalarField& rhs = trhs.ref();
260
261 // Compute (Gs)^(-1)*Λ*A*Dp
262 scalarField GsLADp(cValues_.size(), Zero);
263 forAll(constraintDerivatives_, cI)
264 {
265 const scalarField& cDerivsI = constraintDerivatives_[cI];
266 GsLADp[cI] +=
267 globalSum(scalarField(cDerivsI, activeDesignVars_)*deltaP_);
268 }
269 GsLADp *= lamdas_/gs_;
270
271 // Multiply with A^T
272 forAll(rhs, aI)
273 {
274 const label varI(activeDesignVars_[aI]);
275 forAll(constraintDerivatives_, cI)
276 {
277 rhs[aI] += constraintDerivatives_[cI][varI]*GsLADp[cI];
278 }
279 }
280
281 // Contributions from bounds
282 if (includeBoundConstraints_)
283 {
284 rhs += (lTilda_()/ls_() + uTilda_()/us_())*deltaP_;
285 }
286
288
289 rhs = 0.95*deltaP_ + 0.05*rhs;
290 return trhs;
291}
292
293
295{
296 addProfiling(ISQP, "ISQP::CGforDeltaP");
297 autoPtr<scalarField> precon(nullptr);
298 scalarField r(FDx - matrixVectorProduct(deltaP_));
299 scalarField z(preconVectorProduct(r, precon));
300 scalarField p(z);
301 scalar res(sqrt(globalSum(r*r)));
302 scalar resInit(res);
303 scalar rz(globalSum(r*z));
304 scalar rzOld(rz);
305 label iter(0);
306 do
307 {
308 scalarField Ap(matrixVectorProduct(p));
309 scalar a = rz/globalSum(p*Ap);
310 deltaP_ += a*p;
311 r -= a*Ap;
312 res = sqrt(globalSum(r*r));
313 z = preconVectorProduct(r, precon);
314 rz = globalSum(r*z);
315 scalar beta = rz/rzOld;
316 p = z + beta*p;
317 rzOld = rz;
318 } while (iter++ < maxDxIters_ && res > relTol_*resInit);
320 << "CG, Solving for deltaP, Initial Residual " << resInit
321 << ", Final Residual " << res
322 << ", No Iterations " << iter << endl;
323}
324
325
327(
328 const scalarField& vector
329)
330{
331 addProfiling(ISQP, "ISQP::MatrixVectorProduct");
332 tmp<scalarField> tAp(HessianVectorProduct(vector));
333 scalarField& Ap = tAp.ref();
334 scalarField GsLAv(cValues_.size(), Zero);
335 forAll(constraintDerivatives_, cI)
336 {
337 const scalarField& cDerivsI = constraintDerivatives_[cI];
338 GsLAv[cI] =
339 globalSum(scalarField(cDerivsI, activeDesignVars_)*vector);
340 }
341 scalarField mult(gs_/lamdas_);
342 if (includeExtraVars_)
343 {
344 mult += extraVars_()/z_();
345 }
346 GsLAv /= mult;
347
348 // Multiply with A^T
349 forAll(Ap, aI)
350 {
351 const label varI(activeDesignVars_[aI]);
352 forAll(constraintDerivatives_, cI)
353 {
354 Ap[aI] += constraintDerivatives_[cI][varI]*GsLAv[cI];
355 }
356 }
357
358 // Contributions from bounds
359 if (includeBoundConstraints_)
360 {
361 Ap += (lTilda_()/ls_() + uTilda_()/us_())*vector;
363
364
365 return tAp;
366}
367
368
370(
371 const scalarField& vector,
372 autoPtr<scalarField>& precon
373)
374{
375 addProfiling(ISQP, "ISQP::preconVectorProduct");
376 if (preconType_ == preconditioners::diag)
377 {
378 if (!precon)
379 {
380 precon.reset(diagPreconditioner().ptr());
381 }
382 return precon()*vector;
383 }
384 else if (preconType_ == preconditioners::invHessian)
385 {
386 return invHessianVectorProduct(vector);
387 }
388 else if (preconType_ == preconditioners::ShermanMorrison)
391 }
392 return nullptr;
393}
394
395
397{
398 addProfiling(ISQP, "ISQP::deltaPDiagPreconditioner");
399 tmp<scalarField> tpreconditioner(HessianDiag());
400 scalarField& preconditioner = tpreconditioner.ref();
401
402 // Part related to the constraints
403 forAll(constraintDerivatives_, cI)
404 {
405 scalarField cDerivs(constraintDerivatives_[cI], activeDesignVars_);
406 scalar mult(gs_[cI]/lamdas_[cI]);
407 if (includeExtraVars_)
408 {
409 mult += extraVars_()[cI]/z_()[cI];
410 }
411 preconditioner += sqr(cDerivs)/mult;
412 }
413
414 if (includeBoundConstraints_)
415 {
416 preconditioner += lTilda_()/ls_() + uTilda_()/us_();
417 }
419 preconditioner = 1./preconditioner;
420
421 return tpreconditioner;
422}
423
424
426(
427 const scalarField& vector
428)
429{
430 // Recursively update the inv(LHS)*vector since the LHS consists of the
431 // L-BFGS-based Hessian, computed with rank-2 updates, and the part related
432 // to flow constraints, computed as rank-1 updates. In the inversion
433 // process, the diagonal matrix related to bound constraints is treated as
434 // the initial matrix of the L-BFGS update.
435
436 // Contribution from bound constraints, treated as the seed of the
437 // L-BFGS inverse
438 refPtr<scalarField> tdiag(nullptr);
439 if (includeBoundConstraints_)
440 {
441 tdiag.reset((lTilda_()/ls_() + uTilda_()/us_()).ptr());
442 }
443
444 // List of vectors to be used in the rank-1 updates related to the flow
445 // constraitns
446 PtrList<scalarField> r1Updates(cValues_.size());
447
448 forAll(constraintDerivatives_, cI)
449 {
450 const scalarField& cDerivsI = constraintDerivatives_[cI];
451 r1Updates.set(cI, new scalarField(cDerivsI, activeDesignVars_));
452 }
453 // Multipliers of the rank-1 updates
454 scalarField mult(gs_/lamdas_);
455 if (includeExtraVars_)
456 {
457 mult += extraVars_()/z_();
459
460 return
461 ShermanMorrisonRank1Update(r1Updates, vector, tdiag, mult, mult.size());
462}
463
464
466(
467 const PtrList<scalarField>& r1Updates,
468 const scalarField& p,
469 const refPtr<scalarField>& diag,
470 const scalarField& mult,
471 label n
472)
473{
474 auto tAp(tmp<scalarField>::New(activeDesignVars_.size(), Zero));
475 scalarField& Ap = tAp.ref();
476
477 if (n == 0)
478 {
479 Ap = invHessianVectorProduct(p, counter_, diag.shallowClone());
480 return tAp;
481 }
482
483 do
484 {
485 --n;
486 Ap = ShermanMorrisonRank1Update(r1Updates, p, diag, mult, n);
487 tmp<scalarField> tAv =
488 ShermanMorrisonRank1Update(r1Updates, r1Updates[n], diag, mult, n);
489 scalarField& Av = tAv.ref();
490 scalar yHs = globalSum(r1Updates[n]*Av)/mult[n];
491 Ap -= Av*globalSum(r1Updates[n]*Ap)/(1 + yHs)/mult[n];
492 } while (n > 0);
493 return tAp;
494}
495
496
498{
499 // Zero the updates computed in the previous optimisation cycle
500 //zeroUpdates();
501
502 // Solve equation for deltaP_. The expensive part of the step. Everything
503 // else can be computed based on this
504 addProfiling(ISQP, "ISQP::computeNewtonDirection");
505 solveDeltaPEqn();
506
507 // deltaLamda
508 forAll(constraintDerivatives_, cI)
509 {
510 const scalarField& cDerivsI = constraintDerivatives_[cI];
511 deltaLamda_[cI] =
512 globalSum(scalarField(cDerivsI, activeDesignVars_)*deltaP_);
513 }
514 scalarField mult(gs_/lamdas_);
515 if (includeExtraVars_)
516 {
517 mult += extraVars_()/z_();
518 deltaLamda_ += (resFz() + extraVars_()*resFExtraVars())/z_();
519 }
520 deltaLamda_ += resFGs() - resFlamda()/lamdas_;
521 deltaLamda_ /= mult;
522
523 // deltaGs
524 deltaGs_ = -(gs_*deltaLamda_ + resFlamda())/lamdas_;
525
526 if (includeBoundConstraints_)
527 {
528 // deltaLs
529 deltaLs_() = deltaP_ + resFls();
530
531 // deltaUs
532 deltaUs_() = -deltaP_ + resFus();
533
534 // deltaLTilda
535 deltaLTilda_() = -(lTilda_()*deltaLs_() + resFlTilda())/ls_();
536
537 // deltaUTilda
538 deltaUTilda_() = -(uTilda_()*deltaUs_() + resFuTilda())/us_();
539 }
540
541 if (includeExtraVars_)
544 deltaExtraVars_() = - (extraVars_()*deltaZ_() + resFz())/z_();
545 }
546}
547
548
549Foam::scalar Foam::ISQP::lineSearch()
550{
551 const label n(p_.size());
552 const label m(cValues_.size());
553 scalar step(1.);
554
555 if (includeBoundConstraints_)
556 {
557 for (label i = 0; i < n; ++i)
558 {
559 adjustStep(step, ls_()[i], deltaLs_()[i]);
560 adjustStep(step, us_()[i], deltaUs_()[i]);
561 adjustStep(step, lTilda_()[i], deltaLTilda_()[i]);
562 adjustStep(step, uTilda_()[i], deltaUTilda_()[i]);
563 }
564 }
565
566 // Perform bound checks and adjust step accordingly
567 for (label i = 0; i < m; ++i)
568 {
569 adjustStep(step, lamdas_[i], deltaLamda_[i]);
570 adjustStep(step, gs_[i], deltaGs_[i]);
571 if (includeExtraVars_)
572 {
573 adjustStep(step, extraVars_()[i], deltaExtraVars_()[i]);
574 adjustStep(step, z_()[i], deltaZ_()[i]);
575 }
576 }
577
578 // Each processor might have computed a different step, if design variables
579 // are distributed. Get the global minimum
580 if (globalSum_)
581 {
582 reduce(step, minOp<scalar>());
583 }
584
585 step = min(1, step);
586
587 if (debug > 1)
588 {
589 Info<< "Step before line search is " << step << endl;
590 }
591
592 // Old residual
593 scalar normResOld = sqrt(globalSum(magSqr(computeResiduals())));
594 scalar maxRes(GREAT);
595
596 for (label i = 0; i < maxLineSearchIters_ ; ++i)
597 {
598 // Update the solution with given step
599 updateSolution(step);
600
601 // Compute new residuals and their max value
602 scalarField resNew(computeResiduals());
603 scalar normResNew = sqrt(globalSum(magSqr(resNew)));
604 maxRes = gMax(mag(resNew));
605
606 if (normResNew < normResOld)
607 {
609 << "Initial residual = " << normResOld << ", "
610 << "Final residual = " << normResNew << ", "
611 << "No of LineSearch Iterations = " << i + 1
612 << endl;
613 break;
614 }
615 else
616 {
617 // Return solution to previous and reduce step
618 if (i != maxLineSearchIters_ - 1)
619 {
620 updateSolution(-step);
621 step *= 0.5;
622 }
623 else
624 {
625 Info<< tab << "Line search did not converge. "
626 << "Procceding with the last compute step" << endl;
627 }
628 }
629 }
630
631 if (debug > 1)
632 {
633 Info<< "Step after line search is " << step << nl << endl;
634 }
635
636 return maxRes;
637}
638
639
641(
642 scalar& step,
643 const scalar value,
644 const scalar update
645)
646{
647 if (0.99*value + step*update < scalar(0))
648 {
649 step = -0.99*value/update;
650 }
651}
652
653
654void Foam::ISQP::updateSolution(const scalar step)
655{
656 p_ += step*deltaP_;
657 lamdas_ += step*deltaLamda_;
658 gs_ += step*deltaGs_;
659 if (includeBoundConstraints_)
660 {
661 lTilda_() += step*deltaLTilda_();
662 ls_() += step*deltaLs_();
663 uTilda_() += step*deltaUTilda_();
664 us_() += step*deltaUs_();
665 }
666 if (includeExtraVars_)
668 extraVars_() += step*deltaExtraVars_();
669 z_() += step*deltaZ_();
670 }
671}
672
673
675{
676 const label n(activeDesignVars_.size());
677 const label m(cValues_.size());
678 label size(includeBoundConstraints_ ? 5*n + 2*m : n + 2*m);
679 if (includeExtraVars_)
680 {
681 size += 2*m;
682 }
684 scalarField& res = tres.ref();
685
686 label iRes(0);
687
688 // Gradient of the Lagrangian
689 res.rmap(resFL()(), identity(n));
690 iRes = n;
691
692 // Inequality constraints slacks
693 res.rmap(resFGs()(), identity(m, iRes));
694 iRes += m;
695
696 // Inequality constraints complementarity slackness
697 res.rmap(resFlamda()(), identity(m, iRes));
698 iRes += m;
699
700 if (includeBoundConstraints_)
701 {
702 // Lower bounds slacks
703 res.rmap(resFls()(), identity(n, iRes));
704 iRes += n;
705
706 // Upper bounds slacks
707 res.rmap(resFus()(), identity(n, iRes));
708 iRes += n;
709
710 // Lower bounds complementarity slackness
711 res.rmap(resFlTilda()(), identity(n, iRes));
712 iRes += n;
713
714 // Upper bounds complementarity slackness
715 res.rmap(resFuTilda()(), identity(n, iRes));
716 iRes += n;
717 }
718
719 if (includeExtraVars_)
720 {
721 // Lagragian derivative wrt the extra variables
722 res.rmap(resFExtraVars()(), identity(m, iRes));
723 iRes += m;
724
725 // Lagrange multipliers for the extra variables positiveness
726 res.rmap(resFz(), identity(m, iRes));
727 iRes += m;
728 }
729
730 return tres;
731}
732
733
735{
736 tmp<scalarField> tgradL
737 (tmp<scalarField>::New(objectiveDerivatives_, activeDesignVars_));
738 scalarField& gradL = tgradL.ref();
739
740 scalarField Hp(HessianVectorProduct(p_));
741 //scalarField Hp = SR1HessianVectorProduct(p_);
742 gradL += Hp;
743
744 if (debug > 2)
745 {
746 scalarField H1Hp(invHessianVectorProduct(Hp));
747 Info << "Diff H1Hp - p " << gSum(mag(H1Hp - p_)) << endl;
748 }
749
750 forAll(constraintDerivatives_, cI)
751 {
752 gradL +=
753 lamdas_[cI]
754 *scalarField(constraintDerivatives_[cI], activeDesignVars_);
755 }
756
757 if (includeBoundConstraints_)
758 {
759 gradL += uTilda_() - lTilda_();
760 }
761
762 return tgradL;
763}
764
765
767{
768 tmp<scalarField> tinvHFL
769 (tmp<scalarField>::New(objectiveDerivatives_, activeDesignVars_));
770 scalarField& invHFL = tinvHFL.ref();
771
772 forAll(constraintDerivatives_, cI)
773 {
774 invHFL +=
775 lamdas_[cI]
776 *scalarField(constraintDerivatives_[cI], activeDesignVars_);
777 }
778
779 if (includeBoundConstraints_)
780 {
781 invHFL += uTilda_() - lTilda_();
782 }
783
785 invHFL += p_;
786
787 return tinvHFL;
788}
789
790
792{
794 (
795 tmp<scalarField>::New(gs_ + cValues_ - max((1 - cRed_)*cValues_, Zero))
796 );
797 scalarField& FGs = tFGs.ref();
798
799 forAll(constraintDerivatives_, cI)
800 {
801 FGs[cI] +=
802 globalSum
803 (
804 scalarField(constraintDerivatives_[cI], activeDesignVars_)*p_
805 );
806 }
807
808 if (includeExtraVars_)
809 {
810 FGs -= extraVars_();
811 }
812
813 return tFGs;
814}
815
818{
819 return (lamdas_*gs_ - eps_);
820}
821
822
824{
825 if (includeBoundConstraints_)
826 {
827 const scalarField x(designVars_().getVars(), activeDesignVars_);
828 const scalarField xMin
829 (designVars_().lowerBounds()(), activeDesignVars_);
831 return (x + p_ - xMin - ls_());
832 }
833 return nullptr;
834}
835
836
838{
839 if (includeBoundConstraints_)
840 {
841 const scalarField x(designVars_().getVars(), activeDesignVars_);
842 const scalarField xMax
843 (designVars_().upperBounds()(), activeDesignVars_);
845 return (xMax - x - p_ - us_());
846 }
847 return nullptr;
848}
849
850
852{
853 if (includeBoundConstraints_)
855 return (lTilda_()*ls_() - eps_);
856 }
857 return nullptr;
858}
859
860
862{
863 if (includeBoundConstraints_)
865 return (uTilda_()*us_() - eps_);
866 }
867 return nullptr;
868}
869
870
872{
873 if (includeExtraVars_)
875 return (c_->value(mesh_.time().timeOutputValue()) - lamdas_ - z_());
876 }
877 return nullptr;
878}
879
880
882{
883 if (includeExtraVars_)
885 return (z_()*extraVars_() - eps_);
886 }
887 return nullptr;
888}
889
890
892{
893 addProfiling(ISQP, "ISQP::solveSubproblem");
894 zeroUpdates();
895 if (includeBoundConstraints_ || !cValues_.empty())
896 {
897 scalar resMax(gMax(mag(computeResiduals())));
898 label iter(0);
899 do
900 {
902 << "Newton iter " << iter << nl << endl;
903
904 // Decrease eps
905 if (resMax < 0.9*eps_)
906 {
907 eps_ *= 0.1;
908 }
909
910 // Computes Newton direction for the subproblem
911 computeNewtonDirection();
912
913 // Perform line search and return max residual of the solution
914 // satisfying the bound constraints and the residual reduction.
915 // Upates solution.
916 resMax = lineSearch();
918 << "max residual = " << resMax << ", "
919 << "eps = " << eps_ << nl << endl;
920 } while
921 (
922 iter++ < maxNewtonIters_ && (eps_ > epsMin_ || resMax > 0.9*eps_)
923 );
924 Info<< "Finished solving the QP subproblem" << nl << endl;
925 if (iter == maxNewtonIters_)
926 {
928 << "Iterative solution of the QP problem did not converge"
929 << endl;
930 }
931 if (debug > 2)
932 {
933 scalarField vars(designVars_().getVars(), activeDesignVars_);
934 scalarField newVars(vars + p_);
935 Info<< "Min of updated vars " << gMin(newVars) << endl;
936 Info<< "Max of updated vars " << gMax(newVars) << endl;
937
938 Info<< "Min of lamda " << gMin(lamdas_) << endl;
939 Info<< "Max of gs " << gMax(gs_) << endl;
940 Info<< "Max of lamda*gs " << gMax(lamdas_*gs_) << endl;
941
942 if (includeBoundConstraints_)
943 {
944 Info<< "Min of lTilda " << gMin(lTilda_()) << endl;
945 Info<< "Min of uTilda " << gMin(uTilda_()) << endl;
946 Info<< "Min of ls " << gMin(ls_()) << endl;
947 Info<< "Min of us " << gMin(us_()) << endl;
948 Info<< "Max of lTilda*ls " << gMax(lTilda_()*ls_()) << endl;
949 Info<< "Max of uTilda*us " << gMax(uTilda_()*us_()) << endl;
950 }
951 if (includeExtraVars_)
952 {
953 Info<< "Min of extraVars " << gMin(extraVars_()) << endl;
954 Info<< "Max of extraVars*z " << gMax(extraVars_()*z_()) << endl;
955 }
956 }
957 if (includeExtraVars_)
958 {
960 << "Constraint penalization variables (y) " << extraVars_()
961 << endl;
962 }
963 }
964 else
965 {
966 computeNewtonDirection();
967 lineSearch();
968 }
969
970 // Pass update to correction field
971 correction_.rmap(p_, activeDesignVars_);
972 if (!counter_)
973 {
974 correction_ *= eta_;
975 }
976 else
977 {
979 }
980}
981
982
984{
985 derivativesOld_ = objectiveDerivatives_;
986 if (constraintDerivativesOld_.empty())
987 {
988 constraintDerivativesOld_.setSize(constraintDerivatives_.size());
989 }
990 forAll(constraintDerivativesOld_, cI)
993 }
994 correctionOld_ = correction_;
995}
996
997
999{
1000 // Assumes that all constraints are known by all processors
1001 // What about constraints directly imposed on distributed design variables?
1002 // These should be met in each iteration of the algorithm, so,
1003 // most probably, there is no problem
1004 return sum(pos(cValues_)*cValues_);
1005}
1006
1007
1008// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1009
1010Foam::ISQP::ISQP
1011(
1012 const fvMesh& mesh,
1013 const dictionary& dict,
1014 autoPtr<designVariables>& designVars,
1015 const label nConstraints,
1016 const word& type
1017)
1018:
1019 LBFGS(mesh, dict, designVars, nConstraints, type),
1020 SQPBase(mesh, dict, designVars, *this, type),
1021
1022 doAllocateLagrangeMultipliers_(true),
1023 includeBoundConstraints_
1024 (
1025 designVars->upperBounds() && designVars->lowerBounds()
1026 ),
1027 includeExtraVars_
1028 (
1029 coeffsDict(type).getOrDefault<bool>("includeExtraVars", false)
1030 ),
1031 p_(activeDesignVars_.size(), Zero),
1032 gs_(nConstraints_, 1),
1033 lTilda_
1034 (
1035 includeBoundConstraints_ && found("lTilda") ?
1036 new scalarField("lTilda", *this, activeDesignVars_.size()) :
1037 nullptr
1038 ),
1039 ls_(nullptr),
1040 uTilda_
1041 (
1042 includeBoundConstraints_ && found("uTilda") ?
1043 new scalarField("uTilda", *this, activeDesignVars_.size()) :
1044 nullptr
1045 ),
1046 us_(nullptr),
1047 extraVars_(nullptr),
1048 z_(nullptr),
1049 c_(nullptr),
1050 deltaP_(activeDesignVars_.size(), Zero),
1051 deltaLamda_(nConstraints_, Zero),
1052 deltaGs_(nConstraints_, Zero),
1053 deltaLTilda_(nullptr),
1054 deltaLs_(nullptr),
1055 deltaUTilda_(nullptr),
1056 deltaUs_(nullptr),
1057 deltaExtraVars_(nullptr),
1058 deltaZ_(nullptr),
1059 eps_(1),
1060 epsMin_(coeffsDict(type).getOrDefault<scalar>("epsMin", 1.e-07)),
1061 maxNewtonIters_(coeffsDict(type).getOrDefault<label>("maxIters", 1000)),
1062 maxLineSearchIters_
1063 (
1064 coeffsDict(type).getOrDefault<label>("maxLineSearchIters", 10)
1065 ),
1066 maxDxIters_(coeffsDict(type).getOrDefault<label>("maxDpIters", 1000)),
1067 relTol_(coeffsDict(type).getOrDefault<scalar>("relTol", 0.01)),
1068 preconType_
1069 (
1070 preconditionerNames.getOrDefault
1071 (
1072 "preconditioner", coeffsDict(type), preconditioners::diag
1073 )
1074 ),
1075 cRed_
1076 (coeffsDict(type).getOrDefault<scalar>("targetConstraintReduction", 1)),
1077 disableDamping_
1078 (coeffsDict(type).getOrDefault<bool>("disableDamping", false)),
1079 meritFunctionFile_(nullptr)
1080{
1081 Info<< "Preconditioner type of the SQP subproblem is ::"
1083 << endl;
1084 if (!disableDamping_)
1085 {
1086 // Always apply damping of y in ISQP
1087 useYDamping_ = true;
1088 useSDamping_ = false;
1089 }
1090
1091 // Determine c if necessary
1093 {
1094 if (coeffsDict(type).found("c"))
1095 {
1096 Info<< "Reading constraint penalty function type from dict" << endl;
1097 c_.reset(Function1<scalar>::New("c", coeffsDict(type)));
1098 }
1099 else
1100 {
1101 Info<< "Setting constant penalty factor" << endl;
1102 c_.reset(new Function1Types::Constant<scalar>("c", 100));
1103 }
1104 }
1105
1106 // Allocate multipliers and slack variables for the bound constraints
1109}
1110
1111
1112// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1113
1115{
1116 // Update sizes of fields related to the active design variables
1117 updateSizes();
1118
1119 // The first iteration uses a unitary Hessian. No need to update
1120 LagrangianDerivatives_ = objectiveDerivatives_;
1121 if (counter_)
1122 {
1123 updateYS();
1124 }
1125
1126 // Initiaze variables
1127 initialize();
1128
1129 // Solve subproblem using a Newton optimiser
1130 solveSubproblem();
1131
1132 // Store fields for the next iteration and write them to file
1134
1135 // Increase counter
1136 ++counter_;
1137}
1138
1139
1141{
1144
1145 return L;
1146}
1147
1148
1150{
1151 scalar deriv =
1154
1155 return deriv;
1156}
1157
1158
1160{
1161 updateMethod::updateOldCorrection(oldCorrection);
1162 correctionOld_ = oldCorrection;
1163}
1164
1165
1166bool Foam::ISQP::writeData(Ostream& os) const
1167{
1168 if (includeBoundConstraints_)
1169 {
1170 uTilda_().writeEntry("uTilda", os);
1171 lTilda_().writeEntry("lTilda", os);
1172 }
1173
1175}
1176
1177
1179{
1180 return SQPBase::writeMeritFunction(*this);
1181}
1182
1183
1184// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
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
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
const List< word > & names() const noexcept
The list of enum names, in construction order. Same as toc().
Definition EnumI.H:39
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
Templated function that returns a constant value.
Definition Constant.H:71
An L-BFGS-based SQP algorithm for computing the update of the design variables in the presence of ine...
Definition ISQP.H:69
bool doAllocateLagrangeMultipliers_
Should Lagrange multipliers be allocated.
Definition ISQP.H:93
scalar lineSearch()
Perform line search and return max residual corresponding to the updated solution.
Definition ISQP.C:542
tmp< scalarField > matrixVectorProduct(const scalarField &vector)
Procudt of the LHS of the deltaP eqn with a vector.
Definition ISQP.C:320
void updateSolution(const scalar step)
Update the current solution using the known direction and the given step length.
Definition ISQP.C:647
scalar cRed_
Percentage reduction for the constraints.
Definition ISQP.H:219
tmp< scalarField > resFls()
Residual of the lower bounds slackness.
Definition ISQP.C:816
void CGforDeltaP(const scalarField &FDx)
CG algorithm for the solution of the deltaP eqn.
Definition ISQP.C:287
void solveSubproblem()
Solve subproblem using a Newton optimiser.
Definition ISQP.C:884
void computeCorrection()
Compute design variables correction.
Definition ISQP.C:1107
virtual scalar meritFunctionDirectionalDerivative()
Derivative of the merit function. Could be different than the objective derivative in the presence of...
Definition ISQP.C:1142
bool includeExtraVars_
Are additional design variables included?
Definition ISQP.H:106
tmp< scalarField > resFL()
Residual of the gradient of the Lagrangian.
Definition ISQP.C:727
void initialize()
Allocate fields related to constraints.
Definition ISQP.C:162
virtual bool writeAuxiliaryData()
Write merit function information.
Definition ISQP.C:1171
autoPtr< Function1< scalar > > c_
Multiplier of the additional variables y in the Lagrangian, to make them 'expensive'.
Definition ISQP.H:159
tmp< scalarField > invHFL()
Product of the inverse Hessian with the residual of the gradient of the Lagrangian.
Definition ISQP.C:759
tmp< scalarField > resFz()
Residual of the complementarity slackness for the extra variables.
Definition ISQP.C:874
void adjustStep(scalar &step, const scalar value, const scalar update)
Adjust step to satisfy cireteria.
Definition ISQP.C:634
autoPtr< scalarField > z_
Lagrange multipliers for positive extra variables.
Definition ISQP.H:153
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition ISQP.C:1159
tmp< scalarField > resFExtraVars()
Residual of the Lagrangian derivative wrt the extra variables.
Definition ISQP.C:864
scalarField deltaLamda_
Definition ISQP.H:165
virtual scalar meritFunctionConstraintPart() const
Get the part the merit function that depends on the constraints.
Definition ISQP.C:991
tmp< scalarField > preconVectorProduct(const scalarField &vector, autoPtr< scalarField > &precon)
Preconditioner-vector product for CG.
Definition ISQP.C:363
static const Enum< preconditioners > preconditionerNames
Names of preconditioners for the subproblem.
Definition ISQP.H:83
scalarField p_
The set of design variables being updated during the subproblem.
Definition ISQP.H:115
void computeNewtonDirection()
Compute direction for the Newton problem.
Definition ISQP.C:490
autoPtr< scalarField > us_
Slack variables the upper bound constraints.
Definition ISQP.H:142
scalarField deltaGs_
Definition ISQP.H:166
scalar epsMin_
Final eps quantity to be reached during the solution of the subproblem.
Definition ISQP.H:186
autoPtr< scalarField > deltaUs_
Definition ISQP.H:170
tmp< scalarField > resFuTilda()
Residual of the complementarity slackness for the upper bound constraints.
Definition ISQP.C:854
autoPtr< scalarField > deltaUTilda_
Definition ISQP.H:169
autoPtr< OFstream > meritFunctionFile_
File including the l1 merit function.
Definition ISQP.H:232
autoPtr< scalarField > deltaLs_
Definition ISQP.H:168
void zeroUpdates()
Zero the updates computed in the previous optimisation cycle.
Definition ISQP.C:193
virtual scalar computeMeritFunction()
Compute merit function. Could be different than the objective in the presence of constraints.
Definition ISQP.C:1133
tmp< scalarField > ShermanMorrisonRank1Update(const PtrList< scalarField > &r1Updates, const scalarField &p, const refPtr< scalarField > &diag, const scalarField &mult, label n)
Recursive (naive) implementation of the rank-1 update.
Definition ISQP.C:459
tmp< scalarField > resFGs()
Residual of the inequality constraint slackness.
Definition ISQP.C:784
scalar eps_
Infinitesimal quantity.
Definition ISQP.H:180
scalarField deltaP_
Definition ISQP.H:164
preconditioners
Enumeration of preconditioners for the subproblem.
Definition ISQP.H:76
@ diag
Definition ISQP.H:77
@ invHessian
Definition ISQP.H:77
@ ShermanMorrison
Definition ISQP.H:77
autoPtr< scalarField > deltaExtraVars_
Definition ISQP.H:171
label maxDxIters_
Maxmimum number of iterations for the deltaX equation.
Definition ISQP.H:202
void allocateBoundMultipliers()
Allocate multipliers for the bound constraints.
Definition ISQP.C:91
tmp< scalarField > resFus()
Residual of the upper bounds slackness.
Definition ISQP.C:830
label maxLineSearchIters_
Maxmimum number of line search iterations for each iteration of the subproblem.
Definition ISQP.H:197
autoPtr< scalarField > deltaZ_
Definition ISQP.H:172
scalar relTol_
Relative tolerance of the CG solver, solving for deltaP_.
Definition ISQP.H:207
autoPtr< scalarField > extraVars_
Extra variables for finding solutions even in infeasible problems.
Definition ISQP.H:148
autoPtr< scalarField > ls_
Slack variables the lower bound constraints.
Definition ISQP.H:132
tmp< scalarField > diagPreconditioner()
Diagonal preconditioner of the deltaP eqn.
Definition ISQP.C:389
tmp< scalarField > resFlamda()
Residual of the complementarity slackness for the inequality constraints.
Definition ISQP.C:810
bool includeBoundConstraints_
Are bound constraints included?
Definition ISQP.H:98
bool disableDamping_
Disable damping.
Definition ISQP.H:227
autoPtr< scalarField > lTilda_
Lagrange multipliers for the lower bound constraints.
Definition ISQP.H:127
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-Newton methods coupled with line search.
Definition ISQP.C:1152
void updateYS()
Update the vectors accosiated with the Hessian matrix.
Definition ISQP.C:136
tmp< scalarField > resFlTilda()
Residual of the complementarity slackness for the lower bound constraints.
Definition ISQP.C:844
autoPtr< scalarField > uTilda_
Lagrange multipliers for the upper bound constraints.
Definition ISQP.H:137
void storeOldFields()
Store old fields needed for the next iter.
Definition ISQP.C:976
tmp< scalarField > computeResiduals()
Compute and return residuals based on the current solution.
Definition ISQP.C:667
void updateSizes()
Update sizes of fields related to the active design variables.
Definition ISQP.C:58
void solveDeltaPEqn()
Solve the equation for deltaX, which is the expensive part of the Newtopn step.
Definition ISQP.C:215
autoPtr< scalarField > deltaLTilda_
Definition ISQP.H:167
scalarField gs_
Slack variables for the inequality constraints.
Definition ISQP.H:122
label maxNewtonIters_
Maxmimum number of Newton iterations for the subproblem.
Definition ISQP.H:191
label preconType_
Which preconditioner to use for the solution of the subproblem.
Definition ISQP.H:212
tmp< scalarField > ShermanMorrisonPrecon(const scalarField &vector)
Use the Sherman-Morrison formula to compute the matrix-free product of the approximate inverse of the...
Definition ISQP.C:419
void allocateLagrangeMultipliers()
Allocate Lagrange multipliers for the inequality constraints.
Definition ISQP.C:117
tmp< scalarField > computeRHSForDeltaX(const scalarField &FDx)
Compute the RHS for the deltaX equation.
Definition ISQP.C:247
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 updateVectors(const scalarField &derivatives, const scalarField &derivativesOld)
Update y and s vectors.
Definition LBFGS.C:91
virtual tmp< scalarField > HessianVectorProduct(const scalarField &vector)
Compute the Hessian - vector product.
Definition LBFGS.C:242
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
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition List.H:469
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
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
Base class for Sequantial Quadratic Programming (SQP) methods.
Definition SQPBase.H:53
scalarField lamdas_
Lagrange multipliers.
Definition SQPBase.H:71
virtual bool addToFile(Ostream &os) const
Write continuation info.
Definition SQPBase.C:102
scalar mu_
Penalty value for the merit function.
Definition SQPBase.H:86
virtual bool writeMeritFunction(const updateMethod &UpdateMethod)
Write info about the merit function.
Definition SQPBase.C:115
List< scalarField > constraintDerivativesOld_
The previous constraint derivatives.
Definition SQPBase.H:66
scalar delta_
Safety factor.
Definition SQPBase.H:91
scalarField LagrangianDerivatives_
Derivatives of the Lagrangian function.
Definition SQPBase.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition autoPtrI.H:37
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Abstract base class for optimisation methods supporting constraints. Does not add functionality to up...
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...
const Type & value() const noexcept
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for line search methods.
Definition lineSearch.H:54
scalarField correctionOld_
The previous correction.
Definition quasiNewton.H:88
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
void reset(T *p=nullptr) noexcept
Delete managed pointer and set to new given pointer.
Definition refPtrI.H:314
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.
label nConstraints_
Number of constraints.
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.
PtrList< scalarField > constraintDerivatives_
Derivatives of the constraints.
const labelList & activeDesignVars_
Map to active design variables.
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
virtual void updateOldCorrection(const scalarField &oldCorrection)
Update old correction. Useful for quasi-newton methods coupled with line search.
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
volScalarField & p
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
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
Type gMin(const FieldField< Field, Type > &f)
Vector< scalar > vector
Definition vector.H:57
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
dictionary dict
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const scalar xMin
const scalar xMax
const vector L(dict.get< vector >("L"))