ShiftScaleTransformer#
- class ShiftScaleTransformer(centering: bool = False, scaling: str = None, byrow: bool = False, name: str = None, verbose: bool = False)[source]#
Process snapshots by centering and/or scaling (in that order).
Transformations with this class are notated below as
\[\Q \mapsto \Q' ~\text{(centered)}~ \mapsto \Q'' ~\text{(centered/scaled)},\]where \(\Q\in\RR^{n \times k}\) is the snapshot matrix to be transformed and \(\Q''\in\RR^{n \times k}\) is the transformed snapshot matrix.
All transformations with this class are affine and hence can be written componentwise as \(\Q_{i,j}'' = \alpha_{i,j} \Q_{i,j} + \beta_{i,j}\) for some choice of \(\alpha_{i,j},\beta_{i,j}\in\RR\).
- Parameters
- centeringbool
If
True
, shift the snapshots by the mean training snapshot, i.e.,\[\Q'_{:,j} = \Q_{:,j} - \frac{1}{k}\sum_{j=0}^{k-1}\Q_{:,j}.\]Otherwise, \(\Q' = \Q\) (default).
- scalingstr or None
If given, scale (non-dimensionalize) the centered snapshot entries. Otherwise, \(\Q'' = \Q'\) (default).
Options:
'standard'
: standardize to zero mean and unit standard deviation.\[\Q'' = \Q' - \frac{\mean(\Q')}{\std(\Q')}.\]Guarantees \(\mean(\Q'') = 0\) and \(\std(\Q'') = 1\).
If
byrow=True
, then \(\mean_{j}(\Q_{i,j}'') = 0\) and \(\std_j(\Q_{i,j}'') = 1\) for each row index \(i\).'minmax'
: minmax scaling to \([0, 1]\).\[\Q'' = \frac{\Q' - \min(\Q')}{\max(\Q') - \min(\Q')}.\]Guarantees \(\min(\Q'') = 0\) and \(\max(\Q'') = 1\).
If
byrow=True
, then \(\min_{j}(\Q_{i,j}'') = 0\) and \(\max_{j}(\Q_{i,j}'') = 1\) for each row index \(i\).'minmaxsym'
: minmax scaling to \([-1, 1]\).\[\Q'' = 2\frac{\Q' - \min(\Q')}{\max(\Q') - \min(\Q')} - 1.\]Guarantees \(\min(\Q'') = -1\) and \(\max(\Q'') = 1\).
If
byrow=True
, then \(\min_{j}(\Q_{i,j}'') = -1\) and \(\max_{j}(\Q_{i,j}'') = 1\) for each row index \(i\).'maxabs'
: maximum absolute scaling to \([-1, 1]\) without scalar mean shift.\[\Q'' = \frac{1}{\max(\text{abs}(\Q'))}\Q'.\]Guarantees \(\mean(\Q'') = \frac{\mean(\Q')}{\max(\text{abs}(\Q'))}\) and \(\max(\text{abs}(\Q'')) = 1\).
If
byrow=True
, then \(\mean_{j}(\Q_{i,j}'') = \frac{\mean_j(\Q_{i,j}')}{\max_j(\text{abs}(\Q_{i,j}'))}\) and \(\max_{j}(\text{abs}(\Q_{i,j}'')) = 1\) for each row index \(i\).'maxabssym'
: maximum absolute scaling to \([-1, 1]\) with scalar mean shift.\[\Q'' = \frac{\Q' - \mean(\Q')}{\max(\text{abs}(\Q' - \mean(\Q')))}.\]Guarantees \(\mean(\Q'') = 0\) and \(\max(\text{abs}(\Q'')) = 1\).
If
byrow=True
, then \(\mean_j(\Q_{i,j}'') = 0\) and \(\max_j(\text{abs}(\Q_{i,j}'')) = 1\) for each row index \(i\).
- byrowbool
If
True
, scale each row of the snapshot matrix separately when a scaling is specified. Otherwise, scale the entire matrix at once.- verbosebool
If
True
, print information upon learning a transformation.
Notes
A custom shifting vector (i.e., the mean snapshot) can be specified by setting the
mean_
attribute. Similarly, the scaling \(q'\mapsto q'' = \alpha q' + \beta\) can be adjusted by setting thescale_
(\(\alpha\)) andshift_
(\(\beta\)) attributes. However, callingfit()
orfit_transform()
will overwrite all three attributes.Properties
byrow
If
True
, scale each row of the snapshot matrix separately.centering
If
True
, center the snapshots by the mean training snapshot.mean_
Mean training snapshot.
name
Label for the state variable that this transformer acts on.
scale_
Multiplicative factor of the scaling, the \(\alpha\) of \(q'' = \alpha q' + \beta\).
scaling
Type of scaling (non-dimensionalization).
shift_
Additive factor of the scaling, the \(\beta\) of \(q'' = \alpha q' + \beta\).
state_dimension
Dimension \(n\) of the state.
verbose
If
True
, print information upon learning a transformation.Methods
Learn (but do not apply) the transformation.
Learn and apply the transformation.
Apply the inverse of the learned transformation.
Load a previously saved transformer from an HDF5 file.
Save the transformer to an HDF5 file.
Apply the learned transformation.
Apply the learned transformation to snapshot time derivatives.
Verify that
transform()
andinverse_transform()
are consistent and thattransform_ddts()
, if implemented, is consistent withtransform()
.