flowtorch.rom

flowtorch.rom.base

Definition of a common interface for all reduced-order models (ROMs).

class flowtorch.rom.base.Encoder[source]

Bases: ABC

Abstract base class for dimensionality reduction algorithms.

This base class should be used when defineing new algorithms for dimensionality reduction.

Parameters:

trained (bool) – True if the encoder was trained/set up

abstract decode(reduced_state: Tensor) Tensor[source]

Map the reduced state back to the full state.

Parameters:

reduced_state – snapshot or sequence of snapshots in reduced state space; if there is one more dimension than in reduced_state_shape, the last dimension is considered as time/batch dimension

Returns:

snapshot or sequence of snapshots in full state space

Return type:

Tensor

abstract encode(full_state: Tensor) Tensor[source]

Map the full to the reduced state.

Parameters:

data – snapshot or sequence of snapshots; if the input has one more dimension as the state (state_shape), the last dimension is considered as time/batch dimension

Returns:

snapshot or sequence of snapshots in reduced state space

Return type:

Tensor

abstract property reduced_state_size: int

Size of the reduced state vector.

abstract property state_shape: Size

Shape of the full state tensor.

abstract train(full_state: Tensor) dict[source]

Create a mapping from the full to the reduced state space.

Parameters:

full_state – time series data; the size of the last dimension equals the number of snapshots (batch dimension)

Returns:

information about the training process

Return type:

dict

property trained: bool

Get the training state

Returns:

True if training was completed

Return type:

bool

class flowtorch.rom.base.ROM(reduced_state: Tensor, encoder: Encoder)[source]

Bases: ABC

Abstract base class for reduced-order models.

This base class should be used when defining new ROMs.

property encoder: Encoder

Return encoder instance.

predict(initial_state: Tensor, end_time: float, step_size: float) Tensor[source]

Predict the evolution of a given initial full state vector.

Parameters:
  • initial_state (Tensor) – state from which to start

  • end_time (float) – when to stop the simulation; the corresponding start time is always assumed to be zero

  • step_size (float) – time step size

Returns:

evolution of the full state vector; the last dimension corresponds to the time/batch dimension

Return type:

Tensor

abstract predict_reduced(initial_state: Tensor, end_time: float, step_size: float) Tensor[source]

Predict the evolution of a given initial reduced state vector.

Parameters:
  • initial_state (Tensor) – initial reduced state vector

  • end_time (float) – when to stop the simulation; the corresponding start time is always assumed to be zero

  • step_size (float) – time step size

Returns:

evolution of the reduced state vector; the last dimension corresponds to the time/batch dimension

Return type:

Tensor

flowtorch.rom.cnm

Implementation of cluster-based network modeling (CNM).

A reference implementation by Daniel Fernex is available on Github. Theoretical concepts are covered in the accompanying publication.

class flowtorch.rom.cnm.CNM(reduced_state: Tensor, encoder: Encoder, dt: float, n_clusters: int = 10, model_order: int = 1, cluster_config: dict = {})[source]

Bases: ROM

Cluster-based network modeling implementation.

In contrast to the original implementation, the clustering, the computation of transition probabilities, and the propagation in time are condensed in a single class. The class also relies on Numpy and Scikit-Learn for K-means, KD-tree, and interpolation. These components might be replaced with PyTorch implementations in future releases.

The cluster attribute has a member called labels. Each label corresponds to the cluster id associated with each sample of the input data. Since the input data are sorted according to the time at which they were sampled, we can interpret the labels as a sequence of visited clusters over time However, if two or more consecutive snapshots belong to the same cluster, we are only interested in the next state in the sequence that is different from the current cluster because we want to model the temporal behavior as a transition between cluster centroids. Therefore, sequential duplicates must be removed.

Parameters:
  • dt (float) – time step between two snapshots

  • n_clusters (int, optional) – number of clusters; must be larger than one

  • model_order (int, optional) – number of past cluster to consider when predicting the next cluster

  • Q (Dict[str, np.ndarray]) – transition probabilities

  • T (Dict[str, float]) – transition times

  • times (List[float]) – times at which cluster centroids were visited

  • visited_clusters (List[int]) – clusters visited during the temporal evolution

  • cluster_centers (np.ndarray) – centroids of clusters

Examples

>>> from flowtorch.rom import SVDEncoder, CNM
... assemble data matrix
>>> encoder = SVDEncoder(rank=20)
>>> info = encoder.train(data_matrix)
>>> reduced_state = encoder.encode(data_matrix)
>>> cnm = CNM(reduced_state, encoder, n_clusters=20, model_order=4)
>>> prediction = cnm.predict(data_matrix[:, :5], end_time=1.0, step_size=0.1)
>>> reduced_prediction = cnm.predict_reduced(reduced_state[:, :5], 1.0, 0.1)
>>> full_prediction = encoder.decode(reduced_prediction)
property Q
property T
property cluster_centers
property dt: float
property model_order: int
property n_clusters: int
predict_reduced(initial_state: Tensor, end_time: float, step_size: float) Tensor[source]

Advance given reduced state in time.

Parameters:
  • initial_state (pt.Tensor) – initial reduced state vector or sequence of reduced state vectors; the state vectors must form the columns of the input tensor is a sequence is given

  • end_time (float) – time at which to step the simulation; the simulation time always starts at zero

  • step_size (float) – time step size to advance simulation time; note that the sampling of the next cluster is independent of the step size; the step size is only used to interpolate states inbetween the sampled clusters

Returns:

temporal evolution of the given reduced state vector; each column forms a temporal snapshot

Return type:

pt.Tensor

property times: list
property visited_clusters: list

flowtorch.rom.svd_encoder

Encoder based on the singular value decomposition (SVD).

class flowtorch.rom.svd_encoder.SVDEncoder(rank: int | None = None)[source]

Bases: Encoder

SVD-based dimensionality reduction for ROM.

Examples

>>> from flowtorch import DATASETS
>>> from flowtorch.data import FOAMDataloader
>>> from flowtorch.rom import SVDEncoder
>>> loader = FOAMDataloader(DATASETS["of_cylinder2D_binary"])
>>> data = loader.load_snapshot("p", loader.write_times[1:11])
>>> data.shape
torch.Size([13678, 10]
>>> encoder = SVDEncoder(rank=10)
>>> info = encoder.train(data)
>>> reduced_state = encoder.encode(data)
>>> reduced_state.shape
torch.Size([10, 10]
>>> full_state = encoder.decode(reduced_state)
>>> full_state.shape
torch.Size([13678, 10]
decode(reduced_state: Tensor) Tensor[source]

Compute the full projection onto the POD modes.

Parameters:

reduced_state (pt.Tensor) – 1D or 2D tensor, in which each column holds the mode coefficients of a given state; if the input has two dimensions, the second dimension is considered as the batch dimension

Returns:

full state vector in the subspace spanned by the POD modes

Return type:

pt.Tensor

encode(full_state: Tensor) Tensor[source]

Project one or multiple state vectors onto the POD modes.

This function computes the scalar projections of one or more state vectors onto each POD mode. The result is vector (single state) or a matrix (sequence of states) of scalar projections, in which the row index corresponds to the associated POD mode and the column index corresponds to the associated snapshot in the sequence.

Parameters:

full_state (pt.Tensor) – [description]

Raises:

ValueError – [description]

Returns:

[description]

Return type:

pt.Tensor

property reduced_state_size: int

Get the size of the reduced state.

Returns:

size of the reduced state.

Return type:

int

property state_shape: Size

Get the size of the full state.

Returns:

size of the full state

Return type:

pt.Size

train(data: Tensor) dict[source]

Compute the POD modes of a given data matrix.

Parameters:

data (pt.Tensor) – data matrix containing a sequence of snapshots, where each snapshot corresponds to a column vector of the data matrix

Returns:

empty dictionary since there is no real training process

Return type:

dict

flowtorch.rom.utils

Collection of utilities for reduced-order modeling.

flowtorch.rom.utils.check_int_larger_than(value: int, limit: int, name: str)[source]

Check if input is an integer larger than a given lower limit.

Parameters:
  • value (int) – input value to check

  • limit (int) – the value must be larger than the limit

  • name (str) – name of the parameter

Raises:

ValueError – if the argument is not an integer

flowtorch.rom.utils.check_larger_than(value: int | float, limit: int | float, name: str)[source]

Check if a scalar value is larger than a given lower limit.

Parameters:
  • value (Union[int, float]) – scalar value to check

  • value – lower limit to check against

  • name (str) – name of the parameter

Raises:

ValueError – if the argument is less than or equal to the lower limit

flowtorch.rom.utils.log_time(func) dict[source]

Measure and log a function’s execution time.

Parameters:

func (Callable) – function to be executed; the function is expected to return a dictionary

Returns:

dictionary returned by the wrapped function with additional entry for execution time

Return type:

dict

flowtorch.rom.utils.remove_sequential_duplicates(sequence: ndarray) ndarray[source]

Get sequence of integers without sequential duplicates.

Parameters:

sequence (np.ndarray) – input sequence to check

Returns:

sequence without sequential duplicates

Return type:

np.ndarray