Loading...
Searching...
No Matches
STDMD.H
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-2025 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
26Class
27 Foam::DMDModels::STDMD
28
29Description
30 Streaming Total Dynamic Mode Decomposition (i.e. \c STDMD)
31 is a variant of dynamic mode decomposition.
32
33 Among other Dynamic Mode Decomposition (DMD) variants, \c STDMD is
34 presumed to provide the general DMD method capabilities alongside
35 economised and feasible memory and CPU usage.
36
37 The code implementation corresponds to Figs. 15-16 of the first citation
38 below, more broadly to Section 2.4.
39
40 References:
41 \verbatim
42 DMD and mode-sorting algorithms (tags:K, HRDC, KZ, HWR):
43 Kiewat, M. (2019).
44 Streaming modal decomposition approaches for vehicle aerodynamics.
45 PhD thesis. Munich: Technical University of Munich.
46 URL:mediatum.ub.tum.de/doc/1482652/1482652.pdf
47
48 Hemati, M. S., Rowley, C. W.,
49 Deem, E. A., & Cattafesta, L. N. (2017).
50 De-biasing the dynamic mode decomposition
51 for applied Koopman spectral analysis of noisy datasets.
52 Theoretical and Computational Fluid Dynamics, 31(4), 349-368.
53 DOI:10.1007/s00162-017-0432-2
54
55 Kou, J., & Zhang, W. (2017).
56 An improved criterion to select
57 dominant modes from dynamic mode decomposition.
58 European Journal of Mechanics-B/Fluids, 62, 109-129.
59 DOI:10.1016/j.euromechflu.2016.11.015
60
61 Hemati, M. S., Williams, M. O., & Rowley, C. W. (2014).
62 Dynamic mode decomposition for large and streaming datasets.
63 Physics of Fluids, 26(11), 111701.
64 DOI:10.1063/1.4901016
65
66 Parallel classical Gram-Schmidt process (tag:Ka):
67 Katagiri, T. (2003).
68 Performance evaluation of parallel
69 Gram-Schmidt re-orthogonalization methods.
70 In: Palma J. M. L. M., Sousa A. A., Dongarra J., Hernández V. (eds)
71 High Performance Computing for Computational Science — VECPAR 2002.
72 Lecture Notes in Computer Science, vol 2565, p. 302-314.
73 Berlin, Heidelberg: Springer.
74 DOI:10.1007/3-540-36569-9_19
75
76 Parallel direct tall-skinny QR decomposition (tags:BGD, DGHL):
77 Benson, A. R., Gleich, D. F., & Demmel, J. (2013).
78 Direct QR factorizations for
79 tall-and-skinny matrices in MapReduce architectures.
80 2013 IEEE International Conference on Big Data.
81 DOI:10.1109/bigdata.2013.6691583
82
83 Demmel, J., Grigori, L., Hoemmen, M., & Langou, J. (2012).
84 Communication-optimal parallel
85 and sequential QR and LU factorizations.
86 SIAM Journal on Scientific Computing, 34(1), A206-A239.
87 DOI:10.1137/080731992
88 \endverbatim
89
90Usage
91 Minimal example by using \c system/controlDict.functions:
92 \verbatim
93 DMDFO
94 {
95 // Mandatory entries
96 DMDModel STDMD;
97
98 // Conditional mandatory entries
99
100 // Option-1
101 interval 5.5;
102
103 // Option-2
104 executeInterval 10;
105
106 // Optional entries
107 modeSorter kiewat;
108 nGramSchmidt 5;
109 maxRank 50;
110 nModes 50;
111 fMin 0;
112 fMax 1000000000;
113 nAgglomerationProcs 20;
114
115 // Optional entries (not recommended to change)
116 minBasis 0.00000001;
117 minEVal 0.00000001;
118 sortLimiter 500.0;
119
120 // Inherited entries
121 ...
122 }
123 \endverbatim
124
125 where the entries mean:
126 \table
127 Property | Description | Type | Reqd | Deflt
128 DMDModel | Type: STDMD | word | yes | -
129 interval | STDMD time-step size [s] | scalar | cndtnl <!--
130 --> | executeInterval*(current time-step of the simulation)
131 modeSorter | Mode-sorting algorithm model | word | no <!--
132 --> | firstSnapshot
133 nModes | Number of output modes in input frequency range <!--
134 --> | label | no | GREAT
135 maxRank | Max columns in accumulated matrices | label | no | GREAT
136 nGramSchmidt | Number of Gram-Schmidt iterations | label | no | 5
137 fMin | Min (non-negative) output frequency | label | no | 0
138 fMax | Max output frequency | label | no | GREAT
139 nAgglomerationProcs | Number of processors at each agglomeration <!--
140 --> unit during the computation of reduced Koopman <!--
141 --> operator | label | no | 20
142 minBasis | Orthogonal basis expansion threshold | scalar| no | 1e-8
143 minEVal | Min eigenvalue for below eigenvalues are omitted <!--
144 --> | scalar| no | 1e-8
145 sortLimiter | Max allowable magnitude for <!--
146 --> mag(eigenvalue)*(number of DMD steps) to be used in <!--
147 --> modeSorter={kiewat,kouZhang} to avoid overflow errors <!--
148 --> | scalar | no | 500.0
149 \endtable
150
151 Options for the \c modeSorter entry:
152 \verbatim
153 kiewat | Modified weighted-amplitude scaling method
154 kouZhang | Weighted-amplitude scaling method
155 firstSnapshot | First-snapshot amplitude magnitude method
156 \endverbatim
157
158Note
159 - To specify the STDMD time-step size (not necessarily equal to the
160 time step of the simulation), entries of either \c interval or
161 \c executeInterval must be available (highest to lowest precedence). While
162 \c interval allows users to directly specify the STDMD time-step size
163 in seconds, in absence of \c interval, for convenience,
164 \c executeInterval allows users to compute the STDMD time-step internally
165 by multiplying itself with the current time-step size of the simulation.
166 - Limitation: Restart is currently not available since intermediate writing
167 of STDMD matrices are not supported.
168 - Limitation: Non-physical input (e.g. full-zero fields) can upset STDMD.
169 - Warning: In \c modeSorter algorithms of \c kiewat and \c kouZhang, a
170 function of \f$power(x,y)\f$ exists where \c x is the magnitude of an
171 eigenvalue, and \c y is the number of STDMD snapshots. This function poses
172 a risk of overflow errors since, for example, if \c x ends up above 1.5 with
173 500 snapshots or more, this function automatically throws the floating point
174 error meaning overflow. Therefore, the domain-range of this function is
175 heuristically constrained by the optional entry \c sortLimiter which skips
176 the evaluation of this function for a given
177 mag(eigenvalue)*(no. of STDMD steps), i.e. x*y, whose magnitude is larger
178 than \c sortLimiter.
179
180SourceFiles
181 STDMD.C
182 STDMDTemplates.C
183
184\*---------------------------------------------------------------------------*/
185
186#ifndef Foam_DMDModels_STDMD_H
187#define Foam_DMDModels_STDMD_H
188
189#include "DMDModel.H"
190#include "SquareMatrix.H"
191#include "RectangularMatrix.H"
192
193// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194
195namespace Foam
196{
197namespace DMDModels
198{
199
200/*---------------------------------------------------------------------------*\
201 Class STDMD Declaration
202\*---------------------------------------------------------------------------*/
203
204class STDMD
205:
206 public DMDModel
207{
208 typedef SquareMatrix<scalar> SMatrix;
209 typedef RectangularMatrix<scalar> RMatrix;
210 typedef RectangularMatrix<complex> RCMatrix;
211
212
213 // Private Enumerations
214
215 //- Options for the mode-sorting algorithm
216 enum modeSorterType : char
217 {
218 FIRST_SNAPSHOT = 0,
219 KIEWAT,
220 KOU_ZHANG
221 };
222
223 //- Names for modeSorterType
224 static const Enum<modeSorterType> modeSorterTypeNames;
225
226
227 // Private Data
228
229 //- Mode-sorting algorithm
230 enum modeSorterType modeSorter_;
231
232 //- Accumulated-in-time unitary similarity matrix produced by the
233 //- orthonormal decomposition of the augmented snapshot matrix 'z'
234 // (K:Eq. 60)
235 RMatrix Q_;
236
237 //- Accumulated-in-time (squared) upper triangular matrix produced by
238 //- the orthonormal decomposition of the augmented snapshot matrix 'z'
239 // (K:Eq. 64)
240 SMatrix G_;
241
242 //- Upper half of 'Q'
243 RMatrix Qupper_;
244
245 //- Lower half of 'Q'
246 RMatrix Qlower_;
247
248 //- Moore-Penrose pseudo-inverse of 'R' produced by
249 //- the QR decomposition of the last time-step 'Q'
250 RMatrix RxInv_;
251
252 //- Eigenvalues of modes
253 List<complex> evals_;
254
255 //- Eigenvectors of modes
256 RCMatrix evecs_;
257
258 //- (Non-negative) frequencies of modes
259 List<scalar> freqs_;
260
261 //- Indices of 'freqs' where frequencies are
262 //- non-negative and within ['fMin', 'fMax']
263 DynamicList<label> freqsi_;
264
265 //- Amplitudes of modes
266 List<complex> amps_;
267
268 //- (Descending) magnitudes of (complex) amplitudes of modes
269 List<scalar> mags_;
270
271 //- Indices of 'mags'
272 List<label> magsi_;
273
274 //- Names of operand patches
275 const wordRes patches_;
276
277 //- Name of operand field
278 const word fieldName_;
279
280 //- First-processed snapshot required by the mode-sorting
281 //- algorithms at the final output computations (K:p. 43)
282 word timeName0_;
283
284 //- Min value to execute orthogonal basis expansion of 'Q' and 'G'
285 scalar minBasis_;
286
287 //- Min eigenvalue magnitude below where 'evals' are omitted
288 scalar minEval_;
289
290 //- STDMD time-step size that equals to
291 //- (executeInterval of DMD)*(deltaT of simulation) [s]
292 scalar dt_;
293
294 //- Maximum allowable magnitude for mag(eigenvalue)*(number of
295 //- STDMD steps) to be used in modeSorter={kiewat,kouZhang} to
296 //- avoid overflow errors
297 scalar sortLimiter_;
298
299 //- Min frequency: Output only entries corresponding to freqs >= 'fMin'
300 label fMin_;
301
302 //- Max frequency: Output only entries corresponding to freqs <= 'fMax'
303 label fMax_;
304
305 //- Number of maximum iterations of the classical
306 //- Gram-Schmidt procedure for the orthonormalisation
307 label nGramSchmidt_;
308
309 //- Maximum allowable matrix column for 'Q' and 'G'
310 // 'Q' is assumed to always have full-rank, thus 'Q.n() = rank'
311 label maxRank_;
312
313 //- Current execution-step index of STDMD,
314 //- not necessarily that of the simulation
315 label step_;
316
317 //- Number of output modes within input frequency range
318 //- starting from the most energetic mode
319 label nModes_;
320
321 //- Number of processors at each agglomeration unit
322 //- during the computation of reduced Koopman operator
323 label nAgglomerationProcs_;
324
325 //- (Internal) Flag to tag snapshots which are effectively empty
326 bool empty_;
327
328
329 // Private Member Functions
330
331 // Evaluation
332
333 //- Return (parallel) L2-norm of a given column vector
334 scalar L2norm(const RMatrix& z) const;
335
336 //- Execute (parallel) classical Gram-Schmidt
337 //- process to orthonormalise 'ez' (Ka:Fig. 5)
338 RMatrix orthonormalise(RMatrix ez) const;
339
340 //- Expand orthonormal bases 'Q' and 'G' by stacking a column
341 //- '(ez/ezNorm)' to 'Q', and a row (Zero) and column (Zero)
342 //- to 'G' if '(ezNorm/zNorm > minBasis)'
343 void expand(const RMatrix& ez, const scalar ezNorm);
344
345 //- Update 'G' before the potential orthonormal basis compression
346 void updateG(const RMatrix& z);
347
348 //- Compress orthonormal basis for 'Q' and 'G' if '(Q.n()>maxRank)'
349 void compress();
350
351 //- Return reduced Koopman operator 'Atilde' (K:Eq. 78)
352 // Also fills 'RxInv'.
353 // The function was not divided into subsections to ensure
354 // minimal usage of memory, hence the complexity of the function
355 SMatrix reducedKoopmanOperator();
356
357 //- Compute eigenvalues and eigenvectors of 'Atilde'
358 // Remove any eigenvalues whose magnitude is smaller than
359 // 'minEVal' while keeping the order of elements the same
360 bool eigendecomposition(SMatrix& Atilde);
361
362 //- Compute and filter frequencies and its indices (K:Eq. 81)
363 // Indices of 'freqs' where frequencies are
364 // non-negative and within ['fMin', 'fMax']
365 void frequencies();
366
367 //- Compute amplitudes
368 void amplitudes();
369
370 //- Compute magnitudes and ordered magnitude indices
371 // Includes advanced sorting methods using
372 // the chosen weighted amplitude scaling method
373 void magnitudes();
374
375 //- Eigenvalue weighted amplitude scaling (KZ:Eq. 33)
376 //- Modified eigenvalue weighted amplitude scaling (K)
377 scalar sorter
378 (
379 const List<scalar>& weight,
380 const complex& amplitude,
381 const complex& eigenvalue,
382 const scalar modeNorm
383 ) const;
384
385 //- Compute and write mode dynamics
386 virtual bool dynamics();
387
388 //- Compute and write modes
389 virtual bool modes();
390
391 //- Select field type for modes
392 template<class Type>
393 bool modes();
394
395 //- Compute modes based on input field type
396 template<class GeoFieldType>
397 bool calcModes();
398
399 //- Compute a mode for a given scalar-based or other input field
400 template<class GeoFieldType>
401 void calcMode
402 (
403 GeoFieldType& modeRe,
404 GeoFieldType& modeIm,
405 const RMatrix& primitiveMode,
406 const label magi,
407 const label rowi = 0
408 );
409
410
411 // I-O
412
413 //- Write objects of dynamics
414 void writeToFile(const word& fileName) const;
415
416 // Notifying the compiler that we want both writeToFile functions
418
419 //- Write file header information
420 void writeFileHeader(Ostream& os) const;
421
422 //- Filter objects of dynamics according to 'freqsi' and 'magsi'
423 void filter();
424
425 //- Return a new List containing elems of List at 'indices'
426 template<class Type>
427 void filterIndexed
428 (
429 List<Type>& lst,
430 const UList<label>& indices
431 );
432
433 //- Return a new Matrix containing columns of Matrix at 'indices'
434 template<class MatrixType>
435 void filterIndexed
436 (
437 MatrixType& lst,
438 const UList<label>& indices
439 );
440
441
442public:
443
444 //- Runtime type information
445 TypeName("STDMD");
446
447
448 // Constructors
449
450 //- Construct from components
451 STDMD
452 (
453 const fvMesh& mesh,
454 const word& name,
455 const dictionary& dict
456 );
457
458 //- No copy construct
459 STDMD(const STDMD&) = delete;
460
461 //- No copy assignment
462 void operator=(const STDMD&) = delete;
463
464
465 //- Destructor
466 virtual ~STDMD() = default;
467
468
469 // Member Functions
470
471 // Evaluation
472
473 //- Initialise 'Q' and 'G' (both require the first two snapshots)
474 virtual bool initialise(const RMatrix& z);
475
476 //- Incremental orthonormal basis update (K:Fig. 15)
477 virtual bool update(const RMatrix& z);
478
479 //- Compute and write modes and
480 //- mode dynamics of model data members
481 virtual bool fit();
482
483
484 // I-O
485
486 //- Read the function-object dictionary
487 virtual bool read(const dictionary& dict);
488};
489
490
491// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
492
493} // End namespace DMDModels
494} // End namespace Foam
495
496// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
497
498#ifdef NoRepository
499 #include "STDMDTemplates.C"
500#endif
501
502// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503
504#endif
505
506// ************************************************************************* //
virtual bool modes()=0
Compute and write modes.
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.
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 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
TypeName("STDMD")
Runtime type information.
void operator=(const STDMD &)=delete
No copy assignment.
virtual bool fit()
Compute and write modes and mode dynamics of model data members.
Definition STDMD.C:1041
STDMD(const STDMD &)=delete
No copy construct.
virtual ~STDMD()=default
Destructor.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A templated (M x N) rectangular matrix of objects of <Type>, containing M*N elements,...
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
A complex number, similar to the C++ complex type.
Definition complex.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
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 List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh update()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
A namespace for various dynamic mode decomposition (DMD) model implementations.
Definition STDMD.C:33
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68