STDMD Class Reference

Streaming Total Dynamic Mode Decomposition (i.e. STDMD) is a variant of dynamic mode decomposition. More...

Inheritance diagram for STDMD:
[legend]
Collaboration diagram for STDMD:
[legend]

Public Member Functions

 TypeName ("STDMD")
 Runtime type information. More...
 
 STDMD (const fvMesh &mesh, const word &name, const dictionary &dict)
 Construct from components. More...
 
 STDMD (const STDMD &)=delete
 No copy construct. More...
 
void operator= (const STDMD &)=delete
 No copy assignment. More...
 
virtual ~STDMD ()=default
 Destructor. More...
 
virtual bool initialise (const RMatrix &z)
 Initialise 'Q' and 'G' (both require the first two snapshots) More...
 
virtual bool update (const RMatrix &z)
 Incremental orthonormal basis update (K:Fig. 15) More...
 
virtual bool fit ()
 
virtual bool read (const dictionary &dict)
 Read STDMD settings. More...
 
- Public Member Functions inherited from DMDModel
 TypeName ("DMDModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, DMDModel, dictionary,(const fvMesh &mesh, const word &name, const dictionary &dict),(mesh, name, dict))
 
 DMDModel (const fvMesh &mesh, const word &name, const dictionary &dict)
 Construct from components. More...
 
 DMDModel (const DMDModel &)=delete
 No copy construct. More...
 
void operator= (const DMDModel &)=delete
 No copy assignment. More...
 
virtual ~DMDModel ()=default
 Destructor. More...
 
virtual bool initialise (const RMatrix &snapshot)=0
 Initialise model data members with a given snapshot. More...
 
virtual bool update (const RMatrix &snapshot)=0
 Update model data members with a given snapshot. More...
 
virtual bool fit ()=0
 
virtual void reconstruct (const wordList modes)
 
virtual bool read (const dictionary &dict)=0
 Read model settings. More...
 
- Public Member Functions inherited from writeFile
 writeFile (const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true)
 Construct from objectRegistry, prefix, fileName. More...
 
 writeFile (const objectRegistry &obr, const fileName &prefix, const word &name, const dictionary &dict, const bool writeToFile=true)
 
 writeFile (const writeFile &wf)
 Construct copy. More...
 
virtual ~writeFile ()=default
 Destructor. More...
 
virtual bool read (const dictionary &dict)
 Read. More...
 
virtual OFstreamfile ()
 Return access to the file (if only 1) More...
 
virtual bool writeToFile () const
 Flag to allow writing to file. More...
 
virtual bool canWriteHeader () const
 Flag to allow writing the header. More...
 
virtual label charWidth () const
 Return width of character stream output. More...
 
virtual void writeCommented (Ostream &os, const string &str) const
 Write a commented string to stream. More...
 
virtual void writeTabbed (Ostream &os, const string &str) const
 Write a tabbed string to stream. More...
 
virtual void writeHeader (Ostream &os, const string &str) const
 Write a commented header to stream. More...
 
virtual void writeCurrentTime (Ostream &os) const
 Write the current time to stream. More...
 
virtual void writeBreak (Ostream &os) const
 Write a break marker to the stream. More...
 
template<class Type >
void writeHeaderValue (Ostream &os, const string &property, const Type &value) const
 Write a (commented) header property and value pair. More...
 
template<class Type >
void writeValue (Ostream &os, const Type &val) const
 Write a given value to stream with the space delimiter. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from DMDModel
static autoPtr< DMDModelNew (const fvMesh &mesh, const word &name, const dictionary &dict)
 Return a reference to the selected DMD model. More...
 
- Static Public Attributes inherited from writeFile
static label addChars = 8
 Additional characters for writing. More...
 
virtual bool dynamics ()=0
 Compute and write mode dynamics. More...
 
virtual bool modes ()=0
 Compute and write modes. More...
 
- Protected Member Functions inherited from writeFile
void initStream (Ostream &os) const
 Initialise the output stream for writing. More...
 
fileName baseFileDir () const
 Return the base directory for output. More...
 
fileName baseTimeDir () const
 Return the base directory for the current time value. More...
 
virtual autoPtr< OFstreamcreateFile (const word &name, scalar timeValue) const
 Return autoPtr to a new file for a given time. More...
 
virtual autoPtr< OFstreamcreateFile (const word &name) const
 Return autoPtr to a new file using the simulation start time. More...
 
virtual void resetFile (const word &name)
 Reset internal file pointer to new file with new name. More...
 
Omanip< int > valueWidth (const label offset=0) const
 Return the value width when writing to stream with optional offset. More...
 
void operator= (const writeFile &)=delete
 No copy assignment. More...
 
- Protected Attributes inherited from DMDModel
const fvMeshmesh_
 Reference to the mesh. More...
 
const word name_
 Name of operand function object. More...
 
- Protected Attributes inherited from writeFile
const objectRegistryfileObr_
 Reference to the region objectRegistry. More...
 
const fileName prefix_
 Prefix. More...
 
word fileName_
 Name of file. More...
 
autoPtr< OFstreamfilePtr_
 File pointer. More...
 
label writePrecision_
 Write precision. More...
 
bool writeToFile_
 Flag to enable/disable writing to file. More...
 
bool updateHeader_
 
bool writtenHeader_
 Flag to identify whether the header has been written. More...
 
bool useUserTime_
 
scalar startTime_
 Start time value. More...
 

Detailed Description

Streaming Total Dynamic Mode Decomposition (i.e. STDMD) is a variant of dynamic mode decomposition.

Among other Dynamic Mode Decomposition (DMD) variants, STDMD is presumed to provide the general DMD method capabilities alongside economised and feasible memory and CPU usage.

The code implementation corresponds to Figs. 15-16 of the first citation below, more broadly to Section 2.4.

References:

    DMD and mode-sorting algorithms (tags:K, HRDC, KZ, HWR):
        Kiewat, M. (2019).
        Streaming modal decomposition approaches for vehicle aerodynamics.
        PhD thesis. Munich: Technical University of Munich.
        URL:mediatum.ub.tum.de/doc/1482652/1482652.pdf

        Hemati, M. S., Rowley, C. W.,
        Deem, E. A., & Cattafesta, L. N. (2017).
        De-biasing the dynamic mode decomposition
        for applied Koopman spectral analysis of noisy datasets.
        Theoretical and Computational Fluid Dynamics, 31(4), 349-368.
        DOI:10.1007/s00162-017-0432-2

        Kou, J., & Zhang, W. (2017).
        An improved criterion to select
        dominant modes from dynamic mode decomposition.
        European Journal of Mechanics-B/Fluids, 62, 109-129.
        DOI:10.1016/j.euromechflu.2016.11.015

        Hemati, M. S., Williams, M. O., & Rowley, C. W. (2014).
        Dynamic mode decomposition for large and streaming datasets.
        Physics of Fluids, 26(11), 111701.
        DOI:10.1063/1.4901016

    Parallel classical Gram-Schmidt process (tag:Ka):
        Katagiri, T. (2003).
        Performance evaluation of parallel
        Gram-Schmidt re-orthogonalization methods.
        In: Palma J. M. L. M., Sousa A. A., Dongarra J., Hernández V. (eds)
        High Performance Computing for Computational Science — VECPAR 2002.
        Lecture Notes in Computer Science, vol 2565, p. 302-314.
        Berlin, Heidelberg: Springer.
        DOI:10.1007/3-540-36569-9_19

    Parallel direct tall-skinny QR decomposition (tags:BGD, DGHL):
        Benson, A. R., Gleich, D. F., & Demmel, J. (2013).
        Direct QR factorizations for
        tall-and-skinny matrices in MapReduce architectures.
        2013 IEEE International Conference on Big Data.
        DOI:10.1109/bigdata.2013.6691583

        Demmel, J., Grigori, L., Hoemmen, M., & Langou, J. (2012).
        Communication-optimal parallel
        and sequential QR and LU factorizations.
        SIAM Journal on Scientific Computing, 34(1), A206-A239.
        DOI:10.1137/080731992
Usage
Minimal example by using system/controlDict.functions:
DMD1
{
    // Mandatory entries
    DMDModel            STDMD;

    // Conditional mandatory entries

        // Option-1
        interval            5.5;

        // Option-2
        executeInterval     10;

    // Optional entries
    modeSorter          kiewat;
    nGramSchmidt        5;
    maxRank             50;
    nModes              50;
    fMin                0;
    fMax                1000000000;
    nAgglomerationProcs 20;

    // Optional entries (not recommended to change)
    minBasis            0.00000001;
    minEVal             0.00000001;
    sortLimiter         500.0;

    // Inherited entries
    ...
}

where the entries mean:

Property Description Type Reqd Deflt
DMDModel Type: STDMD word yes -
interval STDMD time-step size [s] scalar

cndtnl

executeInterval*(current time-step of the simulation)
modeSorter Mode-sorting algorithm model word

no

firstSnapshot
nModes

Number of output modes in input frequency range

label no GREAT
maxRank Max columns in accumulated matrices label no GREAT
nGramSchmidt Number of Gram-Schmidt iterations label no 5
fMin Min (non-negative) output frequency label no 0
fMax Max output frequency label no GREAT
nAgglomerationProcs

Number of processors at each agglomeration

unit during the computation of reduced Koopman

operator

label no 20
minBasis Orthogonal basis expansion threshold scalar no 1e-8
minEVal

Min eigenvalue for below eigenvalues are omitted

scalar no 1e-8
sortLimiter

Max allowable magnitude for

mag(eigenvalue)*(number of DMD steps) to be used in

modeSorter={kiewat,kouZhang} to avoid overflow errors

scalar no 500.0

Options for the modeSorter entry:

      kiewat        | Modified weighted-amplitude scaling method
      kouZhang      | Weighted-amplitude scaling method
      firstSnapshot | First-snapshot amplitude magnitude method
Note
  • To specify the STDMD time-step size (not necessarily equal to the time step of the simulation), entries of either interval or executeInterval must be available (highest to lowest precedence). While interval allows users to directly specify the STDMD time-step size in seconds, in absence of interval, for convenience, executeInterval allows users to compute the STDMD time-step internally by multiplying itself with the current time-step size of the simulation.
  • Limitation: Restart is currently not available since intermediate writing of STDMD matrices are not supported.
  • Limitation: Non-physical input (e.g. full-zero fields) can upset STDMD.
  • Warning: In modeSorter algorithms of kiewat and kouZhang, a function of \(power(x,y)\) exists where x is the magnitude of an eigenvalue, and y is the number of STDMD snapshots. This function poses a risk of overflow errors since, for example, if x ends up above 1.5 with 500 snapshots or more, this function automatically throws the floating point error meaning overflow. Therefore, the domain-range of this function is heuristically constrained by the optional entry sortLimiter which skips the evaluation of this function for a given mag(eigenvalue)*(no. of STDMD steps), i.e. x*y, whose magnitude is larger than sortLimiter.
See also
Source files

Definition at line 303 of file STDMD.H.

Constructor & Destructor Documentation

◆ STDMD() [1/2]

STDMD ( const fvMesh mesh,
const word name,
const dictionary dict 
)

Construct from components.

Definition at line 808 of file STDMD.C.

◆ STDMD() [2/2]

STDMD ( const STDMD )
delete

No copy construct.

◆ ~STDMD()

virtual ~STDMD ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "STDMD"  )

Runtime type information.

◆ operator=()

void operator= ( const STDMD )
delete

No copy assignment.

◆ initialise()

bool initialise ( const RMatrix z)
virtual

Initialise 'Q' and 'G' (both require the first two snapshots)

Implements DMDModel.

Definition at line 951 of file STDMD.C.

References UList< T >::begin(), Matrix< Form, Type >::cbegin(), Foam::fileHandler(), Matrix< Form, Type >::m(), Foam::mag(), IOobject::NO_READ, IOobject::NO_WRITE, Foam::sqr(), and fileOperation::writeObject().

Here is the call graph for this function:

◆ update()

bool update ( const RMatrix z)
virtual

Incremental orthonormal basis update (K:Fig. 15)

Working copy of the augmented snapshot matrix "z" being used in the classical Gram-Schmidt process

Implements DMDModel.

Definition at line 1012 of file STDMD.C.

◆ fit()

bool fit ( )
virtual

Compute and write modes and mode dynamics of model data members

Implements DMDModel.

Definition at line 1044 of file STDMD.C.

References splitCell::master(), and Foam::max().

Here is the call graph for this function:

◆ read()

bool read ( const dictionary dict)
virtual

Read STDMD settings.

Implements DMDModel.

Definition at line 864 of file STDMD.C.

References dict, e, Foam::endl(), MinMax< scalar >::ge(), MinMax< T >::ge(), Foam::Info, Foam::nl, and Foam::tab.

Here is the call graph for this function:

The documentation for this class was generated from the following files: