Loading...
Searching...
No Matches
MMA.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 "MMA.H"
31#include "scalarField.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
39 (
41 MMA,
43 );
45 (
47 MMA,
49 );
50}
51
52
53// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54
56{
57 const label n = activeDesignVars_.size();
58 if (n != xNew_.size())
59 {
60 x0_.setSize(n, Zero);
65 a_.setSize(n, Zero);
66 b_.setSize(n, Zero);
69 deltaX_.setSize(n, Zero);
70 deltaEta_.setSize(n, Zero);
71 deltaKsi_.setSize(n, Zero);
72 }
73}
74
75
77{
78 rho_.setSize(cValues_.size() + 1, raa0_);
79
80 // Store old objective and constraint values for GCMMA
81 oldObjectiveValue_ = objectiveValue_;
82 oldCValues_ = cValues_;
83
84 if (variableRho_)
85 {
86 const scalarField x(designVars_().getVars(), activeDesignVars_);
87 const scalarField xMin
88 (designVars_().lowerBounds()(), activeDesignVars_);
89 const scalarField xMax
90 (designVars_().upperBounds()(), activeDesignVars_);
91 const scalarField span(xMax - xMin);
92 const scalarField activeObjDerivs
93 (objectiveDerivatives_, activeDesignVars_);
94
95 // Objective r
96 rho_[0] =
97 maxInitRhoMult_
98 /globalSum(x.size())*globalSum(mag(activeObjDerivs)*span);
99
100 // Constraints r
101 forAll(constraintDerivatives_, i)
102 {
103 const scalarField activeDerivs
104 (
105 constraintDerivatives_[i], activeDesignVars_
106 );
107 rho_[i + 1] =
108 maxInitRhoMult_
109 /globalSum(x.size())*globalSum(mag(activeDerivs)*span);
110
111 }
112
113 // Find a rho value that produces an approximate objective/constraint
114 // combination, evaluated at the new point, that is lower than the
115 // current one
116 if (dynamicRhoInitialisation_)
117 {
118 Info<< "-----------------------------------------" << endl;
119 Info<< "Solving sub problems for initializing rho" << endl;
120 Info<< "-----------------------------------------" << endl;
121 // Backup bounds and initialisation, to be restored after the rho
122 // initialisation phase
123 scalarField lowerBck = lower_;
124 scalarField upperBck = upper_;
125 scalarField aBck = a_;
126 scalarField bBck = b_;
127 scalarField x0Bck = x0_;
128 scalarField x00Bck = x00_;
129
130 // First, do a peudo-update of the the design variables
131 designVars_().storeDesignVariables();
132 x0_.map(designVars_().getVars(), activeDesignVars_);
133
134 Info<< "Initial pseudo update of the design variables" << endl;
135 updateBounds();
136 initialize();
137 solveSubproblem();
138
139 designVars_().update(correction_);
140
141 // Check if the approximate function is higher than the CURRENT
142 // function and, if not, update the rho value until this happens.
143 // Used as an inexpensive way to obtain a decent rho initialisation
144 label iter(0);
145 while (!converged() && iter++ < 10)
146 {
147 Info<< nl << "Dynamic rho initialisation iteration " << iter
148 << nl << endl;
149 designVars_().resetDesignVariables();
150 updateRho();
151 solveSubproblem();
152 designVars_().update(correction_);
153 }
154
155 Info<< "-----------------------------------------" << endl;
156 Info<< "Dynamic rho initialisation converged in " << iter
157 << " iterations " << endl;
158 Info<< "-----------------------------------------" << endl;
159
160 // Restore bounds and design variables
161 lower_ = lowerBck;
162 upper_ = upperBck;
163 a_ = aBck;
164 b_ = bBck;
165 designVars_().resetDesignVariables();
166 x0_ = x0Bck;
167 x00_ = x00Bck;
168 }
169
170 rho_ *= dynamicRhoMult_;
171
172 // Bound too small values
173 if (boundRho_)
174 {
175 rho_ = max(rho_, scalar(1.e-6));
178 << "Computed r values " << rho_ << endl;
179 }
180}
181
182
184{
185 const scalarField x(designVars_().getVars(), activeDesignVars_);
186 const scalarField xMin(designVars_().lowerBoundsRef(), activeDesignVars_);
187 const scalarField xMax(designVars_().upperBoundsRef(), activeDesignVars_);
188 const scalarField span(xMax - xMin);
189
190 if (counter_ > 1)
191 {
192 forAll(x, i)
193 {
194 scalar gamma
195 (
196 (x[i] - x0_[i])*(x0_[i] - x00_[i]) < scalar(0)
197 ? asymDecr_ : asymIncr_
198 );
199 lower_[i] = x[i] - gamma*(x0_[i] - lower_[i]);
200 upper_[i] = x[i] + gamma*(upper_[i] - x0_[i]);
201 }
202
203 lower_ = min(lower_, x - 0.01*(xMax - xMin));
204 lower_ = max(lower_, x - 10.0*(xMax - xMin));
205
206 upper_ = max(upper_, x + 0.01*(xMax - xMin));
207 upper_ = min(upper_, x + 10.0*(xMax - xMin));
208 }
209 else
210 {
211 lower_ = x - sInit_*span;
212 upper_ = x + sInit_*span;
213 }
214
215 a_ = max(xMin, lower_ + 0.1*(x - lower_));
216 a_ = max(a_, x - move_*span);
217 b_ = min(xMax, upper_ - 0.1*(upper_ - x));
218 b_ = min(b_, x + move_*span);
219}
220
221
223{
224 const label m(cValues_.size());
225 // Allocate correct sizes for array depending on the number of constraints
226 if (!c_.good())
227 {
228 alpha_.setSize(m, Zero);
229 if (coeffsDict(typeName).found("c"))
230 {
231 Info<< "Reading constraint penalty function type from dict" << endl;
232 c_.reset(Function1<scalar>::New("c", coeffsDict(typeName)));
233 }
234 else
235 {
236 Info<< "Setting constant penalty factor" << endl;
237 c_.reset(new Function1Types::Constant<scalar>("c", 100));
238 }
239 d_.setSize(m, coeffsDict(typeName).getOrDefault<scalar>("d", 1));
240 deltaLamda_.setSize(m, Zero);
241 deltaY_.setSize(m, Zero);
242 deltaS_.setSize(m, Zero);
243 deltaMu_.setSize(m, Zero);
244 }
245 const scalar t = mesh_.time().timeOutputValue();
246 Info<< "Penalty multiplier (c):: " << c_->value(t) << endl;
247
248 // Scalar quantities
249 eps_ = 1.;
250 z_ = 1.;
251 zeta_ = 1.;
252
253 // Quantities related to design variables
254 xNew_ = 0.5*(a_ + b_);
255 ksi_ = max(scalar(1), 1./(xNew_ - a_));
256 Eta_ = max(scalar(1), 1./(b_ - xNew_));
257
258 // Quantities related to constraints
259 y_.setSize(m, scalar(1));
260 lamda_.setSize(m, scalar(1));
261 s_.setSize(m, scalar(1));
262 mu_.setSize(m, Zero);
263 mu_ = max(scalar(1), 0.5*c_->value(t));
264}
265
266
268{
269 // Store old solutions
270 if (counter_ > 0)
272 x00_ = x0_;
273 }
274 x0_.map(designVars_().getVars(), activeDesignVars_);
275}
276
277
279(
280 const scalarField& derivs,
281 const scalar r,
282 const scalarField& x
283)
284{
285 const scalarField lowerBounds
286 (
287 designVars_().lowerBoundsRef(), activeDesignVars_
288 );
289 const scalarField upperBounds
290 (
291 designVars_().upperBoundsRef(), activeDesignVars_
292 );
293 tmp<scalarField> tres(tmp<scalarField>::New(x.size(), Zero));
294 scalarField& res = tres.ref();
295
296 res =
297 sqr(upper_ - x)
298 *(
299 1.001*max(derivs, scalar(0))
300 + 0.001*max(-derivs, scalar(0))
301 + r/(upperBounds - lowerBounds)
302 );
303
304 return tres;
305}
306
308(
309 const scalarField& derivs,
310 const scalar r
312{
313 return
314 p(derivs, r, scalarField(designVars_().getVars(), activeDesignVars_));
315}
316
317
319(
320 const scalarField& derivs,
321 const scalar r,
322 const scalarField& x
323)
324{
325 const scalarField lowerBounds
326 (
327 designVars_().lowerBoundsRef(), activeDesignVars_
328 );
329 const scalarField upperBounds
330 (
331 designVars_().upperBoundsRef(), activeDesignVars_
332 );
334 scalarField& res = tres.ref();
335
336 res =
337 sqr(x - lower_)
338 *(
339 0.001*max(derivs, scalar(0))
340 + 1.001*max(-derivs, scalar(0))
341 + r/(upperBounds - lowerBounds)
342 );
343
344 return tres;
345}
346
348(
349 const scalarField& derivs,
350 const scalar r
352{
353 return
354 q(derivs, r, scalarField(designVars_().getVars(), activeDesignVars_));
355}
356
357
359(
360 const scalarField& vars
361)
362{
363 tmp<scalarField> tres(tmp<scalarField>::New(cValues_.size(), Zero));
364 scalarField& res = tres.ref();
365 forAll(res, i)
366 {
367 const scalarField activeDerivs
368 (
369 constraintDerivatives_[i], activeDesignVars_
370 );
371 tmp<scalarField> pI(p(activeDerivs, rho_[i + 1]));
372 tmp<scalarField> qI(q(activeDerivs, rho_[i + 1]));
373 res[i] = globalSum(pI/(upper_ - vars) + qI/(vars - lower_));
374 }
375
376 return tres;
377}
378
379
381{
382 const scalarField vars(designVars_().getVars(), activeDesignVars_);
383 tmp<scalarField> tres( - cValues_ + gConstr(vars));
384
385 return tres;
386}
387
388
390{
391 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
392 tmp<scalarField> tres(p(activeObjDerivs, rho_[0]));
393 scalarField& res = tres.ref();
394 forAll(constraintDerivatives_, cI)
395 {
396 const scalarField activeDerivs
397 (
398 constraintDerivatives_[cI], activeDesignVars_
399 );
400 res += lamda_[cI]*p(activeDerivs, rho_[cI + 1]);
401 }
402
403 return tres;
404}
405
406
408{
409 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
410 tmp<scalarField> tres(q(activeObjDerivs, rho_[0]));
411 scalarField& res = tres.ref();
412 forAll(constraintDerivatives_, cI)
413 {
414 const scalarField activeDerivs
415 (
416 constraintDerivatives_[cI], activeDesignVars_
417 );
418 res += lamda_[cI]*q(activeDerivs, rho_[cI + 1]);
419 }
420
421 return tres;
422}
423
424
426{
427 deltaLamda_ = Zero;
428 deltaZ_ = Zero;
429 deltaX_ = Zero;
430 deltaY_ = Zero;
431 deltaS_ = Zero;
433 deltaMu_ = Zero;
434 deltaEta_ = Zero;
435 deltaKsi_ = Zero;
436}
437
438
439Foam::scalar Foam::MMA::fTilda
440(
441 const scalarField& xInit,
442 const scalarField& x,
443 const scalar f,
444 const scalarField& dfdx,
445 const scalar rho
446)
447{
448 return scalar
449 (
450 f
451 + globalSum
452 (
453 p(dfdx, rho, xInit)*(1./(upper_ - x) - 1./(upper_ - xInit))
454 + q(dfdx, rho, xInit)*(1./(x - lower_) - 1./(xInit - lower_))
455 )
456 );
457}
458
459
461(
462 const word& name,
463 const label size,
464 const scalar value
465)
466{
467 tmp<scalarField> tfield(tmp<scalarField>::New(size, value));
468 scalarField& field = tfield.ref();
469 if (found(name))
471 field = scalarField(name, *this, field.size());
472 }
473 return tfield;
474}
475
476
478(
480 const word& name,
481 const label size,
482 const scalarField& defaultField
483)
484{
485 if (found(name))
486 {
487 field = scalarField(name, *this, size);
488 }
489 else
490 {
491 field = defaultField;
492 }
493}
494
495
497{
498 // Zero the updates computed in the previous optimisation cycle
499 zeroUpdates();
500
501 // Easy access to number of design variables and constraints
502 const label n(xNew_.size());
503 const label m(cValues_.size());
504 const scalarField x(designVars_().getVars(), activeDesignVars_);
505
506 // Fields needed to compute Psi and its derivatives
507 tmp<scalarField> tpL(pLamda());
508 const scalarField& pL = tpL();
509 tmp<scalarField> tqL(qLamda());
510 const scalarField& qL = tqL();
511
512 // Build parts of the LHS matrices
513 scalarField DYInv(scalar(1.)/(d_ + mu_/y_));
514 scalarField DLamdaY((s_/lamda_) + DYInv);
515 scalarField DxInv
516 (
517 scalar(1.)/
518 (
519 scalar(2)*pL/pow3(upper_ - xNew_)
520 + scalar(2)*qL/pow3(xNew_ - lower_)
521 + ksi_/(xNew_ - a_)
522 + Eta_/(b_ - xNew_)
523 )
524 );
525
526 // Compute G
528 for (label i = 0; i < m; ++i)
529 {
530 const scalarField activeDerivs
531 (
532 constraintDerivatives_[i], activeDesignVars_
533 );
534 tmp<scalarField> tp = p(activeDerivs, rho_[i + 1]);
535 tmp<scalarField> tq = q(activeDerivs, rho_[i + 1]);
536 scalarField row(tp/sqr(upper_ - xNew_) - tq/sqr(xNew_ - lower_));
537 forAll(row, j)
538 {
539 G[i][j] = row[j];
540 }
541 }
542
543 // Compute parts of the RHS
544 scalarField deltaXTilda
545 (
546 pL/sqr(upper_ - xNew_) - qL/sqr(xNew_ - lower_)
547 - eps_/(xNew_ - a_)
548 + eps_/(b_ - xNew_)
549 );
550 const scalar t = mesh_.time().timeOutputValue();
551 scalarField deltaYTilda
552 (
553 c_->value(t) + d_*y_ - lamda_ - eps_/y_
554 );
555 scalar deltaZTilda = alpha0_ - sum(lamda_*alpha_) - eps_/z_;
556 scalarField deltaLamdaYTilda
557 (
558 // part of deltaLamdaTilda
559 gConstr(xNew_)
560 - alpha_*z_
561 - y_
562 - b()
563 + eps_/lamda_
564 // part of deltaYTilda
565 + DYInv*deltaYTilda
566 );
567
568
569 // Solve system for deltaLamda and deltaZ
570 // The deltaZ equation is substituted into
571 // the one of deltaLamda
572 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
573 scalarSquareMatrix lhs(m, m, Zero);
574 scalar zRatio(z_/zeta_);
575 scalarField rhs(m, Zero);
576
577 // Deal with parts of lhs and rhs that might need global reduction first
578 for (label i = 0; i < m; ++i)
579 {
580 // Prepare LHS
581 for (label j = 0; j < m; ++j)
582 {
583 for (label l = 0; l < n; ++l)
584 {
585 // Part from GDxInvGT
586 lhs[i][j] += DxInv[l]*G[i][l]*G[j][l];
587 }
588 }
589 // Prepare RHS
590 for (label l = 0; l < n; ++l)
591 {
592 // Part from GDxInvDeltaXTilda
593 rhs[i] -= DxInv[l]*G[i][l]*deltaXTilda[l];
594 }
595 }
596 // If design variables are distributed to many processors, reduce
597 if (globalSum_)
598 {
601 }
602
603 // Add remaining parts from deltaLamdaYTilda and the deltaZ eqn
604 rhs += deltaLamdaYTilda + alpha_*deltaZTilda*zRatio;
605 for (label i = 0; i < m; ++i)
606 {
607 // Prepare LHS
608 for (label j = 0; j < m; ++j)
609 {
610 // Part from delta z
611 lhs[i][j] += alpha_[i]*alpha_[j]*zRatio;
612 }
613 // Part from DLamdaY
614 lhs[i][i] += DLamdaY[i];
615 }
616 solve(deltaLamda_, lhs, rhs);
617 deltaZ_ = (sum(alpha_*deltaLamda_) - deltaZTilda)*zRatio;
618
619 // Compute the rest of the corrections using backwards substitution
620
621 // deltaX
622 // No need for reduction since DxInv is a diagonal matrix
623 deltaX_ = - DxInv*(G.Tmul(deltaLamda_) + deltaXTilda);
624
625 // deltaY
626 deltaY_ = DYInv*(deltaLamda_ - deltaYTilda);
627
628 // deltaS
629 deltaS_ = (-s_*deltaLamda_ + eps_)/lamda_ - s_;
630
631 // deltaZeta
632 deltaZeta_ = -zeta_/z_*deltaZ_ - zeta_ + eps_/z_;
633
634 // deltaMu
635 deltaMu_ = (-mu_*deltaY_ + eps_)/y_ - mu_;
636
637 // deltaEta
639
640 // deltaKsi
641 deltaKsi_ = (-ksi_*deltaX_ + eps_)/(xNew_ - a_) - ksi_;
642}
643
644
645Foam::scalar Foam::MMA::lineSearch()
646{
647 const label n(xNew_.size());
648 const label m(cValues_.size());
649 scalar step(1.);
650
651 // Perform bound checks and adjust step accordingly
652 for (label i = 0; i < n; ++i)
653 {
654 if
655 (
656 xNew_[i] + step*deltaX_[i] - a_[i] - 0.01*(xNew_[i] - a_[i])
657 < scalar(0)
658 )
659 {
660 step = -0.99*(xNew_[i] - a_[i])/deltaX_[i];
661 }
662
663 if
664 (
665 -xNew_[i] - step*deltaX_[i] + b_[i] - 0.01*(b_[i] - xNew_[i])
666 < scalar(0)
667 )
668 {
669 step = 0.99*(b_[i] - xNew_[i])/deltaX_[i];
670 }
671
672 adjustStep(step, ksi_[i], deltaKsi_[i]);
673 adjustStep(step, Eta_[i], deltaEta_[i]);
674 }
675
676 for (label i = 0; i < m; ++i)
677 {
678 adjustStep(step, y_[i], deltaY_[i]);
679 adjustStep(step, lamda_[i], deltaLamda_[i]);
680 adjustStep(step, mu_[i], deltaMu_[i]);
681 adjustStep(step, s_[i], deltaS_[i]);
682 }
683
684 adjustStep(step, z_, deltaZ_);
685 adjustStep(step, zeta_, deltaZeta_);
686
687 // Each processor might have computed a different step, if design variables
688 // are distributed. Get the global minimum
689 if (globalSum_)
690 {
691 reduce(step, minOp<scalar>());
692 }
693
694 if (debug > 1)
695 {
696 Info<< "Step before line search is " << step << endl;
697 }
698
699 // Old residual
700 scalar normResOld = sqrt(globalSum(magSqr(computeResiduals())));
701 scalar maxRes(GREAT);
702
703 for (label i = 0; i < maxLineSearchIters_ ; ++i)
704 {
705 // Update the solution with given step
706 updateSolution(step);
707
708 // Compute new residuals and their max value
709 scalarField resNew(computeResiduals());
710 scalar normResNew = sqrt(globalSum(magSqr(resNew)));
711 maxRes = gMax(mag(resNew));
712
713 if (normResNew < normResOld)
714 {
716 << "Initial residual = " << normResOld << ", "
717 << "Final residual = " << normResNew << ", "
718 << "No of LineSearch Iterations = " << i + 1
719 << endl;
720 break;
721 }
722 else
723 {
724 // Return solution to previous and reduce step
725 updateSolution(-step);
726 step *= 0.5;
727
728 // If line search could find an appropriate step, increase eps
729 if (i == maxLineSearchIters_ - 1)
730 {
731 eps_ *= 1.5;
732 Info<< "Line search could not find a step that reduced"
733 << " residuals while satisfying the constraints" << nl
734 << "Increasing eps to " << eps_ << endl;
735 }
736 }
737 }
738
739 if (debug > 1)
740 {
741 Info<< "Step after line search is " << step << nl << endl;
742 }
743
744 return maxRes;
745}
746
747
749(
750 scalar& step,
751 const scalar value,
752 const scalar update
753)
754{
755 if (0.99*value + step*update < scalar(0))
756 {
757 step = -0.99*value/update;
758 }
759}
760
761
762void Foam::MMA::updateSolution(const scalar step)
763{
764 xNew_ += step*deltaX_;
765 y_ += step*deltaY_;
766 z_ += step*deltaZ_;
767 lamda_ += step*deltaLamda_;
768 ksi_ += step*deltaKsi_;
770 mu_ += step*deltaMu_;
771 zeta_ += step*deltaZeta_;
772 s_ += step*deltaS_;
773}
774
775
777{
778 const label n(xNew_.size());
779 const label m(cValues_.size());
780 tmp<scalarField> tres(tmp<scalarField>::New(3*n + 4*m + 2, Zero));
781 scalarField& res = tres.ref();
782
783 label iRes(0);
784
785 // dLdx
786 scalarField dPsidx
787 (pLamda()/sqr(upper_ - xNew_) - qLamda()/sqr(xNew_ - lower_));
788 for (label i = 0; i < n; ++i)
789 {
790 res[iRes++] = dPsidx[i] - ksi_[i] + Eta_[i];
791 }
792
793 // dLdy
794 const scalar t = mesh_.time().timeOutputValue();
795 for (label i = 0; i < m; ++i)
796 {
797 res[iRes++] = c_->value(t) + d_[i]*y_[i] - lamda_[i] - mu_[i];
798 }
799
800 // dLdz
801 res[iRes++] = alpha0_ - zeta_ - sum(lamda_*alpha_);
802
803 // Primal feasibility (constraints)
804 scalarField gb(gConstr(xNew_) - b());
805 for (label i = 0; i < m; ++i)
806 {
807 res[iRes++] = gb[i] - alpha_[i]*z_ - y_[i] + s_[i];
808 }
809
810 // ksi slackness
811 for (label i = 0; i < n; ++i)
812 {
813 res[iRes++] = ksi_[i]*(xNew_[i] - a_[i]) - eps_;
814 }
815
816 // Eta slackness
817 for (label i = 0; i < n; ++i)
818 {
819 res[iRes++] = Eta_[i]*(b_[i] - xNew_[i]) - eps_;
820 }
821
822 // y slackness
823 for (label i = 0; i < m; ++i)
824 {
825 res[iRes++] = mu_[i]*y_[i] - eps_;
826 }
827
828 // z slackness
829 res[iRes++] = z_*zeta_ - eps_;
830
831 // Constraint slackness
832 for (label i = 0; i < m; ++i)
833 {
834 res[iRes++] = lamda_[i]*s_[i] - eps_;
835 }
836
837 return tres;
838}
839
840
842{
843 if (normalise_)
844 {
845 // Compute normalisation factors
846 if
847 (
848 !Jnorm_
849 || (continuousNormalisation_ && counter_ < lastNormalisationStep_)
850 )
851 {
852 scalarField activeSens(objectiveDerivatives_, activeDesignVars_);
853 Jnorm_.reset(new scalar(Foam::sqrt(globalSum(sqr(activeSens)))));
854 Cnorm_.reset(new scalarField(cValues_.size(), Zero));
855 scalarField& Cnorm = Cnorm_.ref();
856 forAll(constraintDerivatives_, cI)
857 {
858 scalarField activeConstrSens
859 (constraintDerivatives_[cI], activeDesignVars_);
860 Cnorm[cI] = Foam::sqrt(globalSum(sqr(activeConstrSens)));
861 }
862 Info<< "Computed normalisation factors " << nl
863 << "J: " << Jnorm_() << nl
864 << "C: " << Cnorm_() << endl;
865 }
866
867 // Normalise objective values and gradients
868 objectiveValue_ /= Jnorm_() + SMALL;
869 objectiveDerivatives_ /= Jnorm_() + SMALL;
870
871 // Normalise constraint values and gradients
872 cValues_ *= cw_/(Cnorm_() + SMALL);
873 forAll(constraintDerivatives_, cI)
874 {
875 constraintDerivatives_[cI] *= cw_/(Cnorm_()[cI] + SMALL);
877 }
878}
879
880
881// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
882
883Foam::MMA::MMA
884(
885 const fvMesh& mesh,
886 const dictionary& dict,
887 autoPtr<designVariables>& designVars,
888 const label nConstraints,
889 const word& type
890)
891:
892 constrainedOptimisationMethod(mesh, dict, designVars, nConstraints, type),
893 updateMethod(mesh, dict, designVars, nConstraints, type),
894 x0_(getOrDefaultScalarField("x0", activeDesignVars_.size())),
895 x00_(getOrDefaultScalarField("x00", activeDesignVars_.size())),
896 xNew_(activeDesignVars_.size()),
897 oldObjectiveValue_(getOrDefault<scalar>("J0", Zero)),
898 oldCValues_(getOrDefaultScalarField("C0", nConstraints_)),
899 valsAndApproxs_(2),
900 z_(coeffsDict(type).getOrDefault<scalar>("z", 1.)),
901 alpha0_(coeffsDict(type).getOrDefault<scalar>("alpha0", 1.)),
902 alpha_(0),
903 y_(0),
904 c_(nullptr),
905 d_(0),
906 lower_(xNew_.size(), -GREAT),
907 upper_(xNew_.size(), GREAT),
908 a_(xNew_.size()),
909 b_(xNew_.size()),
910 rho_(0),
911 boundRho_(coeffsDict(type).getOrDefault<bool>("boundRho", false)),
912 correctDVs_(coeffsDict(type).getOrDefault<bool>("correct", true)),
913 lamda_(0),
914 ksi_(xNew_.size()),
915 Eta_(xNew_.size()),
916 mu_(0),
917 zeta_(coeffsDict(type).getOrDefault<scalar>("zeta", 1.)),
918 s_(0),
919 deltaLamda_(0),
920 deltaZ_(Zero),
921 deltaX_(xNew_.size(), Zero),
922 deltaY_(0),
923 deltaS_(0),
924 deltaZeta_(Zero),
925 deltaMu_(0),
926 deltaEta_(xNew_.size(), Zero),
927 deltaKsi_(xNew_.size(), Zero),
928 eps_(1),
929 maxNewtonIters_(coeffsDict(type).getOrDefault<label>("maxIters", 6000)),
930 maxLineSearchIters_
931 (
932 coeffsDict(type).getOrDefault<label>("maxLineSearchIters", 10)
933 ),
934 initializeEverySubproblem_
935 (coeffsDict(type).getOrDefault<bool>("initializeEverySubproblem", false)),
936 asymDecr_(coeffsDict(type).getOrDefault<scalar>("asymptoteDecrease", 0.7)),
937 asymIncr_(coeffsDict(type).getOrDefault<scalar>("asymptoteIncrease", 1.2)),
938 sInit_(coeffsDict(type).getOrDefault<scalar>("sInit", 0.5)),
939 move_(coeffsDict(type).getOrDefault<scalar>("move", 0.5)),
940 raa0_(coeffsDict(type).getOrDefault<scalar>("raa0", 1.e-05)),
941 maxInitRhoMult_(coeffsDict(type).getOrDefault<scalar>("maxInitRhoMult", 0.1)),
942 maxRhoMult_(coeffsDict(type).getOrDefault<scalar>("maxRhoMult", 10)),
943 variableRho_(coeffsDict(type).getOrDefault<bool>("variableRho", false)),
944 dynamicRhoInitialisation_
945 (coeffsDict(type).getOrDefault<bool>("dynamicRhoInitialisation", false)),
946 dynamicRhoMult_
947 (coeffsDict(type).getOrDefault<scalar>("dynamicRhoMult", 0.1)),
948 normalise_(coeffsDict(type).getOrDefault<bool>("normalise", false)),
949 continuousNormalisation_
950 (coeffsDict(type).getOrDefault<bool>("continuousNormalisation", false)),
951 Jnorm_(nullptr),
952 Cnorm_(nullptr),
953 cw_(coeffsDict(type).getOrDefault<scalar>("constraintWeight", 1)),
954 lastNormalisationStep_
955 (
956 coeffsDict(type).getOrDefault<label>("lastNormalisationStep", 20)
957 )
958{
959 // Check that the design variables bounds have been set
960 if (!designVars().lowerBounds() || !designVars().upperBounds())
961 {
963 << "Lower and upper design variable bounds must be set for MMA"
964 << abort(FatalError);
965 }
966
967 // Set initial bounds
969 (
970 lower_,
971 "lower",
972 activeDesignVars_.size(),
973 scalarField(designVars().lowerBoundsRef(), activeDesignVars_)
974 );
976 (
977 upper_,
978 "upper",
979 activeDesignVars_.size(),
980 scalarField(designVars().upperBoundsRef(), activeDesignVars_)
981 );
982}
983
984
985// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
986
988{
989 if (correctDVs_)
990 {
991 // Update sizes of fields related to the active design variables
992 updateSizes();
993
994 // Perform normalisation of the objective and constraint values
995 normalise();
996
997 // Initialize rho values in the p,q functions
998 // These determine how aggressive or conservative the method will be
999 initializeRho();
1000
1001 // Update the bounds associated with the design variables
1002 updateBounds();
1003
1004 // Initiaze variables
1005 initialize();
1006
1007 // Solve subproblem using a Newton optimiser
1008 solveSubproblem();
1009
1010 // Store old solutions, bjective and constaint values
1011 storeOldValues();
1013 // Increase counter
1014 ++counter_;
1015 }
1016}
1017
1018
1020{
1021 // Initialize the solution
1022 if (initializeEverySubproblem_)
1023 {
1024 initialize();
1025 }
1026
1027 scalar resMax(gMax(mag(computeResiduals())));
1028 label iter(0);
1029 do
1030 {
1031 DebugInfo
1032 << "Newton iter " << iter << nl << endl;
1033
1034 // Decrease eps
1035 if (resMax < 0.9*eps_)
1036 {
1037 eps_ *= 0.1;
1038 }
1039
1040 // Computes Newton direction for the subproblem
1041 computeNewtonDirection();
1042
1043 // Perform line search and return max residual of the solution
1044 // satisfying the bound constraints and the residual reduction.
1045 // Upates solution.
1046 resMax = lineSearch();
1047 DebugInfo
1048 << "max residual = " << resMax << ", "
1049 << "eps = " << eps_ << nl << endl;
1050
1051 mesh_.time().printExecutionTime(Info);
1052 } while (iter++ < maxNewtonIters_ && (eps_ > 1.e-07 || resMax > 0.9*eps_));
1053
1054 Info<< "Solved the MMA Newton problem in " << iter << " iterations "
1055 << nl << endl;
1056
1057 DebugInfo
1058 << "Constraint penalization variables (y) " << y_ << endl;
1059
1060 // Pass update to correction field
1061 const scalarField& oldVars = designVars_().getVars();
1062 forAll(activeDesignVars_, avI)
1064 const label vI(activeDesignVars_[avI]);
1065 correction_[vI] = xNew_[avI] - oldVars[vI];
1066 }
1067}
1068
1069
1071{
1072 // Return values, including the objective and constraint values and
1073 // their approximations
1074 label m(cValues_.size());
1075 valsAndApproxs_.set(0, new scalarField(m + 1, Zero));
1076 valsAndApproxs_.set(1, new scalarField(m + 1, Zero));
1077 scalarField& vals = valsAndApproxs_[0];
1078 scalarField& approx = valsAndApproxs_[1];
1079
1080 // Objective value and approximation
1081 const scalarField activeObjDerivs(objectiveDerivatives_, activeDesignVars_);
1082 vals[0] = objectiveValue_;
1083 approx[0] = fTilda(x0_, xNew_, oldObjectiveValue_, activeObjDerivs, rho_[0]);
1084
1085 // Constraint values and approximations
1086 forAll(constraintDerivatives_, i)
1087 {
1088 const scalarField activeDerivs
1089 (
1090 constraintDerivatives_[i], activeDesignVars_
1091 );
1092 vals[i + 1] = cValues_[i];
1093 approx[i + 1] =
1094 fTilda(x0_, xNew_, oldCValues_[i], activeDerivs, rho_[i + 1]);
1095 }
1096}
1097
1098
1101{
1102 return valsAndApproxs_;
1103}
1104
1105
1107{
1108 // Objective/constraint values and approximations
1109 const scalarField& vals = valsAndApproxs_[0];
1110 const scalarField& approx = valsAndApproxs_[1];
1111
1112 const scalarField xMin(designVars_().lowerBounds()(), activeDesignVars_);
1113 const scalarField xMax(designVars_().upperBounds()(), activeDesignVars_);
1114
1115 scalar d = globalSum
1116 (
1117 (upper_ - lower_)*sqr(xNew_ - x0_)
1118 /(upper_ - xNew_)/(xNew_ - lower_)/(xMax - xMin)
1119 );
1120
1121 // Update rho values
1122 forAll(approx, i)
1123 {
1124 const scalar delta = (vals[i] - approx[i])/d;
1125 if (delta > 0)
1126 {
1127 rho_[i] = min(1.1*(rho_[i] + delta), maxRhoMult_*rho_[i]);
1129 }
1130 DebugInfo
1131 << "Updated rho values to " << rho_ << endl;
1132}
1133
1136{
1137 return rho_;
1138}
1139
1141void Foam::MMA::setVariableRho(bool varRho)
1142{
1143 variableRho_ = varRho;
1144}
1145
1146
1148{
1149 updateValuesAndApproximations();
1150 const scalarField& vals = valsAndApproxs_[0];
1151 const scalarField& approx = valsAndApproxs_[1];
1152
1153 bool isConverged(true);
1154 forAll(vals, i)
1155 {
1156 DebugInfo
1157 << nl << "MMA, objective/constraint " << i << nl
1158 << "Approximation " << approx[i]
1159 << " | old value " << vals[i] << nl << endl;
1160 isConverged = isConverged && (approx[i] - vals[i] > 0.);
1161 }
1162
1163 return isConverged;
1164}
1165
1166
1167bool Foam::MMA::writeData(Ostream& os) const
1168{
1169 x0_.writeEntry("x0", os);
1170 x00_.writeEntry("x00", os);
1171 lower_.writeEntry("lower", os);
1172 upper_.writeEntry("upper", os);
1173 os.writeEntry("J0", oldObjectiveValue_);
1174 oldCValues_.writeEntry("C0", os);
1175
1177}
1178
1179
1180// ************************************************************************* //
scalar delta
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
Templated function that returns a constant value.
Definition Constant.H:71
void setSize(label n)
Alias for resize().
Definition List.H:536
Update design variables using the Method of Moving Assymptotes. Can handle inequality constraints.
Definition MMA.H:65
tmp< scalarField > gConstr(const scalarField &vars)
g of all constraint functions
Definition MMA.C:352
scalarField x0_
The previous design variables. Used to compute new bounds.
Definition MMA.H:73
bool normalise_
Perform the normalisation.
Definition MMA.H:319
scalar lineSearch()
Perform line search and return max residual corresponding to the updated solution.
Definition MMA.C:638
void updateSolution(const scalar step)
Update the current solution using the known direction and the given step length.
Definition MMA.C:755
void solveSubproblem()
Solve subproblem using a Newton optimiser.
Definition MMA.C:1012
void computeCorrection()
Compute design variables correction.
Definition MMA.C:980
scalar maxInitRhoMult_
Multiplier of the mean derivative used for the initialisation of rho.
Definition MMA.H:289
bool converged()
Checks convergence of GCMMA.
Definition MMA.C:1140
scalarField deltaX_
Definition MMA.H:213
void initialize()
Allocate fields related to constraints.
Definition MMA.C:215
scalarField deltaEta_
Definition MMA.H:218
scalar zeta_
Lagrange multiplier for the z constraint.
Definition MMA.H:201
PtrList< scalarField > valsAndApproxs_
Objective/Constraint values and approximations in the new point.
Definition MMA.H:107
scalarField mu_
Lagrange multipliers for the y constraints.
Definition MMA.H:196
autoPtr< Function1< scalar > > c_
y multipliers in the objective function
Definition MMA.H:132
scalarField alpha_
Terms multyplying z in the constraint functions.
Definition MMA.H:122
scalar deltaZ_
Definition MMA.H:212
scalar oldObjectiveValue_
Old objective value.
Definition MMA.H:93
void adjustStep(scalar &step, const scalar value, const scalar update)
Adjust step to satisfy cireteria.
Definition MMA.C:742
scalar maxRhoMult_
Multiplier for the max. rho value during its update.
Definition MMA.H:294
virtual bool writeData(Ostream &os) const
Write useful quantities to files.
Definition MMA.C:1160
scalarField deltaLamda_
Definition MMA.H:211
scalarField xNew_
The set of design variables being updated during the subproblem.
Definition MMA.H:86
scalar sInit_
Used to update the assymptotes in the first optimisation cycle.
Definition MMA.H:270
tmp< scalarField > p(const scalarField &derivs, const scalar r, const scalarField &x)
p and q functions, used to approximate the objective and contraints
Definition MMA.C:272
tmp< scalarField > q(const scalarField &derivs, const scalar r, const scalarField &x)
Definition MMA.C:312
bool initializeEverySubproblem_
Initialize every subproblem with x = 0.5*(a + b) and reset Lagrange multipliers.
Definition MMA.H:252
scalarField x00_
The twice previous design variables. Used to compute new bounds.
Definition MMA.H:78
void updateValuesAndApproximations()
Compute objective/constraint values and their approximations.
Definition MMA.C:1063
scalarField deltaS_
Definition MMA.H:215
void updateBounds()
Update the bounds associated with the design variables.
Definition MMA.C:176
bool correctDVs_
Correct the design variables.
Definition MMA.H:174
void computeNewtonDirection()
Compute direction for the Newton problem.
Definition MMA.C:489
scalarField deltaY_
Definition MMA.H:214
scalar move_
Movement of the unatainable upper and lower bounds.
Definition MMA.H:275
scalarField ksi_
Lagrange multipliers for the lower limits constraints.
Definition MMA.H:186
tmp< scalarField > b()
The rhs of the contraint approximations.
Definition MMA.C:373
void initializeRho()
Initialize rho constants in the (p, q) functions.
Definition MMA.C:69
tmp< scalarField > getOrDefaultScalarField(const word &name, const label size, const scalar value=Zero)
Read in scalarField from self (potentially in binary), or allocate field with the size of the active ...
Definition MMA.C:454
scalar z_
Second term in the approximation function.
Definition MMA.H:112
scalar fTilda(const scalarField &xInit, const scalarField &x, const scalar f, const scalarField &dfdx, const scalar rho)
Computation of the approximate function.
Definition MMA.C:433
scalarField deltaMu_
Definition MMA.H:217
bool dynamicRhoInitialisation_
Change rho constants in each iteration?
Definition MMA.H:304
scalarField b_
Upper design variables constraints.
Definition MMA.H:157
bool continuousNormalisation_
Perform the normalisation in each optimisation cycle.
Definition MMA.H:324
void zeroUpdates()
Zero the updates computed in the previous optimisation cycle.
Definition MMA.C:418
tmp< scalarField > pLamda()
p and q functions, using the dual lamda
Definition MMA.C:382
scalar eps_
Infinitesimal quantity.
Definition MMA.H:227
scalarField d_
y^2 multipliers in the objective function
Definition MMA.H:137
scalarField deltaKsi_
Definition MMA.H:219
scalarField lamda_
Constraint Lagrange multipliers.
Definition MMA.H:181
void updateRho()
Update rho value if needed.
Definition MMA.C:1099
label maxLineSearchIters_
Maxmimum number of line search iterations for each iteration of the subproblem.
Definition MMA.H:238
scalarField upper_
Upper design variables asymptotes.
Definition MMA.H:147
scalarField oldCValues_
Old constraint values.
Definition MMA.H:100
scalar asymIncr_
Upper assymprote increase rate.
Definition MMA.H:265
scalarField s_
Slack variables for the inequality constraints.
Definition MMA.H:206
label lastNormalisationStep_
Constaint weight after the normalisation.
Definition MMA.H:344
scalar cw_
Constaint weight after the normalisation.
Definition MMA.H:339
scalar deltaZeta_
Definition MMA.H:216
scalarField a_
Lower design variables constraints.
Definition MMA.H:152
scalar raa0_
Constant in p, q functions.
Definition MMA.H:280
autoPtr< scalar > Jnorm_
Normalisation factor for the objective.
Definition MMA.H:329
scalar asymDecr_
Lower assymprote reduction rate.
Definition MMA.H:260
void setOrDefaultScalarField(scalarField &field, const word &name, const label size, const scalarField &defaultField)
Read in scalarField from self (potentially in binary), or allocate field with the size of the active ...
Definition MMA.C:471
const scalarField & getRho() const
Get rho value if needed.
Definition MMA.C:1128
scalarField lower_
Lower design variables asymptotes.
Definition MMA.H:142
bool variableRho_
Change rho constants in each iteration?
Definition MMA.H:299
tmp< scalarField > qLamda()
Definition MMA.C:400
void normalise()
Normalise the objective and constraints if necessary.
Definition MMA.C:834
scalarField y_
y in the constraints functions
Definition MMA.H:127
void storeOldValues()
Store old objcective and constraint values.
Definition MMA.C:260
bool boundRho_
Bound the rho value with an upper value?
Definition MMA.H:169
scalarField rho_
Constants in the (p,q) functions.
Definition MMA.H:164
scalar dynamicRhoMult_
The rho values obtained by the dynamic rho initialisation might be too conservative....
Definition MMA.H:311
tmp< scalarField > computeResiduals()
Compute and return residuals based on the current solution.
Definition MMA.C:769
void updateSizes()
Update sizes of fields related to the active design variables.
Definition MMA.C:48
tmp< scalarField > Cnorm_
Normalisation factors for the constraints.
Definition MMA.H:334
scalar alpha0_
Term multiplying z in the objective function.
Definition MMA.H:117
void setVariableRho(bool varRho=true)
Set variable rho.
Definition MMA.C:1134
label maxNewtonIters_
Maxmimum number of Newton iterations for the subproblem.
Definition MMA.H:232
scalarField Eta_
Lagrange multipliers for the upper limits constraints.
Definition MMA.H:191
const PtrList< scalarField > & getValuesAndApproximations() const
Get objective/constraint values and their approximations.
Definition MMA.C:1093
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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
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...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for line search methods.
Definition lineSearch.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.
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.
virtual bool writeData(Ostream &os) const
Write continuation data under uniform.
const labelList & activeDesignVars_
Map to active design variables.
label counter_
Optimisation cycle count.
label nConstraints() const
Get the number of constraints.
autoPtr< designVariables > & designVars_
Design variables.
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
rDeltaTY field()
const scalar gamma
Definition EEqn.H:9
mesh update()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define DebugInfo
Report an information message using Foam::Info.
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
RectangularMatrix< scalar > scalarRectangularMatrix
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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
errorManip< error > abort(error &err)
Definition errorManip.H:139
SquareMatrix< scalar > scalarSquareMatrix
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...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
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
labelList f(nPoints)
dictionary dict
CEqn solve()
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const scalar xMin
const scalar xMax