DMD

class modred.dmd.DMDHandles(inner_product=None, get_array=<function load_array_text>, put_array=<function save_array_text>, max_vecs_per_node=None, verbosity=1)[source]

Dynamic Mode Decomposition implemented for large datasets.

Kwargs:

inner_product: Function that computes inner product of two vector objects.

put_array: Function to put an array out of modred, e.g., write it to file.

get_array: Function to get an array into modred, e.g., load it from file.

max_vecs_per_node: Maximum number of vectors that can be stored in memory, per node.

verbosity: 1 prints progress and warnings, 0 prints almost nothing.

Computes DMD modes from vector objects (or handles). It uses vectorspace.VectorSpaceHandles for low level functions.

Usage:

myDMD = DMDHandles(inner_product=my_inner_product)
myDMD.compute_decomp(vec_handles)
myDMD.compute_exact_modes(range(50), mode_handles)

See also compute_DMD_arrays_snaps_method(), compute_DMD_arrays_direct_method(), and vectors.

compute_adjoint_modes(mode_indices, mode_handles, vec_handles=None)[source]

Computes adjoint DMD modes and calls put on them using mode handles.

Args:

mode_indices: List of indices describing which projected modes to compute, e.g. range(10) or [3, 0, 5].

mode_handles: List of handles for projected modes to compute.

Kwargs:
vec_handles: List of handles for vector objects. Optional if when calling compute_decomp().
compute_decomp(vec_handles, adv_vec_handles=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes eigendecomposition of low-order linear map approximating relationship between vector objects, returning various arrays necessary for computing and characterizing DMD modes.

Args:
vec_handles: List of handles for vector objects.
Kwargs:

adv_vec_handles: List of handles for vector objects advanced in time. If not provided, it is assumed that the vector objects describe a sequential time-series. Thus vec_handles becomes vec_handles[:-1] and adv_vec_handles becomes vec_handles[1:].

atol: Level below which DMD eigenvalues are truncated.

rtol: Maximum relative difference between largest and smallest DMD eigenvalues. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

eigvals: 1D array of eigenvalues of low-order linear map, i.e., the DMD eigenvalues.

R_low_order_eigvecs: Array whose columns are right eigenvectors of approximating low-order linear map.

L_low_order_eigvecs: Array whose columns are left eigenvectors of approximating low-order linear map.

correlation_array_eigvals: 1D array of eigenvalues of correlation array.

correlation_array_eigvecs: Array whose columns are eigenvectors of correlation array.

compute_eigendecomp(atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes eigendecompositions of correlation array and approximating low-order linear map.

Kwargs:

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Useful if you already have the correlation array and cross-correlation array and want to avoid recomputing them.

Usage:

DMD.correlation_array = pre_existing_correlation_array
DMD.cross_correlation_array = pre_existing_cross_correlation_array
DMD.compute_eigendecomp()
DMD.compute_exact_modes(
    mode_idx_list, mode_handles, adv_vec_handles=adv_vec_handles)

Another way to use this is to compute a DMD using a truncated basis for the projection of the approximating linear map. Start by either computing a full decomposition or by loading pre-computed correlation and cross-correlation arrays.

Usage:

# Start with a full decomposition
DMD_eigvals, correlation_array_eigvals = DMD.compute_decomp(
    vec_handles)[0, 3]

# Do some processing to determine the truncation level, maybe based
# on the DMD eigenvalues and correlation array eigenvalues
desired_num_eigvals = my_post_processing_func(
    DMD_eigvals, correlation_array_eigvals)

# Do a truncated decomposition
DMD_eigvals_trunc = DMD.compute_eigendecomp(
  max_num_eigvals=desired_num_eigvals)

# Compute modes for truncated decomposition
DMD.compute_exact_modes(
    mode_idx_list, mode_handles, adv_vec_handles=adv_vec_handles)

Since it doesn’t overwrite the correlation and cross-correlation arrays, compute_eigendecomp can be called many times in a row to do computations for different truncation levels. However, the results of the decomposition (e.g., self.eigvals) do get overwritten, so you may want to call a put method to save those results somehow.

compute_exact_modes(mode_indices, mode_handles, adv_vec_handles=None)[source]

Computes exact DMD modes and calls put on them using mode handles.

Args:

mode_indices: List of indices describing which exact modes to compute, e.g. range(10) or [3, 0, 5].

mode_handles: List of handles for exact modes to compute.

Kwargs:
vec_handles: List of handles for vector objects. Optional if when calling compute_decomp().
compute_proj_coeffs()[source]

Computes projection of vector objects onto DMD modes. Note that a biorthogonal projection onto exact DMD modes is analytically equivalent to a least-squares projection onto projected DMD modes.

Returns:

proj_coeffs: Array of projection coefficients for vector objects, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes.

adv_proj_coeffs: Array of projection coefficients for vector objects advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes.

compute_proj_modes(mode_indices, mode_handles, vec_handles=None)[source]

Computes projected DMD modes and calls put on them using mode handles.

Args:

mode_indices: List of indices describing which projected modes to compute, e.g. range(10) or [3, 0, 5].

mode_handles: List of handles for projected modes to compute.

Kwargs:
vec_handles: List of handles for vector objects. Optional if when calling compute_decomp().
compute_spectrum()[source]

Computes DMD spectral coefficients. These coefficients come from a biorthogonal projection of the first vector object onto the exact DMD modes, which is analytically equivalent to doing a least-squares projection onto the projected DMD modes.

Returns:
spectral_coeffs: 1D array of DMD spectral coefficients, calculated as the magnitudes of the projection coefficients of first data vector. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
get_correlation_array(src)[source]

Gets the correlation array from source (memory or file).

Args:
src: Source from which to retrieve correlation array.
get_cross_correlation_array(src)[source]

Gets the cross-correlation array from source (memory or file).

Args:
src: Source from which to retrieve cross-correlation array.
get_decomp(eigvals_src, R_low_order_eigvecs_src, L_low_order_eigvecs_src, correlation_array_eigvals_src, correlation_array_eigvecs_src)[source]

Gets the decomposition arrays from sources (memory or file).

Args:

eigvals_src: Source from which to retrieve eigenvalues of approximating low-order linear map (DMD eigenvalues).

R_low_order_eigvecs_src: Source from which to retrieve right eigenvectors of approximating low-order linear DMD map.

L_low_order_eigvecs_src: Source from which to retrieve left eigenvectors of approximating low-order linear DMD map.

correlation_array_eigvals_src: Source from which to retrieve eigenvalues of correlation array.

correlation_array_eigvecs_src: Source from which to retrieve eigenvectors of correlation array.

get_proj_coeffs(src, adv_src)[source]

Gets the projection coefficients and advanced projection coefficients from sources (memory or file).

Args:

src: Source from which to retrieve projection coefficients.

adv_src: Source from which to retrieve advanced projection coefficients.

get_spectral_coeffs(src)[source]

Gets the spectral coefficients from source (memory or file).

Args:
src: Source from which to retrieve spectral coefficients.
put_L_low_order_eigvecs(dest)[source]

Puts left eigenvectors of approximating low-order linear map to dest.

put_R_low_order_eigvecs(dest)[source]

Puts right eigenvectors of approximating low-order linear map to dest.

put_correlation_array(dest)[source]

Puts correlation array to dest.

put_correlation_array_eigvals(dest)[source]

Puts eigenvalues of correlation array to dest.

put_correlation_array_eigvecs(dest)[source]

Puts eigenvectors of correlation array to dest.

put_cross_correlation_array(dest)[source]

Puts cross-correlation array to dest.

put_decomp(eigvals_dest, R_low_order_eigvecs_dest, L_low_order_eigvecs_dest, correlation_array_eigvals_dest, correlation_array_eigvecs_dest)[source]

Puts the decomposition arrays in destinations (memory or file).

Args:

eigvals_dest: Destination in which to put eigenvalues of approximating low-order linear map (DMD eigenvalues).

R_low_order_eigvecs_dest: Destination in which to put right eigenvectors of approximating low-order linear map.

L_low_order_eigvecs_dest: Destination in which to put left eigenvectors of approximating low-order linear map.

correlation_array_eigvals_dest: Destination in which to put eigenvalues of correlation array.

correlation_array_eigvecs_dest: Destination in which to put eigenvectors of correlation array.

put_eigvals(dest)[source]

Puts eigenvalues of approximating low-order-linear map (DMD eigenvalues) to dest.

put_proj_coeffs(dest, adv_dest)[source]

Puts projection coefficients to dest, advanced projection coefficients to adv_dest.

put_spectral_coeffs(dest)[source]

Puts DMD spectral coefficients to dest.

sanity_check(test_vec_handle)[source]

Checks that user-supplied vector handle and vector satisfy requirements.

Args:
test_vec_handle: A vector handle to test.

See vectorspace.VectorSpaceHandles.sanity_check().

class modred.dmd.TLSqrDMDHandles(inner_product=None, get_array=<function load_array_text>, put_array=<function save_array_text>, max_vecs_per_node=None, verbosity=1)[source]

Total Least Squares Dynamic Mode Decomposition implemented for large datasets.

Kwargs:

inner_product: Function that computes inner product of two vector objects.

put_array: Function to put an array out of modred, e.g., write it to file.

get_array: Function to get an array into modred, e.g., load it from file.

max_vecs_per_node: Maximum number of vectors that can be stored in memory, per node.

verbosity: 1 prints progress and warnings, 0 prints almost nothing.

Computes Total-Least-Squares DMD modes from vector objects (or handles). It uses vectorspace.VectorSpaceHandles for low level functions.

Usage:

myDMD = TLSqrDMDHandles(inner_product=my_inner_product)
myDMD.compute_decomp(vec_handles, max_num_eigvals=100)
myDMD.compute_exact_modes(range(50), mode_handles)

The total-least-squares variant of DMD accounts for an asymmetric treatment of the non-time-advanced and time-advanced data in standard DMD, in which the former is assumed to be perfect information, whereas all noise is limited to the latter. In many cases, since the non-time-advanced and advanced data come from the same source, noise could be present in either.

Note that max_num_eigvals must be set to a value smaller than the rank of the dataset. In other words, if the projection basis for total-least-squares DMD is not truncated, then the algorithm reduces to standard DMD. For over-constrained datasets (number of vector objects is larger than the size of each vector object), this occurs naturally. For under-constrained datasets, (number of vector objects is smaller than size of vector objects), this must be done explicitly by the user. At this time, there is no standard method for choosing a truncation level. One approach is to look at the roll-off of the correlation array eigenvalues, which contains information about the “energy” content of each projection basis vectors.

Also, note that TLSqrDMDHandles inherits from DMDHandles, so certain methods are available, even though they are not implemented/documented here (namely several put functions).

See also compute_TLSqrDMD_arrays_snaps_method(), compute_TLSqrDMD_arrays_direct_method(), and vectors.

compute_decomp(vec_handles, adv_vec_handles=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes eigendecomposition of low-order linear map approximating relationship between vector objects, returning various arrays necessary for computing and characterizing DMD modes.

Args:
vec_handles: List of handles for vector objects.
Kwargs:

adv_vec_handles: List of handles for vector objects advanced in time. If not provided, it is assumed that the vector objects describe a sequential time-series. Thus vec_handles becomes vec_handles[:-1] and adv_vec_handles becomes vec_handles[1:].

atol: Level below which DMD eigenvalues are truncated.

rtol: Maximum relative difference between largest and smallest DMD eigenvalues. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

eigvals: 1D array of eigenvalues of low-order linear map, i.e., the DMD eigenvalues.

R_low_order_eigvecs: Array whose columns are right eigenvectors of approximating low-order linear map.

L_low_order_eigvecs: Array whose columns are left eigenvectors of approximating low-order linear map.

sum_correlation_array_eigvals: 1D array of eigenvalues of sum correlation array.

sum_correlation_array_eigvecs: Array whose columns are eigenvectors of sum correlation array.

proj_correlation_array_eigvals: 1D array of eigenvalues of projected correlation array.

proj_correlation_array_eigvecs: Array whose columns are eigenvectors of projected correlation array.

Note that the truncation level (corresponding to max_num_eigvals) must be set to a value smaller than the rank of the dataset, otherwise total-least-squares DMD reduces to standard DMD. This occurs naturally for over-constrained datasets, but must be enforced by the user for under-constrined datasets.

compute_eigendecomp(atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes eigendecompositions of correlation array and approximating low-order linear map.

Kwargs:

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Useful if you already have the correlation array and cross-correlation array and want to avoid recomputing them.

Usage:

TLSqrDMD.correlation_array = pre_existing_correlation_array
TLSqrDMD.adv_correlation_array = pre_existing_adv_correlation_array
TLSqrDMD.cross_correlation_array =
    pre_existing_cross_correlation_array
TLSqrDMD.compute_eigendecomp()
TLSqrDMD.compute_exact_modes(
    mode_idx_list, mode_handles, adv_vec_handles=adv_vec_handles)

Another way to use this is to compute TLSqrDMD using different basis truncation levels for the projection of the approximating linear map. Start by either computing a full decomposition or by loading pre-computed correlation and cross-correlation arrays.

Usage:

# Start with a full decomposition
DMD_eigvals, correlation_array_eigvals = TLSqrDMD.compute_decomp(
    vec_handles)[0, 3]

# Do some processing to determine the truncation level, maybe based
# on the DMD eigenvalues and correlation array eigenvalues
desired_num_eigvals = my_post_processing_func(
    DMD_eigvals, correlation_array_eigvals)

# Do a truncated decomposition
DMD_eigvals_trunc = TLSqrDMD.compute_eigendecomp(
  max_num_eigvals=desired_num_eigvals)

# Compute modes for truncated decomposition
TLSqrDMD.compute_exact_modes(
    mode_idx_list, mode_handles, adv_vec_handles=adv_vec_handles)

Since it doesn’t overwrite the correlation and cross-correlation arrays, compute_eigendecomp can be called many times in a row to do computations for different truncation levels. However, the results of the decomposition (e.g., self.eigvals) do get overwritten, so you may want to call a put method to save those results somehow.

Note that the truncation level (corresponding to max_num_eigvals) must be set to a value smaller than the rank of the dataset, otherwise total-least-squares DMD reduces to standard DMD. This occurs naturally for over-constrained datasets, but must be enforced by the user for under-constrined datasets.

compute_proj_coeffs()[source]

Computes projection of (de-noised) vector objects onto DMD modes. Note that a biorthogonal projection onto exact DMD modes is analytically equivalent to a least-squares projection onto projected DMD modes.

Returns:

proj_coeffs: Array of projection coefficients for (de-noised)vector objects, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes.

adv_proj_coeffs: Array of projection coefficients for (de-noised) vector objects advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes.

compute_spectrum()[source]

Computes DMD spectral coefficients. These coefficients come from a biorthogonal projection of the first (de-noised) vector object onto the exact DMD modes, which is analytically equivalent to doing a least-squares projection onto the projected DMD modes.

Returns:
spectral_coeffs: 1D array of DMD spectral coefficients.
get_adv_correlation_array(src)[source]

Gets the advanced correlation array from source (memory or file).

Args:
src: Source from which to retrieve advanced correlation array.
get_decomp(eigvals_src, R_low_order_eigvecs_src, L_low_order_eigvecs_src, sum_correlation_array_eigvals_src, sum_correlation_array_eigvecs_src, proj_correlation_array_eigvals_src, proj_correlation_array_eigvecs_src)[source]

Gets the decomposition arrays from sources (memory or file).

Args:

eigvals_src: Source from which to retrieve eigenvalues of approximating low-order linear map (DMD eigenvalues).

R_low_order_eigvecs_src: Source from which to retrieve right eigenvectors of approximating low-order linear DMD map.

L_low_order_eigvecs_src: Source from which to retrieve left eigenvectors of approximating low-order linear DMD map.

sum_correlation_array_eigvals_src: Source from which to retrieve eigenvalues of sum correlation array.

sum_correlation_array_eigvecs_src: Source from which to retrieve eigenvectors of sum correlation array.

proj_correlation_array_eigvals_src: Source from which to retrieve eigenvalues of projected correlation array.

proj_correlation_array_eigvecs_src: Source from which to retrieve eigenvectors of projected correlation array.

get_proj_correlation_array(src)[source]

Gets the projected correlation array from source (memory or file).

Args:
src: Source from which to retrieve projected correlation array.
get_sum_correlation_array(src)[source]

Gets the sum correlation array from source (memory or file).

Args:
src: Source from which to retrieve sum correlation array.
put_adv_correlation_array(dest)[source]

Puts advanced correlation array to dest.

put_correlation_array_eigvals(dest)[source]

This method is not available for total least squares DMD

put_correlation_array_eigvecs(dest)[source]

This method is not available for total least squares DMD

put_decomp(eigvals_dest, R_low_order_eigvecs_dest, L_low_order_eigvecs_dest, sum_correlation_array_eigvals_dest, sum_correlation_array_eigvecs_dest, proj_correlation_array_eigvals_dest, proj_correlation_array_eigvecs_dest)[source]

Puts the decomposition arrays in destinations (file or memory).

Args:

eigvals_dest: Destination in which to put eigenvalues of approximating low-order linear map (DMD eigenvalues).

R_low_order_eigvecs_dest: Destination in which to put right eigenvectors of approximating low-order linear map.

L_low_order_eigvecs_dest: Destination in which to put left eigenvectors of approximating low-order linear map.

sum_correlation_array_eigvals_dest: Destination in which to put eigenvalues of sum correlation array.

sum_correlation_array_eigvecs_dest: Destination in which to put eigenvectors of sum correlation array.

proj_correlation_array_eigvals_dest: Destination in which to put eigenvalues of projected correlation array.

proj_correlation_array_eigvecs_dest: Destination in which to put eigenvectors of projected correlation array.

put_proj_correlation_array(dest)[source]

Puts projected correlation array to dest.

put_proj_correlation_array_eigvals(dest)[source]

Puts eigenvalues of projected correlation array to dest.

put_proj_correlation_array_eigvecs(dest)[source]

Puts eigenvectors of projected correlation array to dest.

put_sum_correlation_array(dest)[source]

Puts sum correlation array to dest.

put_sum_correlation_array_eigvals(dest)[source]

Puts eigenvalues of sum correlation array to dest.

put_sum_correlation_array_eigvecs(dest)[source]

Puts eigenvectors of sum correlation array to dest.

modred.dmd.compute_DMD_arrays_direct_method(vecs, adv_vecs=None, mode_indices=None, inner_product_weights=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes DMD modes using data stored in arrays, using direct method.

Args:
vecs: Array whose columns are data vectors.
Kwargs:

adv_vecs: Array whose columns are data vectors advanced in time. If not provided, then it is assumed that the vectors describe a sequential time-series. Thus vecs becomes vecs[:, :-1] and adv_vecs becomes vecs[:, 1:].

mode_indices: List of indices describing which modes to compute. Examples are range(10) or [3, 0, 6, 8]. If no mode indices are specified, then all modes will be computed.

inner_product_weights: 1D or 2D array of inner product weights. Corresponds to W in inner product v_1^* W v_2.

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

res: Results of DMD computation, stored in a namedtuple with the following attributes:

  • eigvals: 1D array of eigenvalues of approximating low-order linear map (DMD eigenvalues).
  • spectral_coeffs: 1D array of DMD spectral coefficients, calculated as the magnitudes of the projection coefficients of first data vector. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • exact_modes: Array whose columns are exact DMD modes.
  • proj_modes: Array whose columns are projected DMD modes.
  • adjoint_modes: Array whose columns are adjoint DMD modes.
  • proj_coeffs: Array of projection coefficients for data vectors expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • adv_proj_coeffs: Array of projection coefficients for data vectors advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • R_low_order_eigvecs: Array of right eigenvectors of approximating low-order linear map.
  • L_low_order_eigvecs: Array of left eigenvectors of approximating low-order linear map.
  • correlation_array_eigvals: 1D array of eigenvalues of correlation array.
  • correlation_array_eigvecs: Array of eigenvectors of correlation array.

Attributes can be accessed using calls like res.exact_modes. To see all available attributes, use print(res).

This method does not square the array of vectors as in the method of snapshots (compute_DMD_arrays_snaps_method()). It’s slightly more accurate, but slower when the number of elements in a vector is more than the number of vectors, i.e., when vecs has more rows than columns.

modred.dmd.compute_DMD_arrays_snaps_method(vecs, adv_vecs=None, mode_indices=None, inner_product_weights=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes DMD modes using data stored in arrays, using method of snapshots.

Args:
vecs: Array whose columns are data vectors.
Kwargs:

adv_vecs: Array whose columns are data vectors advanced in time. If not provided, then it is assumed that the vectors describe a sequential time-series. Thus vecs becomes vecs[:, :-1] and adv_vecs becomes vecs[:, 1:].

mode_indices: List of indices describing which modes to compute. Examples are range(10) or [3, 0, 6, 8]. If no mode indices are specified, then all modes will be computed.

inner_product_weights: 1D or 2D array of inner product weights. Corresponds to W in inner product v_1^* W v_2.

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

res: Results of DMD computation, stored in a namedtuple with the following attributes:

  • eigvals: 1D array of eigenvalues of approximating low-order linear map (DMD eigenvalues).
  • spectral_coeffs: 1D array of DMD spectral coefficients, calculated as the magnitudes of the projection coefficients of first data vector. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • exact_modes: Array whose columns are exact DMD modes.
  • proj_modes: Array whose columns are projected DMD modes.
  • adjoint_modes: Array whose columns are adjoint DMD modes.
  • proj_coeffs: Array of projection coefficients for data vectors expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • adv_proj_coeffs: Array of projection coefficients for data vectors advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • R_low_order_eigvecs: Array of right eigenvectors of approximating low-order linear map.
  • L_low_order_eigvecs: Array of left eigenvectors of approximating low-order linear map.
  • correlation_array_eigvals: 1D array of eigenvalues of correlation array.
  • correlation_array_eigvecs: Array of eigenvectors of correlation array.
  • correlation_array: Correlation array; elements are inner products of data vectors with each other.
  • cross_correlation_array: Cross-correlation array; elements are inner products of data vectors with data vectors advanced in time. Going down rows, the data vector changes; going across columns the advanced data vector changes.

Attributes can be accessed using calls like res.exact_modes. To see all available attributes, use print(res).

This uses the method of snapshots, which is faster than the direct method (see compute_DMD_arrays_direct_method()) when vecs has more rows than columns, i.e., when there are more elements in a vector than there are vectors. However, it “squares” this array and its singular values, making it slightly less accurate than the direct method.

modred.dmd.compute_TLSqrDMD_arrays_direct_method(vecs, adv_vecs=None, mode_indices=None, inner_product_weights=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes Total Least Squares DMD modes using data stored in arrays, using direct method.

Args:
vecs: Array whose columns are data vectors.
Kwargs:

adv_vecs: Array whose columns are data vectors advanced in time. If not provided, then it is assumed that the vectors describe a sequential time-series. Thus vecs becomes vecs[:, :-1] and adv_vecs becomes vecs[:, 1:].

mode_indices: List of indices describing which modes to compute. Examples are range(10) or [3, 0, 6, 8]. If no mode indices are specified, then all modes will be computed.

inner_product_weights: 1D or 2D array of inner product weights. Corresponds to W in inner product v_1^* W v_2.

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

res: Results of DMD computation, stored in a namedtuple with the following attributes:

  • eigvals: 1D array of eigenvalues of approximating low-order linear map (DMD eigenvalues).
  • spectral_coeffs: 1D array of DMD spectral coefficients, calculated as the magnitudes of the projection coefficients of first (de-noised) data vector. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • exact_modes: Array whose columns are exact DMD modes.
  • proj_modes: Array whose columns are projected DMD modes.
  • adjoint_modes: Array whose columns are adjoint DMD modes.
  • proj_coeffs: Array of projection coefficients for (de-noised) data vectors expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • adv_proj_coeffs: Array of projection coefficients for (de-noised) data vectors advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • R_low_order_eigvecs: Array of right eigenvectors of approximating low-order linear map.
  • L_low_order_eigvecs: Array of left eigenvectors of approximating low-order linear map.
  • sum_correlation_array_eigvals: 1D array of eigenvalues of sum correlation array.
  • sum_correlation_array_eigvecs: Array whose columns are eigenvectors of sum correlation array.
  • proj_correlation_array_eigvals: 1D array of eigenvalues of projected correlation array.
  • proj_correlation_array_eigvecs: Array whose columns are eigenvectors of projected correlation array.

Attributes can be accessed using calls like res.exact_modes. To see all available attributes, use print(res).

This method does not square the array of vectors as in the method of snapshots (compute_DMD_arrays_snaps_method()). It’s slightly more accurate, but slower when the number of elements in a vector is more than the number of vectors, i.e., when vecs has more rows than columns.

Note that max_num_eigvals must be set to a value smaller than the rank of the dataset. In other words, if the projection basis for total-least-squares DMD is not truncated, then the algorithm reduces to standard DMD. For over-constrained datasets (number of columns in data array is larger than the number of rows), this occurs naturally. For under-constrained datasets, (number of vector objects is smaller than size of vector objects), this must be done explicitly by the user. At this time, there is no standard method for choosing a truncation level. One approach is to look at the roll-off of the correlation array eigenvalues, which contains information about the “energy” content of each projection basis vectors.

modred.dmd.compute_TLSqrDMD_arrays_snaps_method(vecs, adv_vecs=None, mode_indices=None, inner_product_weights=None, atol=1e-13, rtol=None, max_num_eigvals=None)[source]

Computes Total Least Squares DMD modes using data stored in arrays, using method of snapshots.

Args:
vecs: Array whose columns are data vectors.
Kwargs:

adv_vecs: Array whose columns are data vectors advanced in time. If not provided, then it is assumed that the vectors describe a sequential time-series. Thus vecs becomes vecs[:, :-1] and adv_vecs becomes vecs[:, 1:].

mode_indices: List of indices describing which modes to compute. Examples are range(10) or [3, 0, 6, 8]. If no mode indices are specified, then all modes will be computed.

inner_product_weights: 1D or 2D array of inner product weights. Corresponds to W in inner product v_1^* W v_2.

atol: Level below which eigenvalues of correlation array are truncated.

rtol: Maximum relative difference between largest and smallest eigenvalues of correlation array. Smaller ones are truncated.

max_num_eigvals: Maximum number of DMD eigenvalues that will be computed. This is enforced by truncating the basis onto which the approximating linear map is projected. Computationally, this corresponds to truncating the eigendecomposition of the correlation array. If set to None, no truncation will be performed, and the maximum possible number of DMD eigenvalues will be computed.

Returns:

res: Results of DMD computation, stored in a namedtuple with the following attributes:

  • eigvals: 1D array of eigenvalues of approximating low-order linear map (DMD eigenvalues).
  • spectral_coeffs: 1D array of DMD spectral coefficients, calculated as the magnitudes of the projection coefficients of first (de-noised) data vector. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • exact_modes: Array whose columns are exact DMD modes.
  • proj_modes: Array whose columns are projected DMD modes.
  • adjoint_modes: Array whose columns are adjoint DMD modes.
  • proj_coeffs: Array of projection coefficients for (de-noised) data vectors expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • adv_proj_coeffs: Array of projection coefficients for (de-noised) data vectors advanced in time, expressed as a linear combination of DMD modes. Columns correspond to vector objects, rows correspond to DMD modes. The projection is onto the span of the DMD modes using the (biorthogonal) adjoint DMD modes. Note that this is the same as a least-squares projection onto the span of the DMD modes.
  • R_low_order_eigvecs: Array of right eigenvectors of approximating low-order linear map.
  • L_low_order_eigvecs: Array of left eigenvectors of approximating low-order linear map.
  • sum_correlation_array_eigvals: 1D array of eigenvalues of sum correlation array.
  • sum_correlation_array_eigvecs: Array whose columns are eigenvectors of sum correlation array.
  • proj_correlation_array_eigvals: 1D array of eigenvalues of projected correlation array.
  • proj_correlation_array_eigvecs: Array whose columns are eigenvectors of projected correlation array.
  • correlation_array: Correlation array; elements are inner products of data vectors with each other.
  • adv_correlation_array: Advanced correlation array; elements are inner products of advanced data vectors with each other.
  • cross_correlation_array: Cross-correlation array; elements are inner products of data vectors with data vectors advanced in time. Going down rows, the data vector changes; going across columns the advanced data vector changes.

Attributes can be accessed using calls like res.exact_modes. To see all available attributes, use print(res).

This uses the method of snapshots, which is faster than the direct method (see compute_TLSqrDMD_arrays_direct_method()) when vecs has more rows than columns, i.e., when there are more elements in a vector than there are vectors. However, it “squares” this array and its singular values, making it slightly less accurate than the direct method.

Note that max_num_eigvals must be set to a value smaller than the rank of the dataset. In other words, if the projection basis for total-least-squares DMD is not truncated, then the algorithm reduces to standard DMD. For over-constrained datasets (number of columns in data array is larger than the number of rows), this occurs naturally. For under-constrained datasets, (number of vector objects is smaller than size of vector objects), this must be done explicitly by the user. At this time, there is no standard method for choosing a truncation level. One approach is to look at the roll-off of the correlation array eigenvalues, which contains information about the “energy” content of each projection basis vectors.