Loading...
Searching...
No Matches
STDMD.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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "STDMD.H"
29#include "EigenMatrix.H"
30#include "QRMatrix.H"
33using namespace Foam::constant::mathematical;
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace DMDModels
40{
43}
44}
45
46const Foam::Enum
47<
48 Foam::DMDModels::STDMD::modeSorterType
49>
50Foam::DMDModels::STDMD::modeSorterTypeNames
51({
52 { modeSorterType::FIRST_SNAPSHOT, "firstSnapshot" },
53 { modeSorterType::KIEWAT , "kiewat" },
54 { modeSorterType::KOU_ZHANG , "kouZhang" }
55});
56
57
58// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59
60Foam::scalar Foam::DMDModels::STDMD::L2norm(const RMatrix& z) const
61{
62 #ifdef FULLDEBUG
63 // Check if the given RectangularMatrix is effectively a column vector
64 if (z.n() != 1)
65 {
67 << "Input matrix is not a column vector."
68 << exit(FatalError);
69 }
70 #endif
71 const bool noSqrt = true;
72 scalar result = z.columnNorm(0, noSqrt);
73 reduce(result, sumOp<scalar>());
74
75 // Heuristic addition to avoid very small or zero norm
76 return max(SMALL, Foam::sqrt(result));
77}
78
79
80Foam::RectangularMatrix<Foam::scalar> Foam::DMDModels::STDMD::orthonormalise
81(
82 RMatrix ez
83) const
84{
85 List<scalar> dz(Q_.n());
86 const label sz0 = ez.m();
87 const label sz1 = dz.size();
88
89 for (label i = 0; i < nGramSchmidt_; ++i)
90 {
91 // dz = Q_ & ez;
92 dz = Zero;
93 for (label i = 0; i < sz0; ++i)
94 {
95 for (label j = 0; j < sz1; ++j)
96 {
97 dz[j] += Q_(i,j)*ez(i,0);
98 }
99 }
100
101 reduce(dz, sumOp<List<scalar>>());
102
103 // ez -= Q_*dz;
104 for (label i = 0; i < sz0; ++i)
105 {
106 scalar t = 0;
107 for (label j = 0; j < sz1; ++j)
108 {
109 t += Q_(i,j)*dz[j];
110 }
111 ez(i,0) -= t;
112 }
113 }
114
115 return ez;
116}
117
118
119void Foam::DMDModels::STDMD::expand(const RMatrix& ez, const scalar ezNorm)
120{
121 Info<< tab << "Expanding orthonormal basis for field: " << fieldName_
122 << endl;
123
124 // Stack a column "(ez/ezNorm)" to "Q"
125 Q_.resize(Q_.m(), Q_.n() + 1);
126 Q_.subColumn(Q_.n() - 1) = ez/ezNorm;
127
128 // Stack a row (Zero) and column (Zero) to "G"
129 G_.resize(G_.m() + 1);
130}
131
132
133void Foam::DMDModels::STDMD::updateG(const RMatrix& z)
134{
135 List<scalar> zTilde(Q_.n(), Zero);
136 const label sz0 = z.m();
137 const label sz1 = zTilde.size();
138
139 // zTilde = Q_ & z;
140 for (label i = 0; i < sz0; ++i)
141 {
142 for (label j = 0; j < sz1; ++j)
143 {
144 zTilde[j] += Q_(i,j)*z(i,0);
145 }
146 }
147
148 reduce(zTilde, sumOp<List<scalar>>());
149
150 // G_ += SMatrix(zTilde^zTilde);
151 for (label i = 0; i < G_.m(); ++i)
152 {
153 for (label j = 0; j < G_.n(); ++j)
154 {
155 G_(i, j) += zTilde[i]*zTilde[j];
156 }
157 }
158}
159
160
161void Foam::DMDModels::STDMD::compress()
162{
163 Info<< tab << "Compressing orthonormal basis for field: " << fieldName_
164 << endl;
165
166 RMatrix q(1, 1, Zero);
167
168 if (Pstream::master())
169 {
170 const bool symmetric = true;
171 const EigenMatrix<scalar> EM(G_, symmetric);
172 const SquareMatrix<scalar>& EVecs = EM.EVecs();
173 DiagonalMatrix<scalar> EVals(EM.EValsRe());
174
175 // Sort eigenvalues in descending order, and track indices
176 const auto descend = [](scalar a, scalar b){ return a > b; };
177 const labelList permutation(EVals.sortPermutation(descend));
178 EVals.applyPermutation(permutation);
179 EVals.resize(EVals.size() - 1);
180
181 // Update "G"
182 G_ = SMatrix(maxRank_, Zero);
183 G_.diag(EVals);
184
185 q.resize(Q_.n(), maxRank_);
186 for (label i = 0; i < maxRank_; ++i)
187 {
188 q.subColumn(i) = EVecs.subColumn(permutation[i]);
189 }
190 }
193
194 // Update "Q"
195 Q_ = Q_*q;
196}
197
198
199Foam::SquareMatrix<Foam::scalar> Foam::DMDModels::STDMD::
200reducedKoopmanOperator()
201{
202 Info<< tab << "Computing Atilde" << endl;
203
204 Info<< tab << "Computing local Rx" << endl;
205
206 RMatrix Rx(1, 1, Zero);
207 {
208 QRMatrix<RMatrix> QRM
209 (
210 Qupper_,
213 );
214 RMatrix& R = QRM.R();
215 Rx.transfer(R);
216 }
217 Rx.round();
218
219 // Convenience objects A1, A2, A3 to compute "Atilde"
220 RMatrix A1(1, 1, Zero);
221
222 if (Pstream::parRun())
223 {
224 // Parallel direct tall-skinny QR decomposition
225 // (BGD:Fig. 5) & (DGHL:Fig. 2)
226 // Tests revealed that the distribution of "Q" does not affect
227 // the final outcome of TSQR decomposition up to sign
228
229 // Don't clear storage on persistent buffer
230 PstreamBuffers pBufs;
231 pBufs.allowClearRecv(false);
232
233 const label myProcNo = Pstream::myProcNo();
234 const label procNoInSubset = myProcNo % nAgglomerationProcs_;
235
236 // Send Rx from sender subset neighbours to receiver subset master
237 if (procNoInSubset != 0)
238 {
239 const label procNoSubsetMaster = myProcNo - procNoInSubset;
240
241 UOPstream toSubsetMaster(procNoSubsetMaster, pBufs);
242 toSubsetMaster << Rx;
243
244 Rx.clear();
245 }
246
247 pBufs.finishedSends();
248
249 // Receive Rx by receiver subset masters
250 if (procNoInSubset == 0)
251 {
252 // Accept only subset masters possessing sender subset neighbours
253 if (myProcNo + 1 < Pstream::nProcs())
254 {
255 for
256 (
257 label nbr = myProcNo + 1;
258 (
259 nbr < myProcNo + nAgglomerationProcs_
260 && nbr < Pstream::nProcs()
261 );
262 ++nbr
263 )
264 {
265 RMatrix recvMtrx;
266
267 UIPstream fromNbr(nbr, pBufs);
268 fromNbr >> recvMtrx;
269
270 // Append received Rx to Rx of receiver subset masters
271 if (recvMtrx.size() > 0)
272 {
273 Rx.resize(Rx.m() + recvMtrx.m(), Rx.n());
274 Rx.subMatrix
275 (
276 Rx.m() - recvMtrx.m(),
277 Rx.n() - recvMtrx.n()
278 ) = recvMtrx;
279 }
280 }
281 }
282
283 {
284 // Apply interim QR decomposition
285 // on Rx of receiver subset masters
286 QRMatrix<RMatrix> QRM
287 (
291 Rx
292 );
293 RMatrix& R = QRM.R();
294 Rx.transfer(R);
295 }
296 Rx.round();
297 }
298
299 pBufs.clear();
300
301 // Send interim Rx from subset masters to the master
302 if (procNoInSubset == 0)
303 {
304 if (!Pstream::master())
305 {
306 UOPstream toMaster(Pstream::masterNo(), pBufs);
307 toMaster << Rx;
308 }
309 }
310
311 pBufs.finishedSends();
312
313
314 // Receive interim Rx by the master, and apply final operations
315 if (Pstream::master())
316 {
317 for
318 (
319 label nbr = Pstream::masterNo() + nAgglomerationProcs_;
320 nbr < Pstream::nProcs();
321 nbr += nAgglomerationProcs_
322 )
323 {
324 RMatrix recvMtrx;
325
326 UIPstream fromSubsetMaster(nbr, pBufs);
327 fromSubsetMaster >> recvMtrx;
328
329 // Append received Rx to Rx of the master
330 if (recvMtrx.size() > 0)
331 {
332 Rx.resize(Rx.m() + recvMtrx.m(), Rx.n());
333 Rx.subMatrix
334 (
335 Rx.m() - recvMtrx.m(),
336 Rx.n() - recvMtrx.n()
337 ) = recvMtrx;
338 }
339 }
340
341 Info<< tab << "Computing TSQR" << endl;
342 {
343 QRMatrix<RMatrix> QRM
344 (
348 Rx
349 );
350 RMatrix& R = QRM.R();
351 Rx.transfer(R);
352 }
353 Rx.round();
354
355 Info<< tab << "Computing RxInv" << endl;
356 RxInv_ = MatrixTools::pinv(Rx);
357
358 Info<< tab << "Computing A1" << endl;
359 A1 = RxInv_*RMatrix(MatrixTools::pinv(Rx*(G_^Rx)));
360 Rx.clear();
361 }
362 Pstream::broadcast(RxInv_);
364
365 Info<< tab << "Computing A2" << endl;
366 SMatrix A2(Qupper_ & Qlower_);
367 reduce(A2, sumOp<SMatrix>());
368 Qlower_.clear();
369
370 Info<< tab << "Computing A3" << endl;
371 SMatrix A3(Qupper_ & Qupper_);
372 reduce(A3, sumOp<SMatrix>());
373
374 Info<< tab << "Computing Atilde" << endl;
375 // by optimized matrix chain multiplication
376 // obtained by dynamic programming memoisation
377 return SMatrix(RxInv_ & ((A2*(G_*A3))*A1));
378 }
379 else
380 {
381 Info<< tab << "Computing RxInv" << endl;
382 RxInv_ = MatrixTools::pinv(Rx);
383
384 Info<< tab << "Computing A1" << endl;
385 A1 = RxInv_*RMatrix(MatrixTools::pinv(Rx*(G_^Rx)));
386 Rx.clear();
387
388 Info<< tab << "Computing A2" << endl;
389 SMatrix A2(Qupper_ & Qlower_);
390 Qlower_.clear();
391
392 Info<< tab << "Computing A3" << endl;
393 SMatrix A3(Qupper_ & Qupper_);
394
395 Info<< tab << "Computing Atilde" << endl;
396 return SMatrix(RxInv_ & ((A2*(G_*A3))*A1));
397 }
398}
399
400
401bool Foam::DMDModels::STDMD::eigendecomposition(SMatrix& Atilde)
402{
403 bool fail = false;
404
405 // Compute eigenvalues, and clip eigenvalues with values < "minEval"
406 if (Pstream::master())
407 {
408 Info<< tab << "Computing eigendecomposition" << endl;
409
410 // Replace "Atilde" by eigenvalues (in-place eigendecomposition)
411 const EigenMatrix<scalar> EM(Atilde);
412 const DiagonalMatrix<scalar>& evalsRe = EM.EValsRe();
413 const DiagonalMatrix<scalar>& evalsIm = EM.EValsIm();
414
415 evals_.resize(evalsRe.size());
416 evecs_ = RCMatrix(EM.complexEVecs());
417
418 // Convert scalar eigenvalue pairs to complex eigenvalues
419 label i = 0;
420 for (auto& eval : evals_)
421 {
422 eval = complex(evalsRe[i], evalsIm[i]);
423 ++i;
424 }
425
426 Info<< tab << "Filtering eigenvalues" << endl;
427
428 List<complex> cp(evals_.size());
429 auto it =
430 std::copy_if
431 (
432 evals_.cbegin(),
433 evals_.cend(),
434 cp.begin(),
435 [&](const complex& x){ return mag(x) > minEval_; }
436 );
437 cp.resize(std::distance(cp.begin(), it));
438
439 if (cp.size() == 0)
440 {
442 << "No eigenvalue with mag(eigenvalue) larger than "
443 << "minEval = " << minEval_ << " was found." << nl
444 << " Input field may contain only zero-value elements," << nl
445 << " e.g. no-slip velocity condition on a given patch." << nl
446 << " Otherwise, please decrease the value of minEval." << nl
447 << " Skipping further dynamics/mode computations." << nl
448 << endl;
449
450 fail = true;
451 }
452 else
453 {
454 evals_ = cp;
455 }
456 }
457 Pstream::broadcast(fail);
458
459 if (fail)
460 {
461 return false;
462 }
463
464 Pstream::broadcast(evals_);
465 Pstream::broadcast(evecs_);
466
467 return true;
468}
469
470
471void Foam::DMDModels::STDMD::frequencies()
472{
473 if (Pstream::master())
474 {
475 Info<< tab << "Computing frequencies" << endl;
476
477 freqs_.resize(evals_.size());
478
479 // Frequency equation (K:Eq. 81)
480 auto frequencyEquation =
481 [&](const complex& eval)
482 {
483 return Foam::log(max(eval, complex(SMALL))).imag()/(twoPi*dt_);
484 };
485
486 // Compute frequencies
487 std::transform
488 (
489 evals_.cbegin(),
490 evals_.cend(),
491 freqs_.begin(),
492 frequencyEquation
493 );
494
495 Info<< tab << "Computing frequency indices" << endl;
496
497 // Initialise iterator by the first search
498 auto margin = [&](const scalar& x){ return (fMin_ <= x && x < fMax_); };
499
500 auto it = std::find_if(freqs_.cbegin(), freqs_.cend(), margin);
501
502 while (it != freqs_.end())
503 {
504 freqsi_.append(std::distance(freqs_.cbegin(), it));
505 it = std::find_if(std::next(it), freqs_.cend(), margin);
506 }
507 }
508 Pstream::broadcast(freqs_);
509 Pstream::broadcast(freqsi_);
510}
511
512
513void Foam::DMDModels::STDMD::amplitudes()
514{
515 IOField<scalar> snapshot0
516 (
517 IOobject
518 (
519 "snapshot0_" + name_ + "_" + fieldName_,
520 timeName0_,
521 mesh_,
525 )
526 );
527
528 // RxInv^T*(Qupper^T*snapshot0)
529 // tmp = (Qupper^T*snapshot0)
530 List<scalar> tmp(Qupper_.n(), Zero);
531 const label sz0 = snapshot0.size();
532 const label sz1 = tmp.size();
533
534 for (label i = 0; i < sz0; ++i)
535 {
536 for (label j = 0; j < sz1; ++j)
537 {
538 tmp[j] += Qupper_(i,j)*snapshot0[i];
539 }
540 }
541 snapshot0.clear();
542
543 // R = RxInv^T*tmp
544 List<scalar> R(Qupper_.n(), Zero);
545 for (label i = 0; i < sz1; ++i)
546 {
547 for (label j = 0; j < R.size(); ++j)
548 {
549 R[j] += RxInv_(i,j)*tmp[i];
550 }
551 }
552 tmp.clear();
553
554 reduce(R, sumOp<List<scalar>>());
555
556 if (Pstream::master())
557 {
558 Info<< tab << "Computing amplitudes" << endl;
559
560 amps_.resize(R.size());
561 const RCMatrix pEvecs(MatrixTools::pinv(evecs_));
562
563 // amps_ = pEvecs*R;
564 for (label i = 0; i < amps_.size(); ++i)
565 {
566 for (label j = 0; j < R.size(); ++j)
567 {
568 amps_[i] += pEvecs(i,j)*R[j];
569 }
570 }
571 }
572 Pstream::broadcast(amps_);
573}
574
575
576void Foam::DMDModels::STDMD::magnitudes()
577{
578 if (Pstream::master())
579 {
580 Info<< tab << "Computing magnitudes" << endl;
581
582 mags_.resize(amps_.size());
583
584 Info<< tab << "Sorting modes with ";
585
586 switch (modeSorter_)
587 {
588 case modeSorterType::FIRST_SNAPSHOT:
589 {
590 Info<< "method of first snapshot" << endl;
591
592 std::transform
593 (
594 amps_.cbegin(),
595 amps_.cend(),
596 mags_.begin(),
597 [&](const complex& val){ return mag(val); }
598 );
599 break;
600 }
601
602 case modeSorterType::KIEWAT:
603 {
604 Info<< "modified weighted amplitude scaling method" << endl;
605
606 // Eigendecomposition returns evecs with
607 // the unity norm, hence "modeNorm = 1"
608 const scalar modeNorm = 1;
609 const scalar pr = 1;
610 List<scalar> w(step_);
611 std::iota(w.begin(), w.end(), 1);
612 w = sin(twoPi/step_*(w - 1 - 0.25*step_))*pr + pr;
613
614 forAll(mags_, i)
615 {
616 mags_[i] = sorter(w, amps_[i], evals_[i], modeNorm);
617 }
618
619 break;
620 }
621
622 case modeSorterType::KOU_ZHANG:
623 {
624 Info<< "weighted amplitude scaling method" << endl;
625
626 const scalar modeNorm = 1;
627 const List<scalar> w(step_, 1);
628
629 forAll(mags_, i)
630 {
631 mags_[i] = sorter(w, amps_[i], evals_[i], modeNorm);
632 }
633
634 break;
635 }
636
637 default:
638 break;
639 }
640
641 Info<< tab << "Computing magnitude indices" << endl;
642
643 magsi_ = freqsi_;
644
645 auto descend =
646 [&](const label i1, const label i2)
647 {
648 return !(mags_[i1] < mags_[i2]);
649 };
650
651 std::sort(magsi_.begin(), magsi_.end(), descend);
652 }
653 Pstream::broadcast(mags_);
654 Pstream::broadcast(magsi_);
655}
656
657
658Foam::scalar Foam::DMDModels::STDMD::sorter
659(
660 const List<scalar>& weight,
661 const complex& amplitude,
662 const complex& eigenvalue,
663 const scalar modeNorm
664) const
665{
666 // Omit eigenvalues with very large or very small mags
667 if (!(mag(eigenvalue) < GREAT && mag(eigenvalue) > VSMALL))
668 {
669 Info<< " Returning zero magnitude for mag(eigenvalue) = "
670 << mag(eigenvalue) << endl;
671
672 return 0;
673 }
674
675 // Omit eigenvalue-STDMD step combinations that pose a risk of overflow
676 if (mag(eigenvalue)*step_ > sortLimiter_)
677 {
678 Info<< " Returning zero magnitude for"
679 << " mag(eigenvalue) = " << mag(eigenvalue)
680 << " current index = " << step_
681 << " sortLimiter = " << sortLimiter_
682 << endl;
683
684 return 0;
685 }
686
687 scalar magnitude = 0;
688
689 for (label j = 0; j < step_; ++j)
690 {
691 magnitude += weight[j]*modeNorm*mag(amplitude*pow(eigenvalue, j + 1));
692 }
693
694 return magnitude;
695}
696
697
699{
700 SMatrix Atilde(reducedKoopmanOperator());
701 G_.clear();
702
703 if(!eigendecomposition(Atilde))
704 {
705 return false;
706 }
707
708 Atilde.clear();
709
710 frequencies();
711
712 amplitudes();
713
714 magnitudes();
715
716 return true;
717}
718
719
721{
722 Info<< tab << "Computing modes" << endl;
723
724 bool processed = false;
725 processed = processed || modes<scalar>();
726 processed = processed || modes<vector>();
727 processed = processed || modes<sphericalTensor>();
728 processed = processed || modes<symmTensor>();
729 processed = processed || modes<tensor>();
730
731 if (!processed)
732 {
733 return false;
734 }
735
736 return true;
737}
738
739
740void Foam::DMDModels::STDMD::writeToFile(const word& fileName) const
741{
742 Info<< tab << "Writing objects of dynamics" << endl;
743
744 // Write objects of dynamics
745 {
746 autoPtr<OFstream> osPtr =
747 newFileAtTime
748 (
749 fileName + "_" + fieldName_,
750 mesh_.time().timeOutputValue()
751 );
752 OFstream& os = osPtr.ref();
753
754 writeFileHeader(os);
755
756 for (const auto& i : labelRange(0, freqs_.size()))
757 {
758 os << freqs_[i] << tab
759 << mags_[i] << tab
760 << amps_[i].real() << tab
761 << amps_[i].imag() << tab
762 << evals_[i].real() << tab
763 << evals_[i].imag() << endl;
764 }
765 }
766
767 Info<< tab << "Ends output processing for field: " << fieldName_
768 << endl;
769}
770
771
772void Foam::DMDModels::STDMD::writeFileHeader(Ostream& os) const
773{
774 writeHeader(os, "DMD output");
775 writeCommented(os, "Frequency");
776 writeTabbed(os, "Magnitude");
777 writeTabbed(os, "Amplitude (real)");
778 writeTabbed(os, "Amplitude (imag)");
779 writeTabbed(os, "Eigenvalue (real)");
780 writeTabbed(os, "Eigenvalue (imag)");
781 os << endl;
782}
783
784
785void Foam::DMDModels::STDMD::filter()
786{
787 Info<< tab << "Filtering objects of dynamics" << endl;
788
789 // Filter objects according to magsi
790 filterIndexed(evals_, magsi_);
791 filterIndexed(evecs_, magsi_);
792 filterIndexed(freqs_, magsi_);
793 filterIndexed(amps_, magsi_);
794 filterIndexed(mags_, magsi_);
795
796 // Clip objects if need be (assuming objects have the same len)
797 if (freqs_.size() > nModes_)
798 {
799 evals_.resize(nModes_);
800 evecs_.resize(evecs_.m(), nModes_);
801 freqs_.resize(nModes_);
802 amps_.resize(nModes_);
803 mags_.resize(nModes_);
804 }
805}
806
807
808// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
809
811(
812 const fvMesh& mesh,
813 const word& name,
814 const dictionary& dict
815)
816:
818 modeSorter_
819 (
820 modeSorterTypeNames.getOrDefault
821 (
822 "modeSorter",
823 dict,
824 modeSorterType::FIRST_SNAPSHOT
825 )
826 ),
827 Q_(),
828 G_(),
829 Qupper_(1, 1, Zero),
830 Qlower_(1, 1, Zero),
831 RxInv_(1, 1, Zero),
832 evals_(Zero),
833 evecs_(1, 1, Zero),
834 freqs_(Zero),
835 freqsi_(),
836 amps_(Zero),
837 mags_(Zero),
838 magsi_(Zero),
839 patches_
840 (
841 dict.getOrDefault<wordRes>
842 (
843 "patches",
844 dict.found("patch") ? wordRes(1,dict.get<word>("patch")) : wordRes()
845 )
846 ),
847 fieldName_(dict.get<word>("field")),
848 timeName0_(),
849 minBasis_(0),
850 minEval_(0),
851 dt_(0),
852 sortLimiter_(500),
853 fMin_(0),
854 fMax_(pTraits<label>::max),
855 nGramSchmidt_(5),
856 maxRank_(pTraits<label>::max),
857 step_(0),
858 nModes_(pTraits<label>::max),
859 nAgglomerationProcs_(20),
860 empty_(false)
861{}
862
863
864// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
865
867{
868 writeToFile_ = dict.getOrDefault<bool>("writeToFile", true);
869
870 modeSorter_ =
871 modeSorterTypeNames.getOrDefault
872 (
873 "modeSorter",
874 dict,
875 modeSorterType::FIRST_SNAPSHOT
876 );
877
878 minBasis_ =
879 dict.getCheckOrDefault<scalar>("minBasis", 1e-8, scalarMinMax::ge(0));
880
881 minEval_ =
882 dict.getCheckOrDefault<scalar>("minEVal", 1e-8, scalarMinMax::ge(0));
883
884 dt_ =
885 dict.getCheckOrDefault<scalar>
886 (
887 "interval",
888 (
889 dict.getCheck<label>("executeInterval", labelMinMax::ge(1))
890 *mesh_.time().deltaTValue()
891 ),
893 );
894
895 sortLimiter_ =
896 dict.getCheckOrDefault<scalar>("sortLimiter", 500, scalarMinMax::ge(1));
897
898 nGramSchmidt_ =
899 dict.getCheckOrDefault<label>("nGramSchmidt", 5, labelMinMax::ge(1));
900
901 maxRank_ =
902 dict.getCheckOrDefault<label>
903 (
904 "maxRank",
907 );
908
909 nModes_ =
910 dict.getCheckOrDefault<label>
911 (
912 "nModes",
915 );
916
917 fMin_ = dict.getCheckOrDefault<label>("fMin", 0, labelMinMax::ge(0));
918
919 fMax_ =
920 dict.getCheckOrDefault<label>
921 (
922 "fMax",
924 labelMinMax::ge(fMin_ + 1)
925 );
926
927 nAgglomerationProcs_ =
928 dict.getCheckOrDefault<label>
929 (
930 "nAgglomerationProcs",
931 20,
933 );
934
935 Info<< tab << "Settings are read for:" << nl
936 << tab << " field: " << fieldName_ << nl
937 << tab << " modeSorter: " << modeSorterTypeNames[modeSorter_] << nl
938 << tab << " nModes: " << nModes_ << nl
939 << tab << " maxRank: " << maxRank_ << nl
940 << tab << " nGramSchmidt: " << nGramSchmidt_ << nl
941 << tab << " fMin: " << fMin_ << nl
942 << tab << " fMax: " << fMax_ << nl
943 << tab << " minBasis: " << minBasis_ << nl
944 << tab << " minEVal: " << minEval_ << nl
945 << tab << " sortLimiter: " << sortLimiter_ << nl
946 << tab << " nAgglomerationProcs: " << nAgglomerationProcs_ << nl
947 << endl;
948
949 return true;
950}
951
952
953bool Foam::DMDModels::STDMD::initialise(const RMatrix& z)
954{
955 const scalar norm = L2norm(z);
956
957 if (mag(norm) > 0)
958 {
959 // First-processed snapshot required by the mode-sorting
960 // algorithms at the final output computations (K:p. 43)
961 {
962 const label nSnap = z.m()/2;
963 timeName0_ =
964 mesh_.time().timeName(mesh_.time().startTime().value());
965
966 if (nSnap == 0)
967 {
968 empty_ = true;
969 }
970
971 IOField<scalar> snapshot0
972 (
974 (
975 "snapshot0_" + name_ + "_" + fieldName_,
976 timeName0_,
977 mesh_,
981 ),
982 nSnap
983 );
984
985 std::copy
986 (
987 z.cbegin(),
988 z.cbegin() + nSnap,
989 snapshot0.begin()
990 );
991
992 const IOstreamOption streamOpt
993 (
994 mesh_.time().writeFormat(),
995 mesh_.time().writeCompression()
996 );
997
998 fileHandler().writeObject(snapshot0, streamOpt, true);
999 }
1000
1001 Q_ = z/norm;
1002 G_ = SMatrix(1);
1003 G_(0,0) = sqr(norm);
1004
1005 ++step_;
1006
1007 return true;
1008 }
1009
1010 return false;
1011}
1012
1013
1014bool Foam::DMDModels::STDMD::update(const RMatrix& z)
1015{
1016 {
1017 //- Working copy of the augmented snapshot matrix "z"
1018 //- being used in the classical Gram-Schmidt process
1019 const RMatrix ez(orthonormalise(z));
1020
1021 // Check basis for "z" and, if necessary, expand "Q" and "G"
1022 const scalar ezNorm = L2norm(ez);
1023 const scalar zNorm = L2norm(z);
1024
1025 if (ezNorm/zNorm > minBasis_)
1026 {
1027 expand(ez, ezNorm);
1028 }
1029 }
1030
1031 // Update "G" before the potential orthonormal basis compression
1032 updateG(z);
1033
1034 // Compress the orthonormal basis if required
1035 if (Q_.n() >= maxRank_)
1036 {
1037 compress();
1038 }
1039
1040 ++step_;
1042 return true;
1043}
1044
1045
1047{
1048 // DMD statistics and mode evaluation (K:Fig. 16)
1049 const label nSnap = Q_.m()/2;
1050
1051 // Move upper and lower halves of "Q" to new containers
1052 Qupper_ = RMatrix(Q_.subMatrix(0, 0, max(nSnap, 1)));
1053 Qlower_ = RMatrix(Q_.subMatrix(nSnap, 0, max(nSnap, 1)));
1054 Q_.clear();
1055
1056 if (!dynamics())
1057 {
1058 return true;
1059 }
1060
1061 modes();
1062
1063 if (Pstream::master() && writeToFile_)
1064 {
1065 writeToFile(word("dynamics"));
1066
1067 filter();
1068
1069 writeToFile(word("filtered_dynamics"));
1070 }
1071
1072 step_ = 0;
1073
1074 return true;
1075}
1076
1077
1078// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Abstract base class for DMD models to handle DMD characteristics for the DMD function object.
Definition DMDModel.H:65
virtual bool modes()=0
Compute and write modes.
const fvMesh & mesh_
Reference to the mesh.
Definition DMDModel.H:77
DMDModel(const fvMesh &mesh, const word &name, const dictionary &dict)
Construct from components.
Definition DMDModel.C:35
virtual bool dynamics()=0
Compute and write mode dynamics.
const word name_
Name of operand function object.
Definition DMDModel.H:82
Streaming Total Dynamic Mode Decomposition (i.e. STDMD) is a variant of dynamic mode decomposition.
Definition STDMD.H:298
STDMD(const fvMesh &mesh, const word &name, const dictionary &dict)
Construct from components.
Definition STDMD.C:804
virtual bool update(const RMatrix &z)
Incremental orthonormal basis update (K:Fig. 15).
Definition STDMD.C:1007
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition STDMD.C:859
virtual bool initialise(const RMatrix &z)
Initialise 'Q' and 'G' (both require the first two snapshots).
Definition STDMD.C:946
virtual bool fit()
Compute and write modes and mode dynamics of model data members.
Definition STDMD.C:1041
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A simple container for options an IOstream can normally have.
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
static direction m() noexcept
The number of rows.
Definition MatrixSpace.H:79
static direction n() noexcept
The number of columns.
Definition MatrixSpace.H:84
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing a constant Matrix.
Definition MatrixI.H:525
label m() const noexcept
The number of rows.
Definition Matrix.H:261
ConstMatrixBlock< mType > subColumn(const label colIndex, const label rowIndex=0, label len=-1) const
Return const column or column's subset of Matrix.
Definition MatrixI.H:263
static MinMax< scalar > ge(const scalar &minVal)
bool allowClearRecv() const noexcept
Is clearStorage of individual receive buffer by external hooks allowed? (default: true).
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
void clear()
Clear all send/recv buffers and reset states.
@ FALSE
switch off column pivoting
Definition QRMatrix.H:242
@ ONLY_R
compute only R
Definition QRMatrix.H:233
@ ECONOMY
compute economy-size decomposition
Definition QRMatrix.H:224
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
T & ref()
Return reference to the managed object without nullptr checking.
Definition autoPtr.H:231
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool writeToFile_
Flag to enable/disable writing to file.
Definition writeFile.H:146
virtual bool writeToFile() const
Flag to allow writing to file.
Definition writeFile.C:286
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
A namespace for various dynamic mode decomposition (DMD) model implementations.
Definition STDMD.C:33
MatrixType pinv(const MatrixType &A, scalar tol=1e-5)
Moore-Penrose inverse of singular/non-singular square/rectangular scalar/complex matrices (KPP:p....
Definition QRMatrix.C:577
Mathematical constants.
constexpr scalar twoPi(2 *M_PI)
@ complex
Definition Roots.H:55
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
tensor Rx(const scalar omega)
Rotational transformation tensor about the x-axis by omega radians.
Definition transform.H:83
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler().
dimensionedScalar sin(const dimensionedScalar &ds)
static void writeHeader(Ostream &os, const word &fieldName)
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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).
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
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
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
const volScalarField & cp
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299