Proper Orthogonal Decomposition

Description

The POD treatment performs a Proper Orthogonal Decomposition (POD) for every zone within the base.

Warning

A zone must contain all data to be decomposed.

So here, a zone is not necessarily a part of a mesh (i.e.: from a structured multibloc mesh). If you want to apply the treatment on the whole mesh, then you have to merge all domains in one zone.

Warning

The data contained in the zone will be concatenated to form a new vector. The treatment is performed on this vector. e.g.: if the instant contains ‘p’ and ‘t’, the new vector will be [‘p’, ‘t’]

Definition of the POD

The Proper Orthogonal Decomposition (POD) is a technique used to decompose a matrix and characterize it by its principal components which are called modes [Chatterjee2000]. To approximate a function \(z(x,t)\), only a finite sum is required:

\[z(x,t) \approx \sum_{k=1}^{m} a_k(t) \phi_k(x).\]

The basis function \(\phi_{k}(x)\) can be chosen among Fourier series or Chebyshev polynomials, etc. For a chosen basis of functions, a set of unique time-functions \(a_k(t)\) arises. In case of the POD, the basis function are chosen orthonormal, meaning that:

\[\begin{split}\int_{x} \phi_{k_1} \phi_{k_2} dx = \left\{\begin{array}{rcl} 1 & \text{if} & k_1 = k_2 \\ 0 & \text{if} & k_1 \neq k_2\end{array} \right. ,\ \alpha_k (t) = \int_{x} z(x,t) \phi_k(x) dx.\end{split}\]

The principe of the POD is to chose \(\phi_k(x)\) such that the approximation of \(z(x,t)\) is the best in a least-squares sense. These orthonormal functions are called the proper orthogonal modes of the function.

When dealing with CFD simulations, the number of modes \(m\) is usually smaller than the number of measures (or snapshots) \(n\). Hence, from the existing decomposition methods, the Singular Value Decomposition (SVD) is used. It is the snapshots methods [Cordier2006].

The Singular Value Decomposition (SVD) is a factorization operation of a matrix expressed as:

\[A = U \Sigma V^T,\]

with \(V\) diagonalizes \(A^TA\), \(U\) diagonalizes \(AA^T\) and \(\Sigma\) is the singular value matrix which diagonal is composed by the singular values of \(A\). Knowing that a singular value is the square root of an eigenvector. \(u_i\) and \(v_i\) are eigenvectors of respectively \(U\) and \(V\). form orthonormal basis. Thus, the initial matrix can be rewritten:

\[A = \sum_{i=1}^{r} \sigma_i u_i v_i^T,\]

\(r\) being the rank of the matrix. If taken \(k < r\), an approximation of the initial matrix can be constructed. This allows to compress the data as only an extract of \(u\) and \(v\) need to be stored.

Parameters

  • base: Base

    The input base that contains data.

  • tolerance: float

    Tolerance for basis modes filtering [0->1].

  • dim_max: int

    Maximum number of basis modes.

  • POD_vectors: bool

    Output POD vectors as base attributes.

  • variables: list(str)

    The variable names to take into account.

Preconditions

The input base must be an unstructured mesh.

Postconditions

A base with the same number and names of zones as the input base.

Each zone contains 2 or 4 instants depending on the value of the parameter POD_vectors.

The first instant contains the modes, the second the mean of the snapshots, the third the left-singular vectors, and the fourth the right-singular vectors.

References

A. Chatterjee: An introduction to the proper orthogonal decomposition. Current Science 78.7. 2000

L. Cordier: Reduction de dynamique par decomposition orthogonale aux valeurs propres (PDO). Ecole de printemps. 2006

Example

import antares
treatment = antares.Treatment('pod')
treatment['variables'] = ['pressure']

Main functions

class antares.treatment.TreatmentPOD.TreatmentPOD
execute()

Perform the POD treatment.

Returns:

A base with the same number and the same names of zones as the input base.

Each zone contains 2 or 4 instants depending on the value of the parameter POD_vectors.

The first instant contains the modes, the second the mean of the snapshots, the third the left-singular vector, and the fourth the right-singular vectors.

Return type:

Base

S

Singular values matrix, ndarray(nb of modes, nb of snapshots), only the diagonal is stored, of length (nb of modes).

U

Columns of matrix U are left-singular vectors of A, ndarray(nb of snapshots, nb of modes).

VT

Columns of matrix V are right-singular vectors of A, ndarray(nb of snapshots, nb of snapshots), after filtering (nb of snapshots, nb of modes).

Example

import os
if not os.path.isdir('OUTPUT'):
    os.makedirs('OUTPUT')


from antares import Reader, Treatment

# ------------------
# Read the files
# ------------------
r = Reader('fmt_tp')
r['filename'] = os.path.join('..', 'data', 'SNAPSHOTS_POD', '<instant>', 'pressure_<zone>.dat')
base = r.read()

# ------------------
# View the data
# ------------------
print('Reading base: ', base)
print('Zone ex: ', base['ex'])
print('Instant ex -> 0: ', base['ex']['0'])
print('Value ex -> 0 -> p', base['ex']['0']['p'])

# ------------------
# Setup
# ------------------
treatment = Treatment('POD')
treatment['base'] = base
treatment['tolerance'] = 0.99
treatment['dim_max'] = 100
treatment['POD_vectors'] = True
treatment['variables'] = ['p']

# ------------------
# Compute the POD
# ------------------
result = treatment.execute()

# ------------------
# Get some results
# ------------------
print("POD modes: ", result[0]['modes'][0])
print("POD modes: ", result[1]['modes'][0])
print("POD parameters: ")
for k, v in result.attrs.items():
    print(k, v)
# print("POD vectors: ", result.attrs['POD_vectors'])