API reference¶
The reference is generated from the in-source NumPy-style docstrings. Only the public entry points are shown; private helpers (leading underscore) are intentionally omitted.
Package root¶
For convenience the package root re-exports the most common entry points so they can
be imported directly as from vane import ...:
Each is documented in its module section below.
I/O¶
Reader for OpenFAST linearization (.lin) files.
This module parses the modern linearization output produced by the OpenFAST
ModGlue engine into structured, in-memory objects. The modern format is
identified by the presence of a Derivative Order column in the
operating-point tables; the legacy format (which lacks that column and instead
emits StateRotation / mass-matrix blocks) is intentionally not supported and
raises LinFileFormatError.
The on-disk layout consumed here is, in order:
A header block beginning at
Simulation information:holding scalar metadata (simulation time, rotor speed in rad/s, azimuth in rad, wind speed, and channel counts) and terminated byJacobians included in this file?.Operating-point tables (
Order of continuous states:and the analogous blocks for state derivatives, inputs, and outputs). Each row carries an index, an operating-point value, a rotating-frame flag, a derivative order, and a free-text description.Matrix blocks (
A,B,C,DunderLinearized state matrices:and, when Jacobians are written,dUdu/dUdy/dXdy/JunderJacobian matrices:).
- class LinFile(path, sim_time, rotor_speed, azimuth, wind_speed, jacobians_included, x, xdot, u, y, a, b, c, d, jacobians)[source]¶
Bases:
objectParsed contents of a single OpenFAST
.linlinearization file.- Parameters:
path (pathlib.Path) – Source file the data was read from.
sim_time (float) – Simulation time of the linearization, in seconds.
rotor_speed (float) – Rotor speed at the operating point, in rad/s (as written by OpenFAST).
azimuth (float) – Rotor azimuth at the operating point, in rad.
wind_speed (float) – Mean hub-height wind speed, in m/s. Zero when InflowWind is absent.
jacobians_included (bool) – Whether the file contains the
Jacobian matricesblock.x (OperatingPointTable) – Operating-point tables for the continuous states, continuous-state derivatives, inputs, and outputs, respectively.
xdot (OperatingPointTable) – Operating-point tables for the continuous states, continuous-state derivatives, inputs, and outputs, respectively.
u (OperatingPointTable) – Operating-point tables for the continuous states, continuous-state derivatives, inputs, and outputs, respectively.
y (OperatingPointTable) – Operating-point tables for the continuous states, continuous-state derivatives, inputs, and outputs, respectively.
a (npt.NDArray[np.float64] or None) – State-space matrices \(A\) (n_x, n_x), \(B\) (n_x, n_u), \(C\) (n_y, n_x), and \(D\) (n_y, n_u).
Nonewhen the corresponding block is absent (e.g.B/Dwith zero inputs).b (npt.NDArray[np.float64] or None) – State-space matrices \(A\) (n_x, n_x), \(B\) (n_x, n_u), \(C\) (n_y, n_x), and \(D\) (n_y, n_u).
Nonewhen the corresponding block is absent (e.g.B/Dwith zero inputs).c (npt.NDArray[np.float64] or None) – State-space matrices \(A\) (n_x, n_x), \(B\) (n_x, n_u), \(C\) (n_y, n_x), and \(D\) (n_y, n_u).
Nonewhen the corresponding block is absent (e.g.B/Dwith zero inputs).d (npt.NDArray[np.float64] or None) – State-space matrices \(A\) (n_x, n_x), \(B\) (n_x, n_u), \(C\) (n_y, n_x), and \(D\) (n_y, n_u).
Nonewhen the corresponding block is absent (e.g.B/Dwith zero inputs).jacobians (dict[str, npt.NDArray[np.float64]]) – Any additional matrices from the
Jacobian matricesblock, keyed by their label (e.g."dUdu","dUdy","dXdy").
- path: Path¶
- xdot: OperatingPointTable¶
- exception LinFileFormatError[source]¶
Bases:
ValueErrorRaised when a
.linfile cannot be parsed as the modern format.
- class OperatingPointTable(values, rotating_frame, derivative_order, descriptions)[source]¶
Bases:
objectOperating-point metadata for one channel group of a
.linfile.A channel group is the set of continuous states, continuous-state derivatives, inputs, or outputs. All array attributes are parallel and share the channel ordering written in the file.
- Parameters:
values (npt.NDArray[np.float64]) – Operating-point value of each channel, shape
(n_channels,). For orientation channels that are written with multiple comma-separated components, only the first component is retained here.rotating_frame (npt.NDArray[np.bool_]) – Per-channel flag, shape
(n_channels,).Truemarks a channel expressed in the blade-rotating frame (the hook used by the MBC3 transform).derivative_order (npt.NDArray[np.int_]) – Per-channel derivative order (0, 1, or 2), shape
(n_channels,).descriptions (list[str]) – Per-channel free-text description, length
n_channels. Each begins with the originating module prefix (e.g."ED ...").
- values: npt.NDArray[np.float64]¶
- rotating_frame: npt.NDArray[np.bool_]¶
- derivative_order: npt.NDArray[np.int_]¶
- read_lin_file(path)[source]¶
Read and parse a modern OpenFAST
.linfile.- Parameters:
path (pathlib.Path) – Path to the
.linfile to read.- Returns:
LinFile – The parsed linearization data.
- Raises:
LinFileFormatError – If the file is empty, truncated, in the unsupported legacy format, or otherwise does not conform to the modern ModGlue layout.
- Return type:
Notes
Numeric fields containing Fortran overflow markers (runs of
*) are parsed asnumpy.nanand a warning is logged; the parse does not fail. Empty matrix blocks are not written by OpenFAST, so the correspondingLinFilematrix attribute isNone.
Readers for OpenFAST model input decks (.fst and ElastoDyn).
These parsers extract only the information VANE needs to drive linearization
post-processing: the linearization configuration block from the primary .fst
file (so the expected set of .lin files is known), the active module
switches and referenced sub-file paths, and the rotor/tower geometry from
ElastoDyn (used to scale mode shapes).
OpenFAST input files use the convention <value(s)> <Keyword> - <description>;
the value(s) precede the keyword on each line. Quoted string values and
multi-value array lines (e.g. LinTimes) are both handled.
- class FastModel(path, comp, files, lin)[source]¶
Bases:
objectSubset of a primary
.fstmodel relevant to linearization.- Parameters:
path (pathlib.Path) – Path to the primary
.fstfile.comp (dict[str, int]) – Active-module switches keyed by name (e.g.
"CompHydro").files (dict[str, pathlib.Path]) – Referenced sub-file paths resolved relative to
path’s directory, keyed by their keyword (e.g."EDFile"). Entries marked"unused"in the deck are omitted.lin (LinearizationConfig) – Parsed linearization configuration.
- path: Path¶
- lin: LinearizationConfig¶
- property is_floating: bool¶
Return a heuristic for a floating platform (has mooring, no SubDyn).
Notes
This is an approximation based on module switches: floating platforms use a mooring system (
CompMooring >= 1) and no SubDyn substructure (CompSub == 0). Detailed system classification (spar/semi/TLP) is the responsibility of the configuration layer, not this reader.
- exception FstFileError[source]¶
Bases:
ValueErrorRaised when an OpenFAST input file is malformed or missing a field.
- class LinearizationConfig(linearize, calc_steady, trim_case, trim_tol, trim_gain, n_lin_times, lin_times, lin_inputs, lin_outputs, lin_out_jac, lin_out_mod)[source]¶
Bases:
objectLinearization settings parsed from a primary
.fstfile.- Parameters:
linearize (bool) – Master linearization switch (
Linearize).calc_steady (bool) – Whether a periodic steady-state trim precedes linearization (
CalcSteady); whenTrue,n_lin_timesazimuth snapshots are produced per operating point.trim_case (int) – Trim target
{1: yaw, 2: torque, 3: pitch}(TrimCase).trim_tol (float) – Rotor-speed convergence tolerance (
TrimTol).trim_gain (float) – Proportional trim gain (
TrimGain).n_lin_times (int) – Number of linearization snapshots (
NLinTimes); the number of.linfiles written per operating point.lin_times (list[float]) – Absolute linearization times in seconds (
LinTimes); meaningful only whencalc_steadyisFalse.lin_inputs (int) – Inputs included
{0: none, 1: standard, 2: all}(LinInputs).lin_outputs (int) – Outputs included
{0: none, 1: from OutList, 2: all}(LinOutputs).lin_out_jac (bool) – Whether full Jacobians are written (
LinOutJac).lin_out_mod (bool) – Whether per-module
.linfiles are written (LinOutMod).
- class TurbineGeometry(num_blades, tip_radius, hub_radius, tower_height, tower_base_height)[source]¶
Bases:
objectRotor and tower geometry parsed from an ElastoDyn input file.
- Parameters:
num_blades (int) – Number of blades (
NumBl).tip_radius (float) – Rotor apex-to-tip distance in metres (
TipRad).hub_radius (float) – Rotor apex-to-blade-root distance in metres (
HubRad).tower_height (float) – Tower-top height above its reference datum in metres (
TowerHt).tower_base_height (float) – Tower-base height above the same datum in metres (
TowerBsHt).
- read_elastodyn_geometry(path)[source]¶
Read rotor and tower geometry from an ElastoDyn input file.
- Parameters:
path (pathlib.Path) – Path to the ElastoDyn
.datfile.- Returns:
TurbineGeometry – The parsed rotor/tower geometry.
- Raises:
FstFileError – If the file cannot be read or a required field is missing.
- Return type:
- read_fst_file(path)[source]¶
Read a primary OpenFAST
.fstfile.- Parameters:
path (pathlib.Path) – Path to the
.fstfile.- Returns:
FastModel – Parsed module switches, referenced file paths, and linearization config.
- Raises:
FstFileError – If the file cannot be read. Missing module switches and linearization settings are tolerated and fall back to their defaults.
- Return type:
Multi-blade coordinate (MBC3 / Coleman) transformation.
OpenFAST writes linearized state matrices in a mixed frame in which the blade degrees of freedom rotate with the rotor, making the system matrices azimuth-dependent. The MBC3 transform projects the three per-blade rotating coordinates onto non-rotating multi-blade coordinates (collective, cosine, and sine), and averaging the transformed matrices over one revolution yields the azimuth-invariant linear time-invariant system used for modal analysis.
This module implements the transform for three-bladed rotors directly from the
parsed LinFile objects of a single operating point
(one file per azimuth). Continuous states are reordered into the canonical
[q2, q2dot, q1] layout (second-order displacements, their velocities, then
first-order states) and blade triplets are detected from the channel
descriptions. The transformed A, B, C, and D matrices are
averaged over azimuth.
Notes
At zero rotor speed the transform reduces to a similarity transform of the state matrix and therefore preserves its eigenvalues; this property is used to verify the implementation.
- class MBCResult(avg_a, avg_b, avg_c, avg_d, state_descriptions, ndof2, ndof1, n_blades, performed_transformation, rotor_speed_rpm, wind_speed, azimuths_deg, mbc_coordinates=<factory>, per_azimuth_a=None)[source]¶
Bases:
objectResult of an MBC3 transform over one operating point’s azimuth sweep.
- Parameters:
avg_a (npt.NDArray[np.float64] or None) – Azimuth-averaged state-space matrices in the non-rotating frame. Each is
Nonewhen the corresponding input matrix was absent. When the transform is performed, rows/columns follow the reordered[q2, q2dot, q1]state layout; otherwise the original order is kept.avg_b (npt.NDArray[np.float64] or None) – Azimuth-averaged state-space matrices in the non-rotating frame. Each is
Nonewhen the corresponding input matrix was absent. When the transform is performed, rows/columns follow the reordered[q2, q2dot, q1]state layout; otherwise the original order is kept.avg_c (npt.NDArray[np.float64] or None) – Azimuth-averaged state-space matrices in the non-rotating frame. Each is
Nonewhen the corresponding input matrix was absent. When the transform is performed, rows/columns follow the reordered[q2, q2dot, q1]state layout; otherwise the original order is kept.avg_d (npt.NDArray[np.float64] or None) – Azimuth-averaged state-space matrices in the non-rotating frame. Each is
Nonewhen the corresponding input matrix was absent. When the transform is performed, rows/columns follow the reordered[q2, q2dot, q1]state layout; otherwise the original order is kept.state_descriptions (list[str]) – State descriptions in the same order as the rows/columns of
avg_a.ndof2 (int) – Number of second-order degrees of freedom.
ndof1 (int) – Number of first-order states.
n_blades (int) – Number of blades inferred from the rotating triplets (3, or 0 when no rotating triplets were found).
performed_transformation (bool) – Whether the MBC3 transform was applied.
Falsewhen no three-blade triplets were detected, in which case the matrices are plain azimuth averages in the original state order.rotor_speed_rpm (float) – Mean rotor speed over the azimuth sweep, in rev/min.
wind_speed (float) – Mean wind speed over the azimuth sweep, in m/s.
azimuths_deg (npt.NDArray[np.float64]) – Sorted rotor azimuths used, in degrees.
mbc_coordinates (list[str]) – Multi-blade coordinate of each state, aligned with
state_descriptions:"collective","cosine", or"sine"for a transformed blade-triplet state, and""for non-rotating states or when the transform was not performed.per_azimuth_a (list[npt.NDArray[np.float64]] or None) – The transformed (non-rotating)
Amatrix at each azimuth, in the same state order asavg_a, retained only whenretain_per_azimuth=True. This physical ensemble underpins the azimuth-spread uncertainty estimate.
- azimuths_deg: npt.NDArray[np.float64]¶
- find_blade_triplets(descriptions, rotating_frame)[source]¶
Group rotating channels into per-blade triplets by description.
For each rotating channel, every blade-index token is located (prose
blade 1, BeamDynBD_1, blade-coordinate/AeroDynb1/B1,PitchBearing1, or a fallback single digit). A sibling matcher is built in which all blade digits are tied to one captured value, so a triplet is confirmed only when the other rotating channels are identical except for a consistent blade number — node, mode, and span indices stay literal and are never mistaken for the blade index, even when a channel carries more than one blade token. The(internal DOF index = ...)parenthetical, which holds a second blade-numbered token, is removed before matching.- Parameters:
- Returns:
list[tuple[int, …]] – Triplets of indices into
descriptions, each ordered by blade number; channels that do not form a complete blade triplet are omitted.- Return type:
list[Triplet]
- mbc3_transform(lin_files, *, omega_dot=0.0, retain_per_azimuth=False)[source]¶
Apply the MBC3 transform and azimuth-average a set of linearizations.
- Parameters:
lin_files (Sequence[LinFile]) – Linearization files for a single operating point, one per azimuth. They must share the same channel layout; they are sorted by azimuth here.
omega_dot (float, optional) – Rotor angular acceleration in rad/s^2, assumed zero (appropriate for a steady-state periodic operating point).
retain_per_azimuth (bool, optional) – Keep the transformed
Amatrix at each azimuth on the result (per_azimuth_a) for azimuth-spread uncertainty analysis. Off by default to save memory.
- Returns:
MBCResult – The azimuth-averaged non-rotating matrices and associated metadata.
- Raises:
ValueError – If
lin_filesis empty, the files have inconsistent state counts, a file is missing theAmatrix, two files share an azimuth, or the rotor speed drifts across the sweep (the files are not a single operating point).- Return type:
Modal analysis¶
Eigenvalue analysis of the azimuth-averaged state matrix.
Given the non-rotating (MBC-averaged) state matrix \(A\) in the canonical
[q2, q2dot, q1] layout, an eigenvalue analysis yields the full-system modes:
their natural frequencies, damping ratios, and (complex) mode shapes.
Each under-damped mode appears as a complex-conjugate eigenvalue pair; only the representative with positive imaginary part is retained. For an eigenvalue \(\lambda\), the modal quantities are
so the natural and damped frequencies in hertz are
\(\omega_n / 2\pi\) and \(\omega_d / 2\pi\). The velocity (q2dot)
rows are dropped from the returned mode shapes, leaving the second-order
displacements and first-order states; the full eigenvectors are also retained.
Notes
OpenFAST exposes no full-system mass or stiffness matrix, so modal frequencies are obtained directly from \(A\) and vary with operating point.
- class ModalSolution(eigenvalues, natural_frequencies_hz, damping_ratios, damped_frequencies_hz, mode_shapes, full_eigenvectors, n_rigid_body_modes, dof_descriptions, condition_numbers=<factory>, is_degenerate=<factory>, n_unstable=0, n_overdamped=0, dof_mbc_coordinates=<factory>)[source]¶
Bases:
objectModal properties extracted from a state matrix.
All per-mode arrays share the same ordering (ascending natural frequency when the solver sorts). Only under-damped modes (positive imaginary eigenvalue) are represented.
- Parameters:
eigenvalues (npt.NDArray[np.complex128]) – Retained eigenvalues \(\\lambda\), shape
(n_modes,).natural_frequencies_hz (npt.NDArray[np.float64]) – Natural frequencies \(|\\lambda| / 2\\pi\), shape
(n_modes,).damping_ratios (npt.NDArray[np.float64]) – Damping ratios \(-\\operatorname{Re}\\lambda / |\\lambda|\).
damped_frequencies_hz (npt.NDArray[np.float64]) – Damped frequencies \(\\operatorname{Im}\\lambda / 2\\pi\).
mode_shapes (npt.NDArray[np.complex128]) – Reduced complex mode shapes (second-order displacements then first-order states), shape
(ndof2 + ndof1, n_modes), phase-normalized so each mode’s dominant component is real and positive.full_eigenvectors (npt.NDArray[np.complex128]) – Full eigenvectors over all states, shape
(n_states, n_modes), sharing the same phase normalization.n_rigid_body_modes (int) – Number of degrees of freedom without an under-damped representative (rigid-body or over-damped modes).
dof_descriptions (list[str]) – Descriptions of the degrees of freedom that index the rows of
mode_shapes(the second-order displacements then first-order states). Empty when descriptions were not supplied to the solver.dof_mbc_coordinates (list[str]) – Multi-blade coordinate (
"collective"/"cosine"/"sine"/"") of eachmode_shapesrow, aligned withdof_descriptions. Empty when not supplied.condition_numbers (npt.NDArray[np.float64]) – Per-mode eigenvalue condition number \(\\kappa = 1 / |y^H x|\) (the reciprocal cosine of the angle between the left and right eigenvectors), shape
(n_modes,).1is perfectly conditioned; large values flag eigenvalues that are numerically sensitive (highly non-normal system).is_degenerate (npt.NDArray[np.bool_]) – Per-mode flag marking modes whose eigenvalue is (near-)repeated by another retained mode, so the individual mode shape is defined only up to a rotation within the shared invariant subspace.
n_unstable (int) – Number of eigenvalues with positive real part (unstable). Non-zero is a red flag for the operating point.
n_overdamped (int) – Number of non-oscillatory stable real eigenvalues (over-damped modes), which are not returned as oscillatory modes.
- eigenvalues: npt.NDArray[np.complex128]¶
- natural_frequencies_hz: npt.NDArray[np.float64]¶
- damping_ratios: npt.NDArray[np.float64]¶
- damped_frequencies_hz: npt.NDArray[np.float64]¶
- mode_shapes: npt.NDArray[np.complex128]¶
- full_eigenvectors: npt.NDArray[np.complex128]¶
- condition_numbers: npt.NDArray[np.float64]¶
- is_degenerate: npt.NDArray[np.bool_]¶
- property natural_frequencies_rad_s: npt.NDArray[np.float64]¶
Return natural frequencies in rad/s.
- compute_modes(a, ndof2, ndof1, *, sort_by_frequency=True, descriptions=None, mbc_coordinates=None)[source]¶
Compute modal properties from a state matrix.
- Parameters:
a (npt.NDArray[np.float64]) – Square state matrix in
[q2, q2dot, q1]order, shape(2 * ndof2 + ndof1, 2 * ndof2 + ndof1).ndof2 (int) – Number of second-order degrees of freedom (the
q2block size).ndof1 (int) – Number of first-order states (the
q1block size).sort_by_frequency (bool, optional) – Sort the returned modes by ascending natural frequency (the default).
descriptions (list[str] or None, optional) – Descriptions of all
2 * ndof2 + ndof1states, in the matrix order. When given, the subset that indexes the mode-shape rows is stored on the result asdof_descriptions.mbc_coordinates (list[str] or None, optional) – Multi-blade coordinate of all states, in the matrix order; the mode-shape subset is stored on the result as
dof_mbc_coordinates.
- Returns:
ModalSolution – The modal frequencies, damping ratios, and mode shapes.
- Raises:
ValueError – If
ais not square or its size does not match2 * ndof2 + ndof1.- Return type:
- modes_from_mbc(result)[source]¶
Compute modal properties from an MBC transform result.
- Parameters:
result (MBCResult) – Output of
vane.io.mbc_transform.mbc3_transform(); itsavg_amust be present.- Returns:
ModalSolution – The modal solution for the averaged state matrix.
- Raises:
ValueError – If
result.avg_aisNone.- Return type:
Modal Assurance Criterion and mode-to-mode matching.
The Modal Assurance Criterion (MAC) measures the correlation between two mode
shapes; a value of 1 indicates perfect correlation and 0 indicates orthogonality
[1]. For complex (state-space) modes the pole-weighted MACXP criterion [2] is
also provided. match_modes performs a one-to-one Hungarian assignment over a MAC
(optionally frequency-penalized) matrix between two operating points; this is a
pairwise building block. Cross-operating-point tracking over a whole sweep is done
separately by vane.modal.identify_modes(), which extracts globally optimal
continuity paths rather than chaining pairwise matches.
References
- compute_mac(phi_ref, phi_test)[source]¶
Compute the Modal Assurance Criterion (MAC) matrix.
- Parameters:
phi_ref (npt.NDArray[np.complex128]) – Reference mode-shape matrix, shape
(n_dof, n_modes_ref).phi_test (npt.NDArray[np.complex128]) – Test mode-shape matrix, shape
(n_dof, n_modes_test).
- Returns:
npt.NDArray[np.float64] – MAC matrix of shape
(n_modes_ref, n_modes_test)with values in[0, 1].- Raises:
ValueError – If the two mode-shape matrices have incompatible first dimensions.
- Return type:
npt.NDArray[np.float64]
Notes
The MAC between modes \(i\) and \(j\) is
\[\text{MAC}_{ij} = \frac{|\phi_i^H \phi_j|^2} {(\phi_i^H \phi_i)(\phi_j^H \phi_j)}.\]
- compute_macxp(phi_ref, lambda_ref, phi_test, lambda_test)[source]¶
Compute the pole-weighted MACXP matrix for complex modes.
MACXP augments the MAC with the transpose product and weights both by the modal poles, making it more robust than MAC for complex state-space modes.
- Parameters:
phi_ref (npt.NDArray[np.complex128]) – Reference and test mode-shape matrices, shape
(n_dof, n_modes).phi_test (npt.NDArray[np.complex128]) – Reference and test mode-shape matrices, shape
(n_dof, n_modes).lambda_ref (npt.NDArray[np.complex128]) – Corresponding modal poles (eigenvalues), shape
(n_modes,).lambda_test (npt.NDArray[np.complex128]) – Corresponding modal poles (eigenvalues), shape
(n_modes,).
- Returns:
npt.NDArray[np.float64] – MACXP matrix of shape
(n_modes_ref, n_modes_test); a mode correlated with itself yields 1.- Raises:
ValueError – If shapes are inconsistent between mode matrices and pole vectors.
- Return type:
npt.NDArray[np.float64]
Notes
The criterion is defined for damped modes (poles with a non-zero real part); an undamped pole makes the pole-weighting denominator zero, so such entries are returned as NaN rather than raising.
- match_modes(mac, ref_frequencies=None, test_frequencies=None, frequency_weight=0.0)[source]¶
Match reference modes to test modes by maximum-weight assignment.
A Hungarian (optimal) assignment maximizes the total correlation between the two mode sets. When frequencies are supplied, the MAC of each pair is scaled by
1 - frequency_weight * |df| / spanso that frequency proximity breaks ties between similar shapes (the strategy used for Campbell-diagram tracking).- Parameters:
mac (npt.NDArray[np.float64]) – Correlation matrix, shape
(n_ref, n_test)(e.g. fromcompute_mac()).ref_frequencies (npt.NDArray[np.float64] or None, optional) – Per-mode natural frequencies for the reference and test sets. Both must be given to apply the frequency penalty.
test_frequencies (npt.NDArray[np.float64] or None, optional) – Per-mode natural frequencies for the reference and test sets. Both must be given to apply the frequency penalty.
frequency_weight (float, optional) – Strength of the frequency penalty in
[0, 1];0(default) uses MAC alone.
- Returns:
list[tuple[int, int]] –
(ref_index, test_index)pairs, one per matched mode (min(n_ref, n_test)pairs), ordered by reference index.- Return type:
Per-degree-of-freedom participation in each mode.
OpenFAST exposes no full-system mass matrix, so a classical mass-normalized modal participation factor is unavailable. Instead, participation is quantified from the right eigenvector: the magnitude of each state’s component in a mode measures how strongly that degree of freedom participates. Each mode is normalized so its dominant component has magnitude one, and a sign is assigned from the phase of each component relative to the dominant one (components more than 90 degrees out of phase are negated), giving the signed contribution used for physical mode identification.
- class ParticipationResult(magnitude, signed_magnitude, phase_deg, dominant_state)[source]¶
Bases:
objectPer-state participation in each mode.
- Parameters:
magnitude (npt.NDArray[np.float64]) – Normalized participation magnitude, shape
(n_dof, n_modes); each column’s maximum is 1.signed_magnitude (npt.NDArray[np.float64]) –
magnitudewith the sign of components more than 90 degrees out of phase with the dominant component flipped negative.phase_deg (npt.NDArray[np.float64]) – Phase of each component in degrees, shape
(n_dof, n_modes).dominant_state (npt.NDArray[np.int_]) – Index of the largest-magnitude state in each mode, shape
(n_modes,).
- magnitude: npt.NDArray[np.float64]¶
- signed_magnitude: npt.NDArray[np.float64]¶
- phase_deg: npt.NDArray[np.float64]¶
- dominant_state: npt.NDArray[np.int_]¶
- compute_participation(mode_shapes, scale_factors=None)[source]¶
Compute per-state participation from complex mode shapes.
- Parameters:
mode_shapes (npt.NDArray[np.complex128]) – Complex mode-shape matrix, shape
(n_dof, n_modes).scale_factors (npt.NDArray[np.float64] or None, optional) – Per-DOF multipliers applied before normalization (e.g.
1/blade_lengthfor translational blade DOFs) so heterogeneous DOFs are comparable. WhenNone, no scaling is applied.
- Returns:
ParticipationResult – Normalized magnitudes, signed magnitudes, phases, and dominant states.
- Raises:
ValueError – If
scale_factorslength does not match the number of DOFs.- Return type:
- participation_from_modes(solution, scale_factors=None)[source]¶
Compute participation from a
ModalSolution’s mode shapes.- Parameters:
solution (ModalSolution) – Modal solution whose
mode_shapesare analysed.scale_factors (npt.NDArray[np.float64] or None, optional) – Per-DOF multipliers, as in
compute_participation().
- Returns:
ParticipationResult – The participation result for the solution’s mode shapes.
- Return type:
Cross-operating-point mode identification (global tracking).
Modes computed independently at successive operating points (rotor speeds or wind speeds) are linked into tracks — the lines of a Campbell diagram. Rather than a greedy adjacent-point match (which can commit early to a locally optimal pairing that is globally wrong), the tracker extracts each mode line as the globally optimal continuity path through the whole sweep: it repeatedly finds the single best-continuity path over the not-yet-assigned modes (a deterministic dynamic program over the modal-assurance-criterion / frequency-continuity affinity) and removes it, until every mode is assigned.
The result is deterministic — no random restarts or preset cluster counts — and each track carries a confidence (its weakest link’s correlation) and an ambiguity flag (raised where a mode had two near-equally correlated continuations, i.e. a crossing or veering hotspot). Every track is given a physical label.
Optimality caveat: each extracted path is the globally optimal single continuity
path over the remaining modes, but the extraction is sequential, not a simultaneous
globally optimal multi-track assignment. In dense spectra or at near-degenerate
crossings an earlier high-scoring path can still claim a node a later path would have
preferred; the is_ambiguous flag marks exactly those hotspots, and confidence
exposes how trustworthy each line is.
- class IdentificationResult(tracks, n_operating_points)[source]¶
Bases:
objectResult of cross-operating-point mode identification.
- Parameters:
- class ModeTrack(operating_points, mode_indices, natural_frequencies_hz, damping_ratios, label, confidence=1.0, is_ambiguous=False)[source]¶
Bases:
objectA mode followed across operating points.
- Parameters:
operating_points (list[int]) – Indices of the operating points this track spans, in order.
mode_indices (list[int]) – The mode index within each corresponding operating point.
natural_frequencies_hz (list[float]) – Natural frequency at each operating point along the track.
damping_ratios (list[float]) – Damping ratio at each operating point along the track.
label (ModeLabel) – Representative physical label (the most confident along the track).
confidence (float) – The weakest link correlation along the track (minimum MAC between consecutive points), in
[0, 1];1.0for a single-point track.is_ambiguous (bool) – Whether any link was a near-tie (a mode with two near-equally correlated continuations), marking a crossing or veering hotspot.
- identify_modes(solutions, *, frequency_weight=0.5, mac_threshold=0.5, ambiguity_margin=0.2)[source]¶
Link modes across operating points into labeled tracks, globally.
- Parameters:
solutions (Sequence[ModalSolution]) – Modal solutions, one per operating point, in operating-point order. All must share the same number of mode-shape degrees of freedom.
frequency_weight (float, optional) – Weight of the frequency-continuity penalty in the link affinity,
[0, 1].mac_threshold (float, optional) – Minimum modal assurance criterion for two modes to be linkable.
ambiguity_margin (float, optional) – A link is flagged ambiguous when the best continuation’s MAC exceeds the second best by less than this margin.
- Returns:
IdentificationResult – The identified tracks and the number of operating points.
- Raises:
ValueError – If
frequency_weight,mac_threshold, orambiguity_marginis not in[0, 1], or the solutions have inconsistent mode-shape DOF layout.- Return type:
Physics-aware labeling of system modes.
Each mode is labeled by the physical degree-of-freedom category that dominates
its shape. Participation magnitudes are aggregated per DofCategory
(via vane.config.dof_map.classify_dof()), the dominant category names the
mode (e.g. “1st Tower Fore-Aft”), and the fraction of participation it carries is
reported as a confidence — a low value flags an ambiguous, mixed mode.
- class ModeLabel(category, label, confidence, dominant_dofs, multiblade=None)[source]¶
Bases:
objectPhysical label for a single mode.
- Parameters:
category (DofCategory) – The dominant DOF category.
label (str) – Human-readable name of the dominant category.
confidence (float) – Fraction of total participation carried by the dominant category, in
[0, 1]; low values indicate a mixed or ambiguous mode.dominant_dofs (list[str]) – Descriptions of the highest-participation degrees of freedom (up to three).
multiblade (str or None) – For a blade mode after the MBC transform, the multi-blade type
"collective","regressive","progressive", or"cyclic"(a standing cyclic mode with no whirl direction);Nonefor non-blade modes or when no multi-blade information is available.
- category: DofCategory¶
- category_to_label(category)[source]¶
Return the human-readable label for a DOF category.
- Parameters:
category (DofCategory) – The category to name.
- Returns:
str – The display label (
"Unidentified"for unknown categories).- Return type:
- label_mode(magnitudes, descriptions, *, multiblade=None)[source]¶
Label one mode from its per-DOF participation magnitudes.
- Parameters:
magnitudes (npt.NDArray[np.float64]) – Non-negative participation magnitude of each DOF, shape
(n_dof,).descriptions (Sequence[str]) – DOF descriptions parallel to
magnitudes.multiblade (str or None, optional) – Multi-blade type (
"collective"/"regressive"/"progressive") for a blade mode; appended to the label and stored on the result.
- Returns:
ModeLabel – The dominant category, its label, the confidence, and the top DOFs.
- Raises:
ValueError – If
magnitudesanddescriptionshave different lengths.- Return type:
- label_modes(participation, descriptions)[source]¶
Label every mode of a participation result.
- Parameters:
participation (ParticipationResult) – Per-DOF participation (its
magnitudecolumns are one mode each).descriptions (Sequence[str]) – DOF descriptions parallel to the participation rows.
- Returns:
list[ModeLabel] – One label per mode (column).
- Return type:
- label_solution(solution)[source]¶
Label every mode of a modal solution using its mode shapes.
- Parameters:
solution (ModalSolution) – A solution whose
dof_descriptionsare populated (i.e. produced with descriptions, as byvane.modal.eigensolver.modes_from_mbc()).- Returns:
list[ModeLabel] – One label per mode.
- Raises:
ValueError – If the solution carries no DOF descriptions.
- Return type:
Uncertainty quantification for modal estimates.
Two physically grounded uncertainty sources are combined into one per-mode confidence, none of which the reference tools provide:
Azimuth spread. The MBC transform produces one non-rotating state matrix per azimuth; their average is the linear-time-invariant model. For an isotropic rotor every per-azimuth matrix is identical (azimuth-invariant to machine precision), so the eigenvalue spread across azimuths is ~0; for an anisotropic rotor the spread grows, quantifying how well the averaged model represents the periodic system.
Numerical conditioning and degeneracy from the eigen-analysis.
The unified confidence multiplies independent reliability factors (each in [0, 1])
so that a mode is trusted only when it is numerically well posed, non-degenerate,
azimuth-stable, and (optionally) confidently tracked. It feeds the covariance
initialization of the state-space export.
- class AzimuthSpread(natural_frequency_std, damping_ratio_std, n_azimuths)[source]¶
Bases:
objectPer-mode variation of the modal estimates across azimuth.
- Parameters:
natural_frequency_std (npt.NDArray[np.float64]) – Standard deviation of each averaged mode’s natural frequency (Hz) across the per-azimuth perturbations of the averaged model.
damping_ratio_std (npt.NDArray[np.float64]) – Standard deviation of each averaged mode’s damping ratio across azimuths.
n_azimuths (int) – Number of azimuth samples used.
- natural_frequency_std: npt.NDArray[np.float64]¶
- damping_ratio_std: npt.NDArray[np.float64]¶
- azimuth_spread(result)[source]¶
Compute the per-mode azimuth spread of a retained MBC result.
The frozen-time matrices
A(psi)are similar across azimuth, so their spectra coincide; what varies is how each azimuth’s matrix perturbs the averaged model. For each averaged mode the first-order eigenvalue perturbationdlambda = y^H (A(psi) - avg_a) x / (y^H x)(left eigenvectory, rightx) is evaluated at every azimuth, and the standard deviation of the resulting natural frequency and damping is returned. It is zero for an isotropic rotor (A(psi) == avg_a) and grows with anisotropy.- Parameters:
result (MBCResult) – An MBC result computed with
retain_per_azimuth=True.- Returns:
AzimuthSpread – The per-mode frequency and damping spread.
- Raises:
ValueError – If the per-azimuth matrices were not retained or
avg_ais absent.- Return type:
- unified_mode_confidence(solution, *, frequency_spread=None, track_confidence=None, spread_scale=0.05)[source]¶
Combine the uncertainty sources into one per-mode confidence in
[0, 1].The confidence is the product of independent reliability factors: numerical conditioning (
1 / kappa), a degeneracy penalty, an azimuth-stability factorexp(-relative_spread / spread_scale), and an optional tracking confidence.- Parameters:
solution (ModalSolution) – The modal solution (provides condition numbers and degeneracy flags).
frequency_spread (npt.NDArray[np.float64] or None, optional) – Per-mode natural-frequency standard deviation (e.g. from
azimuth_spread());nanentries are treated as zero spread.track_confidence (npt.NDArray[np.float64] or None, optional) – Per-mode tracking confidence in
[0, 1]from the cross-operating-point tracker.spread_scale (float, optional) – Relative-spread scale at which the azimuth-stability factor falls to
1/e.
- Returns:
npt.NDArray[np.float64] – The per-mode confidence.
- Return type:
npt.NDArray[np.float64]
Configuration¶
Classification of OpenFAST degree-of-freedom descriptions.
OpenFAST embeds an internal degree-of-freedom (DOF) token in each ElastoDyn and
BeamDyn state description (e.g. DOF_TFA1 for the first tower fore-aft mode or
DOF_BF(1,1) for the first flapwise mode of blade 1). These tokens are stable
across versions and are used here to map a state description onto a physical DOF
category, the originating module, the blade index for rotating DOFs, and whether
the state is a velocity. This classification underpins physics-aware mode
labeling.
- class DofCategory(value)[source]¶
Bases:
EnumPhysical category of a degree of freedom (string-valued, 3.10 compatible).
- PLATFORM_SURGE = 'platform_surge'¶
- PLATFORM_SWAY = 'platform_sway'¶
- PLATFORM_HEAVE = 'platform_heave'¶
- PLATFORM_ROLL = 'platform_roll'¶
- PLATFORM_PITCH = 'platform_pitch'¶
- PLATFORM_YAW = 'platform_yaw'¶
- TOWER_FORE_AFT_1 = 'tower_fore_aft_1'¶
- TOWER_SIDE_SIDE_1 = 'tower_side_side_1'¶
- TOWER_FORE_AFT_2 = 'tower_fore_aft_2'¶
- TOWER_SIDE_SIDE_2 = 'tower_side_side_2'¶
- NACELLE_YAW = 'nacelle_yaw'¶
- GENERATOR_AZIMUTH = 'generator_azimuth'¶
- DRIVETRAIN_TORSION = 'drivetrain_torsion'¶
- ROTOR_FURL = 'rotor_furl'¶
- TAIL_FURL = 'tail_furl'¶
- TEETER = 'teeter'¶
- BLADE_FLAP_1 = 'blade_flap_1'¶
- BLADE_FLAP_2 = 'blade_flap_2'¶
- BLADE_EDGE_1 = 'blade_edge_1'¶
- BLADE_PITCH = 'blade_pitch'¶
- UNKNOWN = 'unknown'¶
- class DofInfo(category, module, blade, is_velocity)[source]¶
Bases:
objectClassification of a single state description.
- Parameters:
category (DofCategory) – The physical DOF category (
UNKNOWNif not recognized).module (str) – The originating module prefix (e.g.
"ED"or"BD_1"), or an empty string if absent.blade (int or None) – The 1-based blade index for rotating blade DOFs, otherwise
None.is_velocity (bool) – Whether the state is the time derivative (velocity) of a DOF.
- category: DofCategory¶
- classify_dof(description)[source]¶
Classify a state description into a
DofInfo.- Parameters:
description (str) – A state description from a linearization file (e.g.
"ED 1st tower fore-aft ... (internal DOF index = DOF_TFA1), m").- Returns:
DofInfo – The classification;
categoryisDofCategory.UNKNOWNwhen no known DOF token is present.- Return type:
Campbell diagram¶
Assembly of Campbell-diagram data from tracked modes.
A Campbell diagram plots modal natural frequency (and damping) against an operating
parameter — rotor speed or wind speed — with one curve per tracked mode. This
module pairs the mode tracks from vane.modal.identify_modes() with the
operating-parameter value of each operating point, producing a structure ready for
plotting or resonance analysis.
- class CampbellDiagram(parameter_values, parameter_name, tracks, n_operating_points)[source]¶
Bases:
objectCampbell-diagram data: tracked modes against an operating parameter.
- Parameters:
parameter_values (npt.NDArray[np.float64]) – Operating-parameter value (e.g. rotor speed in rev/min) at each operating point, indexed by operating-point number, shape
(n_operating_points,).parameter_name (str) – Name of the operating parameter (e.g.
"rotor_speed_rpm").tracks (list[ModeTrack]) – The tracked modes, ordered by first natural frequency.
n_operating_points (int) – Number of operating points.
- parameter_values: npt.NDArray[np.float64]¶
- track_curve(track)[source]¶
Return the operating-parameter, frequency, and damping arrays of a track.
- Parameters:
track (ModeTrack) – A track belonging to this diagram.
- Returns:
tuple of npt.NDArray[np.float64] –
(parameter_values, natural_frequencies_hz, damping_ratios)along the track, each of length equal to the track’s span.- Return type:
tuple[npt.NDArray[np.float64], npt.NDArray[np.float64], npt.NDArray[np.float64]]
- build_campbell(result, parameter_values, *, parameter_name='rotor_speed_rpm')[source]¶
Build a Campbell diagram from identified mode tracks.
- Parameters:
result (IdentificationResult) – The cross-operating-point mode tracks.
parameter_values (Sequence[float]) – Operating-parameter value at each operating point, in operating-point order; its length must equal
result.n_operating_points.parameter_name (str, optional) – Name of the operating parameter (default
"rotor_speed_rpm").
- Returns:
CampbellDiagram – The assembled diagram.
- Raises:
ValueError – If
parameter_valueslength does not match the number of operating points.- Return type:
- campbell_from_solutions(solutions, parameter_values, *, parameter_name='rotor_speed_rpm', frequency_weight=0.5, mac_threshold=0.5)[source]¶
Identify mode tracks from solutions and build a Campbell diagram.
- Parameters:
solutions (Sequence[ModalSolution]) – Modal solutions, one per operating point, in operating-point order.
parameter_values (Sequence[float]) – Operating-parameter value at each operating point.
parameter_name (str, optional) – Name of the operating parameter.
frequency_weight (float, optional) – Passed through to
vane.modal.identify_modes().mac_threshold (float, optional) – Passed through to
vane.modal.identify_modes().
- Returns:
CampbellDiagram – The assembled diagram.
- Return type:
Per-rev excitation lines and resonance detection.
A rotor turning at Omega excites the structure at integer multiples of the
rotational frequency: the nP lines, with frequency n * rpm / 60 Hz. A
potential resonance occurs where a mode’s natural-frequency curve crosses an nP
line. For a three-bladed rotor the dominant excitations are 1P (rotor imbalance)
and the blade-passing 3P, 6P, and 9P lines.
- class ResonanceCrossing(track_index, harmonic, rotor_speed_rpm, frequency_hz, damping_ratio, severity)[source]¶
Bases:
objectA crossing of a mode curve with an
nPexcitation line.- Parameters:
track_index (int) – Index of the crossing mode track within the diagram.
harmonic (int) – The excitation harmonic
n(thenPline).rotor_speed_rpm (float) – Rotor speed at the crossing, interpolated between operating points.
frequency_hz (float) – Frequency at the crossing, in hertz.
damping_ratio (float) – Damping ratio of the resonating mode at the crossing, interpolated between operating points.
severity (ResonanceSeverity) – Damping-aware severity of the crossing.
- severity: ResonanceSeverity¶
- class ResonanceSeverity(value)[source]¶
Bases:
EnumDamping-aware severity of a resonance crossing.
A crossing’s danger depends on how lightly damped the resonating mode is: a near-undamped mode at a crossing can build large response, while a well-damped one is benign.
- HIGH = 'high'¶
- MEDIUM = 'medium'¶
- LOW = 'low'¶
- excitation_frequencies(rotor_speed_rpm, harmonics)[source]¶
Return the
nPexcitation-line frequencies at each rotor speed.- Parameters:
rotor_speed_rpm (npt.NDArray[np.float64]) – Rotor speed in rev/min, shape
(n_points,).harmonics (Sequence[int]) – The harmonics
nto evaluate.
- Returns:
dict[int, npt.NDArray[np.float64]] – Mapping of each harmonic
nto itsnPfrequency (Hz) at every rotor speed,n * rpm / 60.- Raises:
ValueError – If
harmonicscontains a non-positive, fractional, non-finite, or duplicate value.- Return type:
- find_resonances(diagram, harmonics=(1, 3, 6, 9), *, operating_range=None)[source]¶
Find where mode curves cross the
nPexcitation lines.A crossing is detected wherever
natural_frequency - n * rpm / 60changes sign (or is zero) between consecutive operating points along a track; the rotor speed, frequency, and damping ratio at the crossing are linearly interpolated, and each crossing is assigned a damping-awareResonanceSeverity.- Parameters:
diagram (CampbellDiagram) – A Campbell diagram whose operating parameter is rotor speed in rev/min.
harmonics (Sequence[int], optional) – Excitation harmonics to test (default 1P, 3P, 6P, 9P). Any positive-integer family may be supplied.
operating_range (tuple[float, float] or None, optional) – If given, only crossings whose rotor speed falls within this
(min_rpm, max_rpm)operating window are returned.
- Returns:
list[ResonanceCrossing] – The detected crossings, in track-then-harmonic order.
- Raises:
ValueError – If the diagram’s operating parameter is not rotor speed in rev/min,
harmonicscontains a non-positive, fractional, non-finite, or duplicate value, oroperating_rangeis given withmin_rpm > max_rpm.- Return type:
System identification¶
Clean linear-time-invariant state-space representation.
A StateSpace holds the validated A, B, C, D matrices of a
continuous- or discrete-time system together with the state, input, and output
names. The averaged matrices of an MBC transform are exported directly into this
form, and a zero-order-hold discretization is provided for use in a digital twin.
- class StateSpace(a, b=None, c=None, d=None, state_names=<factory>, input_names=<factory>, output_names=<factory>, dt=None)[source]¶
Bases:
objectA linear time-invariant state-space system.
- Parameters:
a (npt.NDArray[np.float64]) – State matrix, shape
(n_states, n_states).b (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices, shapes
(n_states, n_inputs),(n_outputs, n_states),(n_outputs, n_inputs); any may beNone.c (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices, shapes
(n_states, n_inputs),(n_outputs, n_states),(n_outputs, n_inputs); any may beNone.d (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices, shapes
(n_states, n_inputs),(n_outputs, n_states),(n_outputs, n_inputs); any may beNone.state_names (list[str]) – Names of the states, inputs, and outputs; each is either empty or matches the corresponding dimension.
input_names (list[str]) – Names of the states, inputs, and outputs; each is either empty or matches the corresponding dimension.
output_names (list[str]) – Names of the states, inputs, and outputs; each is either empty or matches the corresponding dimension.
dt (float or None) – Sample time;
Nonefor a continuous-time system, otherwise the discrete sample time in seconds.
- Raises:
ValueError – If the matrix shapes or name lengths are inconsistent.
- a: npt.NDArray[np.float64]¶
- is_stable()[source]¶
Return whether the system is asymptotically stable.
- Returns:
bool – For a continuous system, whether every eigenvalue of
Ahas negative real part; for a discrete system, whether every eigenvalue lies strictly inside the unit circle.- Return type:
- discretized(dt)[source]¶
Return a zero-order-hold discretization at sample time
dt.- Parameters:
dt (float) – Sample time in seconds; must be positive.
- Returns:
StateSpace – The discrete-time system;
CandDare unchanged.- Raises:
ValueError – If
dtis not positive or the system is already discrete.- Return type:
- state_space_from_mbc(result)[source]¶
Export the averaged MBC matrices as a continuous-time state-space system.
- Parameters:
result (MBCResult) – An MBC transform result whose
avg_ais present.- Returns:
StateSpace – The continuous-time system with state names from the result.
- Raises:
ValueError – If
result.avg_aisNone.- Return type:
Real modal state-space realization from modal parameters.
A set of under-damped modes is realized as a real, block-diagonal state-space
system. Each mode with eigenvalue \(\lambda = \sigma + j\omega\) contributes a
2x2 block
whose eigenvalues are \(\sigma \pm j\omega\), acting on the real and imaginary parts of the modal coordinate. When mode shapes are supplied, the output matrix maps the modal states to physical responses via \([\,2\operatorname{Re}\phi, \,-2\operatorname{Im}\phi\,]\). The result is a compact reduced-order model of the identified modes.
- modal_state_space(eigenvalues, mode_shapes=None, *, output_names=None)[source]¶
Build a real modal state-space realization.
- Parameters:
eigenvalues (npt.NDArray[np.complex128]) – Modal eigenvalues, shape
(n_modes,).mode_shapes (npt.NDArray[np.complex128] or None, optional) – Complex mode shapes, shape
(n_outputs, n_modes); when given, the output matrix maps modal states to physical responses.output_names (list[str] or None, optional) – Names of the outputs (the rows of
mode_shapes).
- Returns:
StateSpace – A continuous-time system with
2 * n_modesstates.- Raises:
ValueError – If
mode_shapescolumn count does not match the number of eigenvalues.- Return type:
- modal_state_space_from_solution(solution)[source]¶
Build a modal realization from a
ModalSolution.- Parameters:
solution (ModalSolution) – The modal solution whose eigenvalues and mode shapes are realized.
- Returns:
StateSpace – The real modal state-space system, with outputs named by the solution’s DOF descriptions when available.
- Return type:
Kalman covariance initialization, informed by classification uncertainty.
A Kalman filter needs a process-noise covariance Q, a measurement-noise
covariance R, and an initial state covariance P0. When the modes were
classified with uncertainty quantification, that uncertainty is propagated into
P0: a mode identified with low confidence is given a larger initial variance,
so the filter trusts well-identified modes more than ambiguous ones.
- class CovarianceSet(q, r, p0)[source]¶
Bases:
objectKalman-filter covariance matrices.
- Parameters:
q (npt.NDArray[np.float64]) – Process-noise covariance, shape
(n_states, n_states).r (npt.NDArray[np.float64]) – Measurement-noise covariance, shape
(n_outputs, n_outputs).p0 (npt.NDArray[np.float64]) – Initial state covariance, shape
(n_states, n_states).
- q: npt.NDArray[np.float64]¶
- r: npt.NDArray[np.float64]¶
- p0: npt.NDArray[np.float64]¶
- covariances_from_confidence(confidences, n_outputs, *, process_noise=0.0001, measurement_noise=0.01, base_variance=1.0, uncertainty_scale=10.0)[source]¶
Build covariances whose
P0grows with classification uncertainty.Each mode maps to two modal states (real and imaginary parts); the initial variance of those states is
base_variance * (1 + uncertainty_scale * (1 - confidence)), so a low-confidence mode starts with a larger covariance.- Parameters:
confidences (Sequence[float]) – Per-mode classification confidence in
[0, 1]; values are clipped to that range.n_outputs (int) – Number of measured outputs.
process_noise (float, optional) – Diagonal values for
QandR; must be non-negative.measurement_noise (float, optional) – Diagonal values for
QandR; must be non-negative.base_variance (float, optional) – Initial variance of a fully confident mode; must be non-negative.
uncertainty_scale (float, optional) – How strongly uncertainty inflates the initial variance; must be non-negative.
- Returns:
CovarianceSet – The covariance matrices, with
2 * len(confidences)states.- Raises:
ValueError – If a noise value,
base_variance,uncertainty_scale, orn_outputsis negative.- Return type:
- initialize_covariances(n_states, n_outputs, *, process_noise=0.0001, measurement_noise=0.01, initial_variance=1.0)[source]¶
Build uniform diagonal covariance matrices.
- Parameters:
n_states (int) – Number of states.
n_outputs (int) – Number of measured outputs.
process_noise (float, optional) – Diagonal values for
Q,R, andP0respectively; all must be non-negative.measurement_noise (float, optional) – Diagonal values for
Q,R, andP0respectively; all must be non-negative.initial_variance (float, optional) – Diagonal values for
Q,R, andP0respectively; all must be non-negative.
- Returns:
CovarianceSet – The diagonal covariance matrices.
- Raises:
ValueError – If a dimension is negative or a noise value is negative.
- Return type:
Ready-to-use Kalman-filter model assembly.
A KalmanModel bundles a state-space system with its process,
measurement, and initial covariances into the exact matrices a Kalman filter
consumes. The system can be discretized at a chosen sample time during assembly,
and the covariance dimensions are validated against the system so the resulting
model is internally consistent.
- class KalmanModel(a, b, c, d, q, r, p0, dt=None, state_names=<factory>, output_names=<factory>)[source]¶
Bases:
objectA state-space system with Kalman covariances, ready for filtering.
- Parameters:
a (npt.NDArray[np.float64]) – State matrix, shape
(n_states, n_states).b (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices.
c (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices.
d (npt.NDArray[np.float64] or None) – Input, output, and feedthrough matrices.
q (npt.NDArray[np.float64]) – Process-noise, measurement-noise, and initial state covariances.
r (npt.NDArray[np.float64]) – Process-noise, measurement-noise, and initial state covariances.
p0 (npt.NDArray[np.float64]) – Process-noise, measurement-noise, and initial state covariances.
dt (float or None) – Sample time;
Nonefor a continuous-time model.state_names (list[str]) – State and output names carried from the source system.
output_names (list[str]) – State and output names carried from the source system.
- a: npt.NDArray[np.float64]¶
- q: npt.NDArray[np.float64]¶
- r: npt.NDArray[np.float64]¶
- p0: npt.NDArray[np.float64]¶
- build_kalman_model(state_space, covariances, *, dt=None)[source]¶
Assemble a Kalman model from a state-space system and covariances.
- Parameters:
state_space (StateSpace) – The continuous- or discrete-time system.
covariances (CovarianceSet) – The process, measurement, and initial covariances.
dt (float or None, optional) – When given and the system is continuous, discretize it at this sample time before assembly.
- Returns:
KalmanModel – The assembled, dimension-checked model.
- Raises:
ValueError – If
dtis given for an already-discrete system, or the covariance dimensions do not match the system.- Return type:
Visualization¶
Interactive Campbell-diagram plotting.
Renders a CampbellDiagram as an interactive Plotly
figure: one line per tracked mode (natural frequency versus the operating
parameter), the nP excitation lines, and markers at the detected resonance
crossings.
- plot_campbell(diagram, *, harmonics=(1, 3, 6, 9), show_excitation=True, show_resonances=True)[source]¶
Plot a Campbell diagram as an interactive figure.
- Parameters:
diagram (CampbellDiagram) – The diagram to plot.
harmonics (Sequence[int], optional) – Excitation harmonics for the
nPlines and resonance markers (only used for rotor-speed diagrams).show_excitation (bool, optional) – Draw the
nPexcitation lines (rotor-speed diagrams only).show_resonances (bool, optional) – Mark detected resonance crossings (rotor-speed diagrams only).
- Returns:
plotly.graph_objects.Figure – The Campbell-diagram figure.
- Return type:
go.Figure
Damping-ratio plotting against the operating parameter.
Renders the damping ratio of each tracked mode as a function of the operating parameter (rotor speed or wind speed), the companion to the Campbell diagram for assessing stability margins across the operating envelope.
- plot_damping(diagram, *, as_percent=True)[source]¶
Plot damping ratio versus the operating parameter.
- Parameters:
diagram (CampbellDiagram) – The diagram whose tracks are plotted.
as_percent (bool, optional) – Plot damping as a percentage (the default) rather than a fraction.
- Returns:
plotly.graph_objects.Figure – The damping figure, one line per tracked mode.
- Return type:
go.Figure
Static mode-shape plotting.
Since OpenFAST mode shapes are expressed over modal degrees of freedom rather than spatial coordinates, a mode is visualized as a horizontal bar chart of the signed participation of its most active degrees of freedom: bar length is the relative contribution and the sign encodes the in-/out-of-phase relationship to the dominant degree of freedom.
- plot_mode_shape(solution, mode_index, *, top_n=10)[source]¶
Plot the signed participation of a mode’s most active degrees of freedom.
- Parameters:
solution (ModalSolution) – A solution whose
dof_descriptionsare populated.mode_index (int) – Index of the mode to plot.
top_n (int, optional) – Number of highest-participation degrees of freedom to show.
- Returns:
matplotlib.figure.Figure – The bar-chart figure.
- Raises:
ValueError – If the solution has no DOF descriptions,
mode_indexis out of range, ortop_nis not positive.- Return type:
Figure
Three-dimensional complex mode-shape plotting.
A complex mode shape carries both magnitude and phase per degree of freedom. This view renders the most active degrees of freedom as phasors in three dimensions: the horizontal axis orders the degrees of freedom by participation while the other two axes are the real and imaginary parts of each component, so the spread out of the real plane shows how non-proportional (complex) the mode is.
- plot_mode_3d(solution, mode_index, *, top_n=12)[source]¶
Plot a mode’s most active degrees of freedom as 3-D complex phasors.
- Parameters:
solution (ModalSolution) – A solution whose
dof_descriptionsare populated.mode_index (int) – Index of the mode to plot.
top_n (int, optional) – Number of highest-magnitude degrees of freedom to show.
- Returns:
plotly.graph_objects.Figure – The 3-D phasor figure.
- Raises:
ValueError – If the solution has no DOF descriptions,
mode_indexis out of range, ortop_nis not positive.- Return type:
go.Figure
Export¶
Tabular export of modal, Campbell, and uncertainty results.
The analysis products — a per-operating-point modal solution and a Campbell diagram
of tracked modes — are flattened into pandas.DataFrame tables and written to
CSV, JSON, or Excel for downstream tools (spreadsheets, reports, notebooks). The
tables optionally carry the physical labels, per-mode confidence, and azimuth-spread
uncertainty produced by the rest of the library.
- campbell_table(diagram)[source]¶
Flatten a Campbell diagram into a long-format table of tracked modes.
- Parameters:
diagram (CampbellDiagram) – The Campbell diagram to export.
- Returns:
pandas.DataFrame – One row per (track, operating point): the track index and label, the operating-parameter value, natural frequency, damping ratio, track confidence, and ambiguity flag. The operating-parameter column is named after
diagram.parameter_name.- Raises:
ValueError – If
diagram.parameter_namecollides with a fixed export column.- Return type:
pd.DataFrame
- modes_table(solution, *, labels=None, confidence=None, spread=None)[source]¶
Flatten a modal solution into a per-mode table.
- Parameters:
solution (ModalSolution) – The modal solution to export.
labels (Sequence[ModeLabel] or None, optional) – Per-mode physical labels; their name, category, and multi-blade type are added as columns.
confidence (npt.NDArray[np.float64] or Sequence[float] or None, optional) – Per-mode unified confidence, added as a column.
spread (AzimuthSpread or None, optional) – Per-mode azimuth spread; its frequency and damping standard deviations are added as columns.
- Returns:
pandas.DataFrame – One row per mode, ordered as in
solution.- Raises:
ValueError – If any optional argument’s length does not match the number of modes.
- Return type:
pd.DataFrame
- write_table(table, path)[source]¶
Write a table to CSV, JSON, or Excel, inferred from the file suffix.
- Parameters:
table (pandas.DataFrame) – The table to write.
path (str or pathlib.Path) – Destination path; the suffix selects the format (
.csv,.json, or.xlsx).
- Raises:
ValueError – If the suffix is not a supported format, or Excel is requested without the optional
openpyxldependency installed.- Return type:
None
Benchmark¶
Ground-truth benchmarking of modal extraction accuracy and robustness.
A synthetic system is built from prescribed natural frequencies, damping ratios, and mode shapes, so the exact modal answer is known. Running the extraction on it and scoring the result against that ground truth measures accuracy directly; adding controlled noise to the system matrix and re-scoring measures robustness. This is the quantitative evidence that the extraction is accurate and degrades gracefully, which the reference tools do not provide.
A classically damped system x'' + C x' + K x = 0 with M = I is assembled
from a modal matrix Phi and modal parameters via K = Phi diag(wn^2) Phi^-1
and C = Phi diag(2 zeta wn) Phi^-1; its state matrix [[0, I], [-K, -C]] then
has exactly the prescribed eigenvalues and (displacement) mode shapes Phi.
- class GroundTruthSystem(a, natural_frequencies_hz, damping_ratios, mode_shapes, ndof2)[source]¶
Bases:
objectA synthetic system with known modal parameters.
- Parameters:
a (npt.NDArray[np.float64]) – State matrix of shape
(2 n, 2 n)in[displacement, velocity]order.natural_frequencies_hz (npt.NDArray[np.float64]) – The prescribed natural frequencies (Hz), shape
(n,).damping_ratios (npt.NDArray[np.float64]) – The prescribed damping ratios, shape
(n,).mode_shapes (npt.NDArray[np.float64]) – The prescribed displacement mode shapes as columns, shape
(n, n).ndof2 (int) – The number of second-order degrees of freedom
n.
- a: npt.NDArray[np.float64]¶
- natural_frequencies_hz: npt.NDArray[np.float64]¶
- damping_ratios: npt.NDArray[np.float64]¶
- mode_shapes: npt.NDArray[np.float64]¶
- class ModeRecoveryScore(max_frequency_error, max_damping_error, min_mac, mean_mac, n_matched)[source]¶
Bases:
objectAccuracy of a recovered modal solution against a ground truth.
- Parameters:
max_frequency_error (float) – Largest relative natural-frequency error over the matched modes.
max_damping_error (float) – Largest absolute damping-ratio error over the matched modes.
min_mac (float) – Smallest mode-shape MAC over the matched modes (1 = exact);
0when any ground-truth mode was not recovered.mean_mac (float) – Mean mode-shape MAC over all ground-truth modes (unrecovered modes count as
0).n_matched (int) – Number of ground-truth modes matched to a recovered mode (less than the prescribed count when modes are missing).
- robustness_curve(truth, noise_levels, *, seed=0)[source]¶
Score recovery as increasing noise is added to the system dynamics.
One fixed seeded perturbation direction is drawn and, at each level
sigma, scaled bysigma * ||dynamics||and added to the stiffness/damping block of the state matrix (the kinematic[0, I]block is left exact). Reusing a single direction makes the sweep a pure amplitude sweep, so the error grows monotonically withsigma(to first order) rather than varying with an arbitrary fresh direction at each level.- Parameters:
truth (GroundTruthSystem) – The ground-truth system.
noise_levels (Sequence[float]) – Relative noise amplitudes, each
>= 0.seed (int, optional) – Seed for the deterministic perturbation direction.
- Returns:
list[tuple[float, ModeRecoveryScore]] – Each noise level paired with its recovery score.
- Raises:
ValueError – If any noise level is negative or non-finite.
- Return type:
- score_recovery(solution, truth)[source]¶
Score a recovered solution against the ground truth.
Recovered modes are matched to ground-truth modes by maximum total MAC (Hungarian assignment), and the worst-case frequency and damping errors and the mode-shape correlations are reported over the matched pairs.
- Parameters:
solution (ModalSolution) – The extracted modal solution.
truth (GroundTruthSystem) – The ground-truth system.
- Returns:
ModeRecoveryScore – The accuracy metrics.
- Raises:
ValueError – If the solution has no modes to score.
- Return type:
- synthetic_system(natural_frequencies_hz, damping_ratios, *, mode_shapes=None)[source]¶
Build a classically damped system with the prescribed modal parameters.
- Parameters:
natural_frequencies_hz (Sequence[float]) – Natural frequencies (Hz), all strictly positive.
damping_ratios (Sequence[float]) – Damping ratios in
[0, 1), parallel tonatural_frequencies_hz.mode_shapes (npt.NDArray[np.float64] or None, optional) – Invertible modal matrix
Phi(columns are mode shapes), shape(n, n); the identity is used when omitted.
- Returns:
GroundTruthSystem – The state matrix and the exact modal parameters it realizes.
- Raises:
ValueError – If the inputs have mismatched lengths, non-positive frequencies, damping outside
[0, 1), or a non-square / singularmode_shapes.- Return type:
Pipeline¶
High-level modal-analysis pipeline.
ModalPipeline runs the full workflow end to end: for each operating point it
applies the MBC3 transform, solves for modes, then links the modes across
operating points, assembles the Campbell diagram, and flags resonance crossings.
The PipelineResult exposes every intermediate product and convenience
exports to a state-space system, so a complete analysis is a few lines of code.
- class ModalPipeline(*, frequency_weight=0.5, mac_threshold=0.5, harmonics=(1, 3, 6, 9))[source]¶
Bases:
objectEnd-to-end modal analysis from linearization files to a Campbell diagram.
- run(operating_points, *, parameter_name='rotor_speed_rpm')[source]¶
Run the full pipeline over a set of operating points.
- Parameters:
- Returns:
PipelineResult – Every intermediate product of the analysis. Operating points are ordered by ascending operating parameter regardless of input order, so tracking and resonance detection see a monotonic sweep.
- Raises:
ValueError – If
operating_pointsis empty,parameter_nameis unknown, or the operating-parameter values are non-finite or not distinct across points.- Return type:
- class PipelineResult(mbc_results, solutions, identification, campbell, resonances, parameter_name)[source]¶
Bases:
objectOutputs of a full modal-analysis run.
- Parameters:
mbc_results (list[MBCResult]) – The MBC transform result for each operating point.
solutions (list[ModalSolution]) – The modal solution for each operating point.
identification (IdentificationResult) – The cross-operating-point mode tracks.
campbell (CampbellDiagram) – The assembled Campbell diagram.
resonances (list[ResonanceCrossing]) – Detected resonance crossings (empty for non-rotor-speed runs).
parameter_name (str) – Name of the operating parameter used.
- solutions: list[ModalSolution]¶
- identification: IdentificationResult¶
- campbell: CampbellDiagram¶
- resonances: list[ResonanceCrossing]¶
- state_space(operating_point=0)[source]¶
Export the state-space system at one operating point.
- Parameters:
operating_point (int, optional) – Index of the operating point to export (default the first).
- Returns:
StateSpace – The continuous-time system from that operating point’s MBC result.
- Raises:
IndexError – If
operating_pointis negative or out of range.- Return type:
Study and provenance¶
High-level study orchestration with reproducibility provenance.
A Study discovers the linearization files of a campaign, runs the modal
pipeline, and records the provenance of the result — the exact source files (with
content hashes), the azimuth coverage of each operating point, the tuning thresholds,
and the library version — so a Campbell diagram can be tied back to precisely the
inputs and assumptions that produced it. write_bundle serializes a reproducibility
package (a JSON manifest plus the result tables).
- class OperatingPointProvenance(n_azimuths, azimuth_min_deg, azimuth_max_deg, rotor_speed_rpm)[source]¶
Bases:
objectAzimuth-sweep coverage of one operating point.
- Parameters:
- class Provenance(vane_version, created_at, parameter_name, frequency_weight, mac_threshold, harmonics, source_files, operating_points, n_tracks, n_resonances)[source]¶
Bases:
objectA reproducibility record for a study run.
- Parameters:
vane_version (str) – Version of the library that produced the result.
created_at (str) – ISO-8601 timestamp of when the provenance was recorded.
parameter_name (str) – Operating parameter the sweep was run against.
frequency_weight (float) – Tracking tuning thresholds used.
mac_threshold (float) – Tracking tuning thresholds used.
harmonics (tuple[int, ...]) – Excitation harmonics used for resonance detection.
source_files (tuple[SourceFile, ...]) – Every input file with its content hash.
operating_points (tuple[OperatingPointProvenance, ...]) – Per-operating-point azimuth coverage.
n_tracks (int) – Number of identified mode tracks and detected resonance crossings.
n_resonances (int) – Number of identified mode tracks and detected resonance crossings.
- source_files: tuple[SourceFile, ...]¶
- operating_points: tuple[OperatingPointProvenance, ...]¶
- class SourceFile(path, sha256)[source]¶
Bases:
objectA source linearization file and the hash of its contents.
- Parameters:
- class StudyResult(pipeline, provenance)[source]¶
Bases:
objectThe pipeline result of a study together with its provenance.
- Parameters:
pipeline (PipelineResult) – Every intermediate product of the analysis.
provenance (Provenance) – The reproducibility record.
- pipeline: PipelineResult¶
- provenance: Provenance¶
- write_bundle(output_dir)[source]¶
Write a reproducibility bundle to
output_dir.The bundle contains
provenance.json(the manifest of inputs and assumptions) andcampbell.csv(the tracked-mode result table).- Parameters:
output_dir (str or pathlib.Path) – Destination directory; created if it does not exist.
- Return type:
None
- discover_operating_points(directory)[source]¶
Discover and group a directory’s
.linfiles into operating points.Files are named
<case>.<index>.linby OpenFAST; those sharing a<case>root are one operating point’s azimuth sweep.- Parameters:
directory (str or pathlib.Path) – Directory containing
.linlinearization files.- Returns:
list[list[LinFile]] – One list of parsed linearization files per operating point, ordered by case.
- Raises:
FileNotFoundError – If the directory contains no
.linfiles.- Return type:
- run_study(operating_points, *, parameter_name='rotor_speed_rpm', frequency_weight=0.5, mac_threshold=0.5, harmonics=(1, 3, 6, 9), timestamp=None)[source]¶
Run the modal pipeline and record the provenance of the result.
- Parameters:
operating_points (Sequence[Sequence[LinFile]]) – One azimuth sweep per operating point.
parameter_name (str, optional) – Operating parameter to run against.
frequency_weight (float, optional) – Tracking tuning thresholds.
mac_threshold (float, optional) – Tracking tuning thresholds.
harmonics (Sequence[int], optional) – Excitation harmonics for resonance detection.
timestamp (datetime or None, optional) – Timestamp to record (defaults to the current UTC time); injectable for reproducible tests.
- Returns:
StudyResult – The pipeline result and its provenance.
- Return type: