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-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::DMDModels::STDMD
28 
29 Description
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 
90 Usage
91  Minimal example by using \c system/controlDict.functions:
92  \verbatim
93  DMD1
94  {
95  // Mandatory/Optional entries
96  ...
97 
98  // Mandatory entries (unmodifiable)
99  DMDModel STDMD;
100 
101  // Conditional mandatory entries (runtime modifiable)
102 
103  // Option-1
104  interval 5.5;
105 
106  // Option-2
107  executeInterval 10;
108 
109  // Optional entries (runtime modifiable)
110  modeSorter kiewat;
111  nGramSchmidt 5;
112  maxRank 50;
113  nModes 50;
114  fMin 0;
115  fMax 1000000000;
116  nAgglomerationProcs 20;
117 
118  // Optional entries (runtime modifiable, yet not recommended)
119  minBasis 0.00000001;
120  minEVal 0.00000001;
121  sortLimiter 500.0;
122 
123  // Mandatory/Optional (inherited) entries
124  ...
125  }
126  \endverbatim
127 
128  where the entries mean:
129  \table
130  Property | Description | Type | Reqd | Deflt
131  DMDModel | Type: STDMD | word | yes | -
132  interval | STDMD time-step size [s] | scalar | cndtnl <!--
133  --> | executeInterval*(current time-step of the simulation)
134  modeSorter | Mode-sorting algorithm model | word | no <!--
135  --> | firstSnapshot
136  nModes | Number of output modes in input frequency range <!--
137  --> | label | no | GREAT
138  maxRank | Max columns in accumulated matrices | label | no | GREAT
139  nGramSchmidt | Number of Gram-Schmidt iterations | label | no | 5
140  fMin | Min (non-negative) output frequency | label | no | 0
141  fMax | Max output frequency | label | no | GREAT
142  nAgglomerationProcs | Number of processors at each agglomeration <!--
143  --> unit during the computation of reduced Koopman <!--
144  --> operator | label | no | 20
145  minBasis | Orthogonal basis expansion threshold | scalar| no | 1e-8
146  minEVal | Min eigenvalue for below eigenvalues are omitted <!--
147  --> | scalar| no | 1e-8
148  sortLimiter | Max allowable magnitude for <!--
149  --> mag(eigenvalue)*(number of DMD steps) to be used in <!--
150  --> modeSorter={kiewat,kouZhang} to avoid overflow errors <!--
151  --> | scalar | no | 500.0
152  \endtable
153 
154  Options for the \c modeSorter entry:
155  \verbatim
156  kiewat | Modified weighted-amplitude scaling method
157  kouZhang | Weighted-amplitude scaling method
158  firstSnapshot | First-snapshot amplitude magnitude method
159  \endverbatim
160 
161 Note
162  - To specify the STDMD time-step size (not necessarily equal to the
163  time step of the simulation), entries of either \c interval or
164  \c executeInterval must be available (highest to lowest precedence). While
165  \c interval allows users to directly specify the STDMD time-step size
166  in seconds, in absence of \c interval, for convenience,
167  \c executeInterval allows users to compute the STDMD time-step internally
168  by multiplying itself with the current time-step size of the simulation.
169  - Limitation: Restart is currently not available since intermediate writing
170  of STDMD matrices are not supported.
171  - Limitation: Non-physical input (e.g. full-zero fields) can upset STDMD.
172  - Warning: In \c modeSorter algorithms of \c kiewat and \c kouZhang, a
173  function of \f$power(x,y)\f$ exists where \c x is the magnitude of an
174  eigenvalue, and \c y is the number of STDMD snapshots. This function poses
175  a risk of overflow errors since, for example, if \c x ends up above 1.5 with
176  500 snapshots or more, this function automatically throws the floating point
177  error meaning overflow. Therefore, the domain-range of this function is
178  heuristically constrained by the optional entry \c sortLimiter which skips
179  the evaluation of this function for a given
180  mag(eigenvalue)*(no. of STDMD steps), i.e. x*y, whose magnitude is larger
181  than \c sortLimiter.
182 
183 See also
184  - Foam::functionObjects::DMD
185  - Foam::DMDModel
186 
187 SourceFiles
188  STDMD.C
189  STDMDTemplates.C
190 
191 \*---------------------------------------------------------------------------*/
192 
193 #ifndef DMDModels_STDMD_H
194 #define DMDModels_STDMD_H
195 
196 #include "DMDModel.H"
197 #include "SquareMatrix.H"
198 #include "RectangularMatrix.H"
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 namespace Foam
203 {
204 namespace DMDModels
205 {
206 
207 /*---------------------------------------------------------------------------*\
208  Class STDMD Declaration
209 \*---------------------------------------------------------------------------*/
210 
211 class STDMD
212 :
213  public DMDModel
214 {
215  typedef SquareMatrix<scalar> SMatrix;
216  typedef RectangularMatrix<scalar> RMatrix;
217  typedef RectangularMatrix<complex> RCMatrix;
218 
219 
220  // Private Enumerations
221 
222  //- Options for the mode-sorting algorithm
223  enum modeSorterType : char
224  {
225  FIRST_SNAPSHOT = 0,
226  KIEWAT,
227  KOU_ZHANG
228  };
229 
230  //- Names for modeSorterType
231  static const Enum<modeSorterType> modeSorterTypeNames;
232 
233 
234  // Private Data
235 
236  //- Mode-sorting algorithm
237  enum modeSorterType modeSorter_;
238 
239  //- Accumulated-in-time unitary similarity matrix produced by the
240  //- orthonormal decomposition of the augmented snapshot matrix 'z'
241  // (K:Eq. 60)
242  RMatrix Q_;
243 
244  //- Accumulated-in-time (squared) upper triangular matrix produced by
245  //- the orthonormal decomposition of the augmented snapshot matrix 'z'
246  // (K:Eq. 64)
247  SMatrix G_;
248 
249  //- Upper half of 'Q'
250  RMatrix Qupper_;
251 
252  //- Lower half of 'Q'
253  RMatrix Qlower_;
254 
255  //- Moore-Penrose pseudo-inverse of 'R' produced by
256  //- the QR decomposition of the last time-step 'Q'
257  RMatrix RxInv_;
258 
259  //- Eigenvalues of modes
260  List<complex> evals_;
261 
262  //- Eigenvectors of modes
263  RCMatrix evecs_;
264 
265  //- (Non-negative) frequencies of modes
266  List<scalar> freqs_;
267 
268  //- Indices of 'freqs' where frequencies are
269  //- non-negative and within ['fMin', 'fMax']
270  DynamicList<label> freqsi_;
271 
272  //- Amplitudes of modes
273  List<complex> amps_;
274 
275  //- (Descending) magnitudes of (complex) amplitudes of modes
276  List<scalar> mags_;
277 
278  //- Indices of 'mags'
279  List<label> magsi_;
280 
281  //- Name of operand field
282  const word fieldName_;
283 
284  //- Name of operand patch
285  const word patch_;
286 
287  //- First-processed snapshot required by the mode-sorting
288  //- algorithms at the final output computations (K:p. 43)
289  word timeName0_;
290 
291  //- Min value to execute orthogonal basis expansion of 'Q' and 'G'
292  scalar minBasis_;
293 
294  //- Min eigenvalue magnitude below where 'evals' are omitted
295  scalar minEval_;
296 
297  //- STDMD time-step size that equals to
298  //- (executeInterval of DMD)*(deltaT of simulation) [s]
299  scalar dt_;
300 
301  //- Maximum allowable magnitude for mag(eigenvalue)*(number of
302  //- STDMD steps) to be used in modeSorter={kiewat,kouZhang} to
303  //- avoid overflow errors
304  scalar sortLimiter_;
305 
306  //- Min frequency: Output only entries corresponding to freqs >= 'fMin'
307  label fMin_;
308 
309  //- Max frequency: Output only entries corresponding to freqs <= 'fMax'
310  label fMax_;
311 
312  //- Number of maximum iterations of the classical
313  //- Gram-Schmidt procedure for the orthonormalisation
314  label nGramSchmidt_;
315 
316  //- Maximum allowable matrix column for 'Q' and 'G'
317  // 'Q' is assumed to always have full-rank, thus 'Q.n() = rank'
318  label maxRank_;
319 
320  //- Current execution-step index of STDMD,
321  //- not necessarily that of the simulation
322  label step_;
323 
324  //- Number of output modes within input frequency range
325  //- starting from the most energetic mode
326  label nModes_;
327 
328  //- Number of processors at each agglomeration unit
329  //- during the computation of reduced Koopman operator
330  label nAgglomerationProcs_;
331 
332  //- (Internal) Flag to tag snapshots which are effectively empty
333  bool empty_;
334 
335 
336  // Private Member Functions
337 
338  // Process
339 
340  //- Return (parallel) L2-norm of a given column vector
341  scalar L2norm(const RMatrix& z) const;
342 
343  //- Execute (parallel) classical Gram-Schmidt
344  //- process to orthonormalise 'ez' (Ka:Fig. 5)
345  RMatrix orthonormalise(RMatrix ez) const;
346 
347  //- Expand orthonormal bases 'Q' and 'G' by stacking a column
348  //- '(ez/ezNorm)' to 'Q', and a row (Zero) and column (Zero)
349  //- to 'G' if '(ezNorm/zNorm > minBasis)'
350  void expand(const RMatrix& ez, const scalar ezNorm);
351 
352  //- Compress orthonormal basis for 'Q' and 'G' if '(Q.n()>maxRank)'
353  void compress();
354 
355 
356  // Evaluation
357 
358  //- Return reduced Koopman operator 'Atilde' (K:Eq. 78)
359  // Also fills 'RxInv'.
360  // The function was not divided into subsections to ensure
361  // minimal usage of memory, hence the complexity of the function
362  SMatrix reducedKoopmanOperator();
363 
364  //- Compute eigenvalues and eigenvectors of 'Atilde'
365  // Remove any eigenvalues whose magnitude is smaller than
366  // 'minEVal' while keeping the order of elements the same
367  bool eigendecomposition(SMatrix& Atilde);
368 
369  //- Compute and filter frequencies and its indices (K:Eq. 81)
370  // Indices of 'freqs' where frequencies are
371  // non-negative and within ['fMin', 'fMax']
372  void frequencies();
373 
374  //- Compute amplitudes
375  void amplitudes();
376 
377  //- Compute magnitudes and ordered magnitude indices
378  // Includes advanced sorting methods using
379  // the chosen weighted amplitude scaling method
380  void magnitudes();
381 
382  //- Eigenvalue weighted amplitude scaling (KZ:Eq. 33)
383  //- Modified eigenvalue weighted amplitude scaling (K)
384  scalar sorter
385  (
386  const List<scalar>& weight,
387  const complex& amplitude,
388  const complex& eigenvalue,
389  const scalar modeNorm
390  ) const;
391 
392  //- Compute and write mode dynamics
393  virtual bool dynamics();
394 
395  //- Compute and write modes
396  virtual bool modes();
397 
398  //- Select field type for modes
399  template<class Type>
400  bool modes();
401 
402  //- Compute modes based on input field type
403  template<class GeoFieldType>
404  bool calcModes();
405 
406  //- Compute a mode for a given scalar-based input field
407  template<class GeoFieldType>
408  typename std::enable_if
409  <
410  std::is_same<scalar, typename GeoFieldType::value_type>::value,
411  void
412  >::type calcMode
413  (
414  GeoFieldType& modeRe,
415  GeoFieldType& modeIm,
416  const RMatrix& primitiveMode,
417  const label i
418  );
419 
420  //- Compute a mode for a given non-scalar-based input field
421  template<class GeoFieldType>
422  typename std::enable_if
423  <
424  !std::is_same<scalar, typename GeoFieldType::value_type>::value,
425  void
426  >::type calcMode
427  (
428  GeoFieldType& modeRe,
429  GeoFieldType& modeIm,
430  const RMatrix& primitiveMode,
431  const label i
432  );
433 
434 
435  // IO
436 
437  //- Write objects of dynamics
438  void writeToFile(const word& fileName) const;
439 
440  // Notifying the compiler that we want both writeToFile functions
442 
443  //- Write file header information
444  void writeFileHeader(Ostream& os) const;
445 
446  //- Filter objects of dynamics according to 'freqsi' and 'magsi'
447  void filter();
448 
449  //- Return a new List containing elems of List at 'indices'
450  template<class Type>
451  void filterIndexed
452  (
453  List<Type>& lst,
454  const UList<label>& indices
455  );
456 
457  //- Return a new Matrix containing columns of Matrix at 'indices'
458  template<class MatrixType>
459  void filterIndexed
460  (
461  MatrixType& lst,
462  const UList<label>& indices
463  );
464 
465 
466 public:
467 
468  //- Runtime type information
469  TypeName("STDMD");
470 
471 
472  // Constructors
473 
474  //- Construct from components
475  STDMD
476  (
477  const fvMesh& mesh,
478  const word& name,
479  const dictionary& dict
480  );
481 
482  //- No copy construct
483  STDMD(const STDMD&) = delete;
484 
485  //- No copy assignment
486  void operator=(const STDMD&) = delete;
487 
488 
489  //- Destructor
490  virtual ~STDMD() = default;
491 
492 
493  // Member Functions
494 
495  // Evaluation
496 
497  //- Initialise 'Q' and 'G' (both require the first two snapshots)
498  virtual bool initialise(const RMatrix& z);
499 
500  //- Incremental orthonormal basis update (K:Fig. 15)
501  virtual bool update(const RMatrix& z);
502 
503  //- Compute and write modes and
504  //- mode dynamics of model data members
505  virtual bool fit();
506 
507 
508  // IO
509 
510  //- Read STDMD settings
511  virtual bool read(const dictionary& dict);
512 };
513 
514 
515 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
516 
517 } // End namespace DMDModels
518 } // End namespace Foam
519 
520 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
521 
522 #ifdef NoRepository
523  #include "STDMDTemplates.C"
524 #endif
525 
526 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
527 
528 #endif
529 
530 // ************************************************************************* //
Foam::Enum< modeSorterType >
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::DynamicList< label >
Foam::DMDModels::STDMD::TypeName
TypeName("STDMD")
Runtime type information.
Foam::DMDModels::STDMD::fit
virtual bool fit()
Definition: STDMD.C:965
Foam::RectangularMatrix< scalar >
Foam::DMDModels::STDMD::operator=
void operator=(const STDMD &)=delete
No copy assignment.
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::DMDModels::STDMD::~STDMD
virtual ~STDMD()=default
Destructor.
Foam::DMDModels::STDMD::read
virtual bool read(const dictionary &dict)
Read STDMD settings.
Definition: STDMD.C:785
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::DMDModels::STDMD::initialise
virtual bool initialise(const RMatrix &z)
Initialise 'Q' and 'G' (both require the first two snapshots)
Definition: STDMD.C:872
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::complex
A complex number, similar to the C++ complex type.
Definition: complex.H:82
Foam::DMDModel
Abstract base class for DMD models to handle DMD characteristics for the DMD function object.
Definition: DMDModel.H:65
Foam::DMDModels::STDMD::update
virtual bool update(const RMatrix &z)
Incremental orthonormal basis update (K:Fig. 15)
Definition: STDMD.C:928
Foam::SquareMatrix< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::UList< label >
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
STDMDTemplates.C
Foam::functionObjects::writeFile::writeToFile
virtual bool writeToFile() const
Flag to allow writing to file.
Definition: writeFile.C:253
DMDModel.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::DMDModels::STDMD
Streaming Total Dynamic Mode Decomposition (i.e. STDMD) is a variant of dynamic mode decomposition.
Definition: STDMD.H:306
SquareMatrix.H
Foam::DMDModels::STDMD::STDMD
STDMD(const fvMesh &mesh, const word &name, const dictionary &dict)
Construct from components.
Definition: STDMD.C:737
RectangularMatrix.H