4.2.4. Solver
The OpenFAST tight-coupling solver is implemented in
modules/openfast-library/src/FAST_Solver.f90. It integrates the continuous
states and resolves the input-output coupling between modules using a
generalized-alpha scheme with Newton-Raphson convergence iterations.
4.2.4.1. User input parameters
All solver parameters are set in the main OpenFAST input file
(*.fst) under the Feature Switches and Flags and
Tight-Coupling / Solver sections.
Parameter |
Type |
Description |
|---|---|---|
|
real |
Global (solver) time step in seconds. All module time steps must be
equal to or an integer sub-divisor of |
|
integer |
Coupling method.
|
|
real |
Numerical damping parameter ρ∞ for the generalized-alpha integrator.
Range [0, 1]; 1 = no numerical damping (second-order accurate), 0 =
maximum damping (first-order accurate). Typical value: 0.9.
Reducing |
|
integer |
Maximum number of Newton convergence iterations per time step before the
solver declares convergence failure. Typical value: 20.
With |
|
real |
Convergence tolerance. The iteration stops when the average
L2-norm of the Newton update vector falls below this value.
Typical value: |
|
real |
Time interval (seconds) between Jacobian rebuilds when
|
|
real |
Conditioning scale factor applied to load rows and columns of the Jacobian. Force and moment variables are divided by this factor before the linear solve and multiplied back afterwards, equalising the magnitude of load entries relative to displacement/velocity entries. Typical value: 1.0e5 for offshore systems; may need adjustment for very large or very small turbines. |
|
integer |
Select the structural dynamics module: |
|
integer |
Sub-structural module: |
|
integer |
|
|
integer |
|
|
integer |
Aerodynamics module: |
|
integer |
Controller module: |
4.2.4.2. Generalized-alpha integration
The tight-coupling solver integrates second-order ODEs of the form
using the generalized-alpha method (Chung & Hulbert, 1993). Given the
spectral radius ρ∞ specified by RhoInf, the method parameters are:
Two derived coefficients used throughout the convergence loop are:
where h = DT.
State vector layout – the solver maintains a per-module generalized coordinate (q) vector with four columns:
Column |
Meaning |
|---|---|
|
Displacement / orientation states ( |
|
Velocity states ( |
|
Acceleration (physical, from module |
|
Algorithmic acceleration (generalized-alpha internal variable) |
State prediction at the start of each step:
4.2.4.3. Module ordering
During FAST_SolverInit each module is categorised based on ModCoupling
and its own physics type, and assigned to one of the ordered index arrays
in the Glue_TCParam structure:
Array |
Modules (in order) |
|---|---|
|
ElastoDyn, BeamDyn, SubDyn (when |
|
ServoDyn (when StC active), SED, AD (MHK), ExtPtfm, HydroDyn, OrcaFlex,
MoorDyn; ED/BD/SD also appear here when |
|
ServoDyn, SED, ED, BD, SD, InflowWind, SeaState, AeroDyn (land), AeroDisk, ExtLoads, MAP++, FEAMooring, IceDyn, IceFloe, SoilDyn |
|
ServoDyn, ExternalInflow |
|
SED, ED, BD, SD, InflowWind, ExtLoads (Step 0 initialisation only) |
4.2.4.4. Jacobian construction
Two separate Jacobians are assembled:
TC/Option-1 Jacobian (
BuildJacobianTC) — for the main time-stepping convergence loop.IO Jacobian (
BuildJacobianIO) — for the initial and linearization input-output solve.
4.2.4.4.1. Variable selection (VF_Solve flag)
During FAST_SolverInit → SetVarSolveFlags, the VF_Solve flag is set on
the variables that must appear in the Jacobian:
Continuous states of all TC modules (automatically).
Motion mesh inputs/outputs of TC-to-TC mappings (all fields).
Motion mesh input accelerations of TC-to-Option1 or Option1-to-TC mappings.
Load mesh inputs and outputs involved in any TC/Option1 mapping.
Load mesh displacement outputs of the destination module when the mapping carries moments (needed for moment-arm Jacobian terms).
Variable-to-variable mapped inputs/outputs of TC/Option1 modules.
Any variable with
VF_NoLinis excluded fromVF_Solve.
4.2.4.4.2. Jacobian structure (TC Jacobian)
The assembled TC Jacobian J has size NumJ × NumJ, where:
The columns and rows are partitioned as:
where
J₁₁ (
NumQ × NumQ) — derivative of the acceleration residual with respect to TC displacement/velocity states (formed from the moduledXdxsub-Jacobians plus the generalized-alpha tangent).J₁₂ (
NumQ × NumU_T) — derivative of the acceleration residual with respect to TC inputs (fromdXdu).J₂₁ (
NumU_T × NumQ) — derivative of the input residual with respect to TC states (fromdUdx = dUdy · dydx).J₂₂ (
NumU × NumU) — derivative of the input residual with respect to inputs, including load conditioning rows/columns.
The right-hand side (XB) contains the residuals:
State residual (rows
iJX): difference between the predicted velocity derivative and the module-computed accelerations.Input residual (rows
iJU): difference between the inputs computed from mesh mappings (FAST_InputSolve) and the current iterate.
The loads portion (rows iJL) is pre-divided by UJacSclFact before the
factorisation to improve conditioning.
4.2.4.4.3. Jacobian update strategy
ModCoupling = 2(fixed updates)The Jacobian is rebuilt if either of these counters reaches zero:
UJacStepsRemain— steps remaining; initialised toCEILING(DT_UJac/DT)each time the Jacobian is rebuilt.UJacIterRemain— iteration budget; initialised toCEILING(DT_UJac/DT · MaxConvIter)whenDT_UJac < DT.
On convergence failure the solver returns a fatal error immediately.
ModCoupling = 3(adaptive updates)The Jacobian is rebuilt the first time the convergence loop fails. If the step still does not converge after the forced rebuild, a non-fatal warning is issued and the simulation proceeds.
4.2.4.4.4. Per-module Jacobian contributions
The module-level Jacobian sub-matrices are computed by finite differencing
inside BuildJacobianTC and BuildJacobianIO using the MV_Perturb /
MV_ComputeDiff / MV_ComputeCentralDiff utilities from ModVar. For
each variable flagged VF_Solve:
Apply a positive perturbation of magnitude
Var%Perturbto the working state/input array.Call
FAST_CalcOutput(orFAST_GetContStateDeriv).Apply an equal negative perturbation.
Call again.
Compute the central difference:
(y_plus - y_minus) / (2·Perturb).
For orientation variables (FieldOrientation), perturbations are applied by
quaternion composition rather than direct addition (MV_Perturb), and
differences are extracted as rotation vectors (MV_ComputeDiff).
4.2.4.4.5. Linear solve
The LU factorisation of J is computed with LAPACK_getrf and the system
is solved with LAPACK_getrs (packed in NWTC_LAPACK). The same
factored matrix is reused across convergence iterations until the update
strategy triggers a rebuild.
4.2.4.4.6. Convergence check
After each Newton step the convergence error is the average L2-norm of the update vector:
where \(\Delta \mathbf{z}\) combines state and input updates. The loop
exits if e < ConvTol (ErrID_None) or the iteration count reaches
MaxConvIter (ErrID_Fatal or ErrID_Warn depending on
ModCoupling).
4.2.4.5. Output channels from the solver
Three output channels are written to DriverWriteOutput each step and
appear in the output file when enabled:
Index |
Content |
|---|---|
1 |
Total convergence iterations in the step ( |
2 |
Final convergence error ( |
3 |
Number of Jacobian rebuilds in the step ( |