Skip to content
Shane edited this page Feb 11, 2021 · 6 revisions

This page contains code documentation for rom_operator_inference classes and functions. The code itself is also internally documented and can be accessed on the fly with dynamic object introspection, e.g.,

>>> import rom_operator_inference as roi
>>> help(roi.InferredContinuousROM)

Contents

ROM Classes

The core of rom_operator_inference is highly object oriented and defines several ROM classes that serve as the workhorse of the package. The API for these classes adopts some principles from the scikit-learn API: there are fit() and predict() methods, hyperparameters are set in the constructor, estimated attributes end with underscore, and so on.

Each class corresponds to a type of full-order model (continuous vs. discrete, non-parametric vs. parametric) and a strategy for constructing the ROM.

Class Name Problem Statement ROM Strategy
InferredContinuousROM Operator Inference
InferredDiscreteROM Operator Inference
InterpolatedInferredContinuousROM Operator Inference
InterpolatedInferredDiscreteROM Operator Inference
AffineInferredContinuousROM Operator Inference
AffineInferredDiscreteROM Operator Inference
IntrusiveContinuousROM Intrusive Projection
IntrusiveDiscreteROM Intrusive Projection
AffineIntrusiveContinuousROM Intrusive Projection
AffineIntrusiveDiscreteROM Intrusive Projection

The following function may be helpful for selecting an appropriate class.

select_model_class(time, rom_strategy, parametric=False): select the appropriate ROM model class for the situation. Parameters:

  • time: Type of full-order model to be reduced, either "continuous" or "discrete".
  • rom_strategy: Whether to use Operator Inference ("inferred") or intrusive projection ("intrusive") to compute the operators of the reduced model.
  • parametric: Whether or not the model depends on an external parameter, and how to handle the parametric dependence. Options:
    • False (default): the problem is nonparametric.
    • "interpolated": construct individual models for each sample parameter and interpolate them for general parameter inputs. Only valid for rom_strategy="inferred", and only when the parameter is a scalar.
    • "affine": one or more operators in the problem depends affinely on the parameter. Only valid for rom_strategy="intrusive".

The return value is the class type for the situation, e.g., InferredContinuousROM.

Constructor

All ROM classes are instantiated with a single argument, modelform, which is a string denoting the structure of the full-order operator f. Each character in the string corresponds to a single term of the operator, given in the following table.

Character Name Continuous Term Discrete Term
c Constant
A Linear
H Quadratic
G Cubic
B Input

These are all input as a single string. Examples:

modelform Continuous ROM Structure Discrete ROM Structure
"A"
"cA"
"HB"
"cAHB"

Attributes

All ROM classes have the following attributes.

  • Structure of model:

    • modelform: set in the constructor.
    • has_constant: boolean, whether or not there is a constant term c.
    • has_linear: boolean, whether or not there is a linear term Ax.
    • has_quadratic: boolean, whether or not there is a quadratic term H(xx).
    • has_cubic: boolean, whether or not there is a cubic term G(xxx).
    • has_inputs: boolean, whether or not there is an input term Bu.
  • Dimensions, set in fit():

    • n: Dimension of the original model
    • r: Dimension of the learned reduced-order model
    • m: Dimension of the input u, or None if 'B' is not in modelform.
  • Reduced operators c_, A_, H_, G_, and B_, learned in fit(): the NumPy arrays corresponding to the learned parts of the reduced-order model. Set to None if the operator is not included in the prescribed modelform (e.g., if modelform="AHG", then c_ and B_ are None).

  • Basis matrix Vr: the n x r basis defining the mapping between the n-dimensional space of the full-order model and the reduced r-dimensional subspace of the reduced-order model (e.g., POD basis). This is the first input to all fit() methods. To save memory, inferred (but not intrusive) ROM classes allow entering Vr=None, which then assumes that other inputs for training are already projected to the r-dimensional subspace (e.g., VrTX instead of X).

  • Reduced model function f_, learned in fit(): the ROM function, defined by the reduced operators listed above. For continuous models, f_ has the following signature:

def f_(t, x_, u):
    """ROM function for continuous models.

    Parameters
    ----------
    t : float
        Time, a scalar.

    x_ : (r,) ndarray
        Reduced state vector.

    u : func(float) -> (m,)
        Input function that maps time `t` to an input vector of length m.
    """

For discrete models, the signature is the following.

def f_(x_, u):
    """ROM function for discrete models.

    Parameters
    ----------
    x_ : (r,) ndarray
        Reduced state vector.

    u : (m,) ndarray
        Input vector of length m corresponding to the state.
    """

The input argument u is only used if B is in modelform.

Model Persistence

Trained ROM objects can be saved in HDF5 format with the save_model() method, and recovered later with the load_model() function. Such files store metadata for the model class and structure, the reduced-order model operators (c_, A_, etc.), other attributes learned in fit(), and (optionally) the basis Vr.

load_model(loadfile): Load a serialized model from an HDF5 file, created previously from a ROM object's save_model() method.

<ROMclass>.save_model(savefile, save_basis=True, overwrite=False): Serialize the learned model, saving it in HDF5 format. The model can then be loaded with load_model(). Currently implemented for nonparametric classes only. Parameters:

  • savefile: File to save to. If it does not end with '.h5', this extension will be tacked on to the end.
  • savebasis: If True, save the basis Vr as well as the reduced operators. If False, only save reduced operators.
  • overwrite: If True and the specified file already exists, overwrite the file. If False and the specified file already exists, raise an error.
>>> import rom_operator_inference as roi

# Assume model is a trained roi.InferredContinuousROM object.
>>> model.save_model("trained_rom.h5")          # Save a trained model.
>>> model2 = roi.load_model("trained_rom.h5")   # Load a model from file.

InferredContinuousROM

This class constructs a reduced-order model for the continuous, nonparametric system

via Operator Inference [1]. That is, given snapshot data, a basis, and a form for a reduced model, it computes the reduced model operators by solving an ordinary least-squares problem (see Operator-Inference).

InferredContinuousROM.fit(Vr, X, Xdot, U=None, P=0): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order model will be projected (for example, a POD basis matrix; see pre.pod_basis()). Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrix X. If given as None, X is assumed to be the projected snapshot matrix VrTX and Xdot is assumed to be the projected time derivative matrix.
    • X: An n x k snapshot matrix of solutions to the full-order model, or the r x k projected snapshot matrix VrTX. Each column is one snapshot.
    • Xdot: n x k snapshot time derivative matrix for the full-order model, or the r x k projected snapshot time derivative matrix. Each column is the time derivative dx/dt for the corresponding column of X. See the pre submodule for some simple derivative approximation tools.
    • U: m x k input matrix (or a k-vector if m = 1). Each column is the input vector for the corresponding column of X. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns
    • Trained InferredContinuousROM object.

InferredContinuousROM.predict(x0, t, u=None, **options): Simulate the learned reduced-order model with scipy.integrate.solve_ivp().

  • Parameters
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector). If Vr=None in fit(), this must be the projected initial state VrTx0.
    • t: Time domain, an nt-vector, over which to integrate the reduced-order model.
    • u: Input as a function of time, that is, a function mapping a float to an m-vector (or to a scalar if m = 1). Alternatively, the m x nt matrix (or nt-vector if m = 1) where column j is the input vector corresponding to time t[j]. In this case, u(t) is approximated by a cubic spline interpolating the given inputs. This argument is only required if 'B' is in modelform.
    • Other keyword arguments for scipy.integrate.solve_ivp().
  • Returns
    • X_ROM: n x nt matrix of approximate solution to the full-order system over t, or, if Vr=None in fit(), the r x nt solution in the reduced-order space. Each column is one snapshot of the solution.

InferredDiscreteROM

This class constructs a reduced-order model for the discrete, nonparametric system

via Operator Inference.

InferredDiscreteROM.fit(Vr, X, U=None, P=0): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order model will be projected. Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrix X. If given as None, X is assumed to be the projected snapshot matrix VrTX.
    • X: n x k snapshot matrix of solutions to the full-order model, or the r x k projected snapshot matrix VrTX. Each column is one snapshot.
    • U: m x k-1 input matrix (or a (k-1)-vector if m = 1). Each column is the input for the corresponding column of X. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns
    • Trained InferredDiscreteROM object.

InferredDiscreteROM.predict(x0, niters, U=None): Step forward the learned ROM niters steps.

  • Parameters
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector). If Vr=None in fit(), this must be the projected initial state VrTx0.
    • niters: Number of times to step the system forward.
    • U: Inputs for the next niters-1 time steps, as an m x niters-1 matrix (or an (niters-1)-vector if m = 1). This argument is only required if 'B' is in modelform.
  • Returns
    • X_ROM: n x niters matrix of approximate solutions to the full-order system, including the initial condition; or, if Vr=None in fit(), the r x niters solution in the reduced-order space. Each column is one iteration of the solution.

InterpolatedInferredContinuousROM

This class constructs a reduced-order model for the continuous, parametric system

via Operator Inference. The strategy is to take snapshot data for several parameter samples and a global basis, compute a reduced model for each parameter sample via Operator Inference, then construct a general parametric model by interpolating the entries of the inferred operators [1].

InterpolatedInferredContinuousROM.fit(Vr, µs, Xs, Xdots, Us=None, P=0): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters
    • Vr: The (global) n x r basis for the linear reduced space on which the full-order model will be projected. Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrices Xs. If given as None, Xs is assumed to be the list of the projected snapshot matrices VrTXi and Xdots is assumed to be the list of projected time derivative matrices.
    • µs: s parameter values corresponding to the snapshot sets.
    • Xs: List of s snapshot matrices, each n x k (full-order solutions) or r x k (projected solutions). The _i_th array Xs[i] corresponds to the _i_th parameter, µs[i]; each column each of array is one snapshot.
    • Xdots: List of s snapshot time derivative matrices, each n x k (full-order velocities) or r x k (projected velocities). The _i_th array Xdots[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Xdots[i][:,j], is the time derivative dx/dt for the corresponding snapshot column Xs[i][:,j].
    • Us: List of s input matrices, each m x k (or a k-vector if m=1). The _i_th array Us[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Us[i][:,j], is the input for the corresponding snapshot Xs[i][:,j]. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns
    • Trained InterpolatedInferredContinuousROM object.

InterpolatedInferredContinuousROM.predict(µ, x0, t, u=None, **options): Simulate the learned reduced-order model with scipy.integrate.solve_ivp().

  • Parameters
    • µ: Parameter value at which to simulate the ROM.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector). If Vr=None in fit(), this must be the projected initial state VrTx0.
    • t: Time domain, an nt-vector, over which to integrate the reduced-order model.
    • u: Input as a function of time, that is, a function mapping a float to an m-vector (or to a scalar if m = 1). Alternatively, the m x nt matrix (or nt-vector if m = 1) where column j is the input vector corresponding to time t[j]. In this case, u(t) is approximated by a cubic spline interpolating the given inputs. This argument is only required if 'B' is in modelform.
    • Other keyword arguments for scipy.integrate.solve_ivp().
  • Returns
    • X_ROM: n x nt matrix of approximate solution to the full-order system over t, or, if Vr=None in fit(), the r x nt solution in the reduced-order space. Each column is one snapshot of the solution.

InterpolatedInferredDiscreteROM

This class constructs a reduced-order model for the continuous, parametric system

via Operator Inference. The strategy is to take snapshot data for several parameter samples and a global basis, compute a reduced model for each parameter sample via Operator Inference, then construct a general parametric model by interpolating the entries of the inferred operators [1].

InterpolatedInferredDiscreteROM.fit(Vr, µs, Xs, Us=None, P=0): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters
    • Vr: The (global) n x r basis for the linear reduced space on which the full-order model will be projected. Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrices Xs. If given as None, Xs is assumed to be the list of the projected snapshot matrices VrTXi.
    • µs: s parameter values corresponding to the snapshot sets.
    • Xs: List of s snapshot matrices, each n x k (full-order solutions) or r x k (projected solutions). The _i_th array Xs[i] corresponds to the _i_th parameter, µs[i]; each column each of array is one snapshot.
    • Us: List of s input matrices, each m x k (or a k-vector if m=1). The _i_th array Us[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Us[i][:,j], is the input for the corresponding snapshot Xs[i][:,j]. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns
    • Trained InterpolatedInferredDiscreteROM object.

InterpolatedInferredDiscreteROM.predict(µ, x0, niters, U=None): Step forward the learned ROM niters steps.

  • Parameters
    • µ: Parameter value at which to simulate the ROM.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector). If Vr=None in fit(), this must be the projected initial state VrTx0.
    • niters: Number of times to step the system forward.
    • U: Inputs for the next niters-1 time steps, as an m x niters-1 matrix (or an (niters-1)-vector if m = 1). This argument is only required if 'B' is in modelform.
  • Returns
    • X_ROM: n x niters matrix of approximate solutions to the full-order system, including the initial condition; or, if Vr=None in fit(), the r x nt solution in the reduced-order space. Each column is one iteration of the solution.

AffineInferredContinuousROM

This class constructs a reduced-order model via Operator Inference for the continuous, affinely parametric system

where the operators that define f may only depend affinely on the parameter, e.g.,

AffineInferredContinuousROM.fit(Vr, µs, affines, Xs, Xdots, Us=None, P=0): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters
    • Vr: The (global) n x r basis for the linear reduced space on which the full-order model will be projected. Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrices Xs. If given as None, Xs is assumed to be the list of the projected snapshot matrices VrTXi and Xdots is assumed to be the list of projected time derivative matrices.
    • µs: s parameter values corresponding to the snapshot sets.
    • affines A dictionary mapping labels of the operators that depend affinely on the parameter to the list of functions that define that affine dependence. The keys are entries of modelform. For example, if the constant term has the affine structure c(µ) = θ1(µ)c1 + θ2(µ)c2 + θ3(µ)c3, then 'c' -> [θ1, θ2, θ3].
    • Xs: List of s snapshot matrices, each n x k (full-order solutions) or r x k (projected solutions). The _i_th array Xs[i] corresponds to the _i_th parameter, µs[i]; each column each of array is one snapshot.
    • Xdots: List of s snapshot time derivative matrices, each n x k (full-order velocities) or r x k (projected velocities). The _i_th array Xdots[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Xdots[i][:,j], is the time derivative dx/dt for the corresponding snapshot column Xs[i][:,j].
    • Us: List of s input matrices, each m x k (or a k-vector if m=1). The _i_th array Us[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Us[i][:,j], is the input for the corresponding snapshot Xs[i][:,j]. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns:
    • Trained AffineInferredContinuousROM object.

AffineInferredContinuousROM.predict(µ, x0, t, u=None, **options): Simulate the learned reduced-order model at the given parameter value with scipy.integrate.solve_ivp().

  • Parameters
    • µ: Parameter value at which to simulate the model.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • t: Time domain, an nt-vector, over which to integrate the reduced-order model.
    • u: Input as a function of time, that is, a function mapping a float to an m-vector (or to a scalar if m = 1). Alternatively, the m x nt matrix (or nt-vector if m = 1) where column j is the input vector corresponding to time t[j]. In this case, u(t) is approximated by a cubic spline interpolating the given inputs. This argument is only required if 'B' is in modelform.
    • Other keyword arguments for scipy.integrate.solve_ivp().
  • Returns
    • X_ROM: n x nt matrix of approximate solution to the full-order system over t. Each column is one snapshot of the solution.

AffineInferredDiscreteROM

This class constructs a reduced-order model via Operator Inference for the discrete, affinely parametric system

where the operators that define f may only depend affinely on the parameter, e.g.,

AffineInferredDiscreteROM.fit(Vr, µs, affines, Xs, Us, P): Compute the operators of the reduced-order model via Operator Inference.

  • Parameters:
    • Vr: The (global) n x r basis for the linear reduced space on which the full-order model will be projected. Each column is a basis vector. The column space of Vr should be a good approximation of the column space of the full-order snapshot matrices Xs. If given as None, Xs is assumed to be the list of the projected snapshot matrices VrTXi.
    • µs: s parameter values corresponding to the snapshot sets.
    • affines A dictionary mapping labels of the operators that depend affinely on the parameter to the list of functions that define that affine dependence. The keys are entries of modelform. For example, if the constant term has the affine structure c(µ) = θ1(µ)c1 + θ2(µ)c2 + θ3(µ)c3, then 'c' -> [θ1, θ2, θ3].
    • Xs: List of s snapshot matrices, each n x k (full-order solutions) or r x k (projected solutions). The _i_th array Xs[i] corresponds to the _i_th parameter, µs[i]; each column each of array is one snapshot.
    • Us: List of s input matrices, each m x k (or a k-vector if m=1). The _i_th array Us[i] corresponds to the _i_th parameter, µs[i]. The _j_th column of the _i_th array, Us[i][:,j], is the input for the corresponding snapshot Xs[i][:,j]. Only required when 'B' is in modelform.
    • P: Tikhonov regularization matrix for the least-squares problem; see lstsq.
  • Returns:
    • Trained AffineInferredDiscreteROM object.

AffineInferredDiscreteROM.predict(µ, x0, niters, U=None): Step forward the learned ROM niters steps at the given parameter value.

  • Parameters
    • µ: Parameter value at which to simulate the model.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • niters: Number of times to step the system forward.
    • U: Inputs for the next niters-1 time steps, as an m x niters-1 matrix (or an (niters-1)-vector if m = 1). This argument is only required if 'B' is in modelform.
  • Returns
    • X_ROM: n x niters matrix of approximate solutions to the full-order system, including the initial condition. Each column is one iteration of the solution.

IntrusiveContinuousROM

This class constructs a reduced-order model for the continuous, nonparametric system

via intrusive projection, i.e.,

where ⊗ denotes the full matrix Kronecker product. The class requires the actual full-order operators (c, A, H, G, and/or B) that define f.

IntrusiveContinuousROM.fit(Vr, operators): Compute the operators of the reduced-order model by projecting the operators of the full-order model.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order operators will be projected.
    • operators: A dictionary mapping labels to the full-order operators that define f. The operators are indexed by the entries of modelform; for example, if modelform="cHB", then operators={'c':c, 'H':H, 'B':B}.
  • Returns
    • Trained IntrusiveContinuousROM object.

IntrusiveContinuousROM.predict(x0, t, u=None, **options): Simulate the learned reduced-order model with scipy.integrate.solve_ivp().

  • Parameters
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • t: Time domain, an nt-vector, over which to integrate the reduced-order model.
    • u: Input as a function of time, that is, a function mapping a float to an m-vector (or to a scalar if m = 1). Alternatively, the m x nt matrix (or nt-vector if m = 1) where column j is the input vector corresponding to time t[j]. In this case, u(t) is approximated by a cubic spline interpolating the given inputs. This argument is only required if 'B' is in modelform.
    • Other keyword arguments for scipy.integrate.solve_ivp().
  • Returns
    • X_ROM: n x nt matrix of approximate solution to the full-order system over t. Each column is one snapshot of the solution.

IntrusiveDiscreteROM

This class constructs a reduced-order model for the discrete, nonparametric system

via intrusive projection, i.e.,

The class requires the actual full-order operators (c, A, H, G, and/or B) that define f.

IntrusiveDiscreteROM.fit(Vr, operators): Compute the operators of the reduced-order model by projecting the operators of the full-order model.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order operators will be projected.
    • operators: A dictionary mapping labels to the full-order operators that define f. The operators are indexed by the entries of modelform; for example, if modelform="cHB", then operators={'c':c, 'H':H, 'B':B}.
  • Returns
    • Trained IntrusiveDiscreteROM object.

IntrusiveDiscreteROM.predict(x0, niters, U=None): Step forward the learned ROM niters steps.

  • Parameters
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • niters: Number of times to step the system forward.
    • U: Inputs for the next niters-1 time steps, as an m x niters-1 matrix (or an (niters-1)-vector if m = 1). This argument is only required if 'B' is in modelform.
  • Returns
    • X_ROM: n x niters matrix of approximate solutions to the full-order system, including the initial condition. Each column is one iteration of the solution.

AffineIntrusiveContinuousROM

This class constructs a reduced-order model for the continuous, affinely parametric system

where the operators that define f may only depend affinely on the parameter, e.g.,

The reduction is done via intrusive projection, i.e.,

The class requires the actual full-order operators (c, A, H, and/or B) that define f and the functions that define any affine parameter dependencies (i.e., the θ functions).

AffineIntrusiveContinuousROM.fit(Vr, affines, operators): Compute the operators of the reduced-order model by projecting the operators of the full-order model.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order operators will be projected.
    • affines A dictionary mapping labels of the operators that depend affinely on the parameter to the list of functions that define that affine dependence. The keys are entries of modelform. For example, if the constant term has the affine structure c(µ) = θ1(µ)c1 + θ2(µ)c2 + θ3(µ)c3, then 'c' -> [θ1, θ2, θ3].
    • operators: A dictionary mapping labels to the full-order operators that define f. The keys are entries of modelform. Terms with affine structure should be given as a list of the component matrices. For example, suppose modelform="cA". If A has the affine structure A(µ) = θ1(µ)A1 + θ2(µ)A2, then 'A' -> [A1, A2]. If c does not vary with the parameter, then 'c' -> c, the complete full-order order.
  • Returns:
    • Trained AffineIntrusiveContinuousROM object.

AffineIntrusiveContinuousROM.predict(µ, x0, t, u=None, **options): Simulate the learned reduced-order model at the given parameter value with scipy.integrate.solve_ivp().

  • Parameters
    • µ: Parameter value at which to simulate the model.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • t: Time domain, an nt-vector, over which to integrate the reduced-order model.
    • u: Input as a function of time, that is, a function mapping a float to an m-vector (or to a scalar if m = 1). Alternatively, the m x nt matrix (or nt-vector if m = 1) where column j is the input vector corresponding to time t[j]. In this case, u(t) is approximated by a cubic spline interpolating the given inputs. This argument is only required if 'B' is in modelform.
    • Other keyword arguments for scipy.integrate.solve_ivp().
  • Returns
    • X_ROM: n x nt matrix of approximate solution to the full-order system over t. Each column is one snapshot of the solution.

AffineIntrusiveDiscreteROM

This class constructs a reduced-order model for the discrete, affinely parametric system

where the operators that define f may only depend affinely on the parameter, e.g.,

The reduction is done via intrusive projection, i.e.,

The class requires the actual full-order operators (c, A, H, and/or B) that define f and the functions that define any affine parameter dependencies (i.e., the θ functions).

AffineIntrusiveDiscreteROM.fit(Vr, affines, operators): Compute the operators of the reduced-order model by projecting the operators of the full-order model.

  • Parameters
    • Vr: n x r basis for the linear reduced space on which the full-order operators will be projected.
    • affines A dictionary mapping labels of the operators that depend affinely on the parameter to the list of functions that define that affine dependence. The keys are entries of modelform. For example, if the constant term has the affine structure c(µ) = θ1(µ)c1 + θ2(µ)c2 + θ3(µ)c3, then 'c' -> [θ1, θ2, θ3].
    • operators: A dictionary mapping labels to the full-order operators that define f. The keys are entries of modelform. Terms with affine structure should be given as a list of the component matrices. For example, suppose modelform="cA". If A has the affine structure A(µ) = θ1(µ)A1 + θ2(µ)A2, then 'A' -> [A1, A2]. If c does not vary with the parameter, then 'c' -> c, the complete full-order order.
  • Returns:
    • Trained AffineIntrusiveDiscreteROM object.

AffineIntrusiveDiscreteROM.predict(µ, x0, niters, U=None): Step forward the learned ROM niters steps at the given parameter value.

  • Parameters
    • µ: Parameter value at which to simulate the model.
    • x0: Initial state vector, either full order (n-vector) or projected to reduced order (r-vector).
    • niters: Number of times to step the system forward.
    • U: Inputs for the next niters-1 time steps, as an m x niters-1 matrix (or an (niters-1)-vector if m = 1). This argument is only required if 'B' is in modelform.
  • Returns
    • X_ROM: n x niters matrix of approximate solutions to the full-order system, including the initial condition. Each column is one iteration of the solution.

Least-squares Solvers

The lstsq module defines classes to handle the solution of the actual Operator Inference least-squares problems. Each class has fit(A, B)and predict(𝚪) methods, where the arguments A, B = [b1, ..., br], and 𝚪 correspond to the problem

which is the Frobenius-norm least-squares problem decoupled by columns. In the context of Operator Inference, A is the data matrix D, B is the projected time derivative matrix RT, and 𝚪 determines the penalties on the entries of the operator matrix OT. See Operator-Inference for additional mathematical explanation.

lstsq.SolverL2

Solve (LS) with L2 regularization, meaning 𝚪j = λI, j = 1,...,r, for some λ ≥ 0. Setting λ = 0 is equivalent to ordinary, non-regularized least squares. The solution is obtained via the SVD of A:

lstsq.SolverL2.fit(A, B): Take the SVD of A in preparation to solve (LS).

  • Parameters
    • A: k x d data matrix.
    • B: k x r right-hand-side matrix.
  • Returns
    • Initialized lstsq.SolverL2 object.

lstsq.SolverL2.predict(λ): Solve (LS) with regularization hyperparameter λ.

  • Parameters
    • λ: Scalar, non-negative regularization hyperparameter.
  • Returns
    • X: d x r matrix of least-squares solutions.

lstsq.SolverL2Decoupled

Solve (LS) with a different L2 regularization for each column of B, i.e., 𝚪j = λjI where λj ≥ 0, j = 1,...,r. The solution is obtained via the SVD of A:

lstsq.SolverL2Decoupled.fit(): Take the SVD of A in preparation to solve (LS).

  • Parameters
    • A: k x d data matrix.
    • B: k x r right-hand-side matrix.
  • Returns
    • Initialized lstsq.SolverL2Decoupled object.

lstsq.SolverL2Decoupled.predict(λs): Solve (LS) with parameters λs = [λ1, ..., λr].

  • Parameters
    • λs: r non-negative regularization hyperparameters.
  • Returns
    • X: d x r matrix of least-squares solutions.

lstsq.SolverTikhonov

Solve (LS) with a given matrix 𝚪, using the same 𝚪 for each j: 𝚪j = 𝚪 for j = 1, ..., r. The solution is obtained by solving the modified normal equations (ATA + 𝚪T𝚪)xj = ATbj via scipy.linalg.solve() with assume_a="pos" (POSV from LAPACK). If these equations are highly ill-conditioned, the solver uses scipy.linalg.lstsq() (GELSD from LAPACK) as a backup, which is slightly slower but generally more stable.

lstsq.SolverTikhonov.fit(A, B): Compute ATA and ATB in preparation to solve (LS).

  • Parameters
    • A: k x d data matrix.
    • B: k x r right-hand-side matrix.
  • Returns
    • Initialized lstsq.SolverTikhonov object.

lstsq.SolverTikhonov.predict(P):

  • Parameters
    • P: d x d regularization matrix 𝚪 OR a d-vector, in which case 𝚪 is interpreted as diag(P).
  • Returns
    • X: d x r matrix of least-squares solutions.

lstsq.SolverTikhonovDecoupled

Solve (LS) with given matrices 𝚪j, j = 1, ..., r, by solving the modified normal equations (ATA + 𝚪jT𝚪j)xj = ATbj (unless ill conditioned).

lstsq.SolverTikhonovDecoupled.fit(A, B):

  • Parameters
    • A: k x d data matrix.
    • B: k x r right-hand-side matrix.
  • Returns
    • Initialized lstsq.SolverTikhonovDecoupled object.

lstsq.SolverTikhonovDecoupled.predict(Ps):

  • Parameters
    • Ps: sequence of r (d x d) regularization matrices 𝚪1, ..., 𝚪r, OR sequence of r d-vectors, in which case 𝚪j is interpreted as diag(P[j]).
  • Returns
    • X: d x r matrix of least-squares solutions.

The following helper functions interface with the least-squares solver classes.

lstsq.lstsq_size(modelform, r, m=0, affines=None): Calculate the number of columns of the operator matrix O in the Operator Inference least squares problem (called d(r,m) in Operator-Inference). Useful for determining the dimensions of the regularization matrix 𝚪.

  • Parameters
    • modelform: the structure of the desired model.
    • r: Dimension of the reduced order model.
    • m: Dimension of the inputs of the model. Must be zero unless 'B' is in modelform.
    • affines: A dictionary mapping labels of the operators that depend affinely on the parameter to the list of functions that define that affine dependence. The keys are entries of modelform. For example, if the constant term has the affine structure c(µ) = θ1(µ)c1 + θ2(µ)c2 + θ3(µ)c3, then 'c' -> [θ1, θ2, θ3].
  • Returns
    • Number of columns of the unknown matrix in the Operator Inference least squares problem.

lstsq.solver(A, b, P=0): Select and initialize an appropriate solver for the (LS), i.e., pick a solver class and call its fit() method.

lstsq.solve(A, b, P=0): Select, initialize, and execute an appropriate solver for (LS), i.e., pick a solver class and call its fit() and predict() methods.

Preprocessing Tools

The pre submodule is a collection of common routines for preparing data to be used by the ROM classes. None of these routines are novel, but they may be instructive for new Python users.

pre.shift(X, shift_by=None): Shift the columns of X by the vector shift_by. If shift_by=None, shift X by the mean of its columns.

pre.scale(X, scale_to, scale_from=None): Scale the entries of X from the interval [scale_from[0], scale_from[1]] to the interval [scale_to[0], scale_to[1]]. If scale_from=None, learn the scaling by setting scale_from[0] = min(X); scale_from[1] = max(X).

pre.pod_basis(X, r=None, mode="dense", **options): Compute the POD basis of rank r and the associated singular values for a snapshot matrix X. If r = None, compute all singular vectors / values. This function simply wraps a few SVD methods, selected by mode:

Use **options to specify additional parameters for these wrapped functions.

pre.svdval_decay(singular_values, eps, plot=False): Count the number of singular values that are greater than eps. The singular values can be computed with, for example, singular_values = scipy.linalg.svdvals(X) where X is a snapshot matrix. If plot=True, plot the singular values on a log scale.

pre.cumulative_energy(singular_values, thresh, plot=False): Compute the number of singular values needed to surpass the energy threshold thresh; the energy of the first r singular values is defined by

The singular values can be computed with, for example, singular_values = scipy.linalg.svdvals(X) where X is a snapshot matrix. If plot=True, plot the cumulative energy on a log scale.

pre.projection_error(X, Vr): Compute the relative Frobenius-norm projection error on X induced by the basis matrix Vr,

pre.minimal_projection_error(X, V, eps, plot=False): Compute the number of POD basis vectors required to obtain a projection error less than eps, up to the number of columns of V. If plot=True, plot the projection error on a log scale as a function of the basis size.

pre.reproject_continuous(f, Vr, X, U=None): Sample re-projected trajectories [5] of the continuous system of ODEs defined by f.

pre.reproject_discrete(f, Vr, x0, niters, U=None): Sample re-projected trajectories [5] of the discrete dynamical system defined by f.

pre.xdot_uniform(X, dt, order=2): Approximate the first time derivative of a snapshot matrix X in which the snapshots are evenly spaced in time.

pre.xdot_nonuniform(X, t): Approximate the first time derivative of a snapshot matrix X in which the snapshots are not evenly spaced in time.

pre.xdot(X, *args, **kwargs): Call pre.xdot_uniform() or pre.xdot_nonuniform(), depending on the arguments.

Postprocessing Tools

The post submodule is a collection of common routines for computing the absolute and relative errors produced by a ROM approximation. Given a norm, "true" data X, and an approximation Y to X, these errors are defined by

post.frobenius_error(X, Y): Compute the absolute and relative Frobenius-norm errors between snapshot sets X and Y, assuming Y is an approximation to X. The Frobenius matrix norm is defined by

post.lp_error(X, Y, p=2, normalize=False): Compute the absolute and relative lp-norm errors between snapshot sets X and Y, assuming Y is an approximation to X. The lp norm is defined by

With p = 2 this is the usual Euclidean norm. The errors are calculated for each pair of columns of X and Y. If normalize=True, then the normalized absolute error is computed instead of the relative error:

post.Lp_error(X, Y, t=None, p=2): Approximate the absolute and relative Lp-norm errors between snapshot sets X and Y corresponding to times t, assuming Y is an approximation to X. The Lp norm for vector-valued functions is defined by

For finite p, the integrals are approximated by the trapezoidal rule:

The t argument can be omitted if p is infinity (p = np.inf).

Utility Functions

These functions are helper routines for matricized higher-order tensor operations.

utils.kron2c(x): Compute the compact, quadratic, column-wise (Khatri-Rao) Kronecker product of x with itself.

utils.kron3c(x): Compute the compact, cubic, column-wise (Khatri-Rao) Kronecker product of x with itself three times.

utils.compress_H(H): Convert the full r x r2 matricized quadratic operator H to its compact r x (r(r+1)/2) form.

utils.expand_H(H): Convert the compact r x (r(r+1)/2) matricized quadratic operator H to the full, symmetric, r x r2 form.

utils.compress_G(G): Convert the full r x r3 matricized cubic operator G to its compact r x (r(r+1)(r+2)/6) form.

utils.expand_G(G): Convert the compact r x (r(r+1)(r+2)/6) matricized cubic operator G to the full, symmetric, r x r3 form.