NNFoil
Documentation for NNFoil.
NNFoil.KulfanParameters
NNFoil.KulfanParameters
NNFoil.NeuralNetworkOutput
NNFoil.NeuralNetworkParameters
NNFoil.__evaluate_aero
NNFoil.__flip_x!
NNFoil.airfoil_cst
NNFoil.airfoil_cst_zero_trailing_edge
NNFoil.bernstein
NNFoil.class_function
NNFoil.cst
NNFoil.get_aero_from_kulfan_parameters
NNFoil.get_kulfan_parameters
NNFoil.load_network_parameters
NNFoil.net
NNFoil.normalize_coordinates!
NNFoil.shape_function
NNFoil.sigmoid
NNFoil.split_upper_lower_surfaces
NNFoil.squared_mahalanobis_distance
NNFoil.swish
NNFoil.KulfanParameters
— TypeKulfanParameters{T,V}
Parameter container for the Kulfan (CST) airfoil shape parameterization.
This type stores the upper and lower Bernstein polynomial weights, along with leading-edge and trailing-edge scalar parameters, all sharing a common floating-point element type T
. The coordinate arrays are stored with a concrete vector type V<:AbstractVector{T}
(e.g. Vector{T}
, SVector{N,T}
).
Type parameters
T<:Real
: floating-point element type.V<:AbstractVector{T}
: concrete vector type used for the weight arrays.
Fields
upper_weights::V
: weights for the upper surface (commonly length 8).lower_weights::V
: weights for the lower surface (commonly length 8).leading_edge_weight::T
: scalar parameter controlling leading-edge thickness/rounding.trailing_edge_thickness::T
: scalar trailing-edge thickness parameter.
NNFoil.KulfanParameters
— MethodKulfanParameters(upper_weights, lower_weights, leading_edge_weight,
trailing_edge_thickness)
Converting constructor for KulfanParameters
that promotes all inputs to a common floating-point type T
and returns a KulfanParameters{T, V}
where V
matches the concrete vector type of the provided weights arrays (after promotion).
Arguments
upper_weights::AbstractVector{<:Real}
: upper weights.lower_weights::AbstractVector{<:Real}
: lower weights.leading_edge_weight::Real
: LE scalar parameter.trailing_edge_thickness::Real
: TE thickness parameter.
NNFoil.NeuralNetworkOutput
— TypeNeuralNetworkOutput{V}
Stores the aerodynamic coefficients predicted by the neural network.
Type Parameters
V<:AbstractVector{<:Real}
: numeric vector type used for all output quantities.
Fields
analysis_confidence::V
: confidence level of the neural network prediction.CL::V
: lift coefficient values.CD::V
: drag coefficient values.CM::V
: moment coefficient values.Top_Xtr::V
: transition location on the upper surface.Bot_Xtr::V
: transition location on the lower surface.
NNFoil.NeuralNetworkParameters
— TypeNeuralNetworkParameters{R, V, M, W, B}
Stores parameters of the pretrained neural network model.
Type Parameters
R<:Real
: numeric type used for all elements.V<:AbstractVector{R}
: vector type (for input means and biases).M<:AbstractMatrix{R}
: matrix type (for covariances and weight matrices).W<:AbstractVector{M}
: collection of weight matrices, one per layer.B<:AbstractVector{V}
: collection of bias vectors, one per layer.
Fields
mean_inputs_scaled::V
: mean values of the scaled input features.cov_inputs_scaled::M
: covariance matrix of the scaled inputs.inv_cov_inputs_scaled::M
: inverse of the covariance matrix.weights::W
: vector of weight matrices for each layer.biases::B
: vector of bias vectors for each layer.
NNFoil.__evaluate_aero
— Method__evaluate_aero(network_parameters, x) -> NeuralNetworkOutput
Evaluate the neural network for an input x
using the pretrained parameters network_parameters
, performing symmetry fusion and applying post-processing transformations to produce physically meaningful aerodynamic coefficients.
Arguments
network_parameters
: Pretrained parameters of the neural network.x::AbstractMatrix
: Input features for the original flow condition. Each column corresponds to one input sample.
Returns
A NeuralNetworkOutput
struct with fields:
analysis_confidence
: confidence level (0–1)CL
: lift coefficientCD
: drag coefficientCM
: moment coefficientTop_Xtr
: upper-surface transition location (0–1)Bot_Xtr
: lower-surface transition location (0–1)
NNFoil.__flip_x!
— Method__flip_x!(x)
Flip the input vector or matrix x
in-place, creating a geometrically mirrored version of the input features.
Arguments
x::AbstractVecOrMat{<:Real}
: Input vector or matrix where each column represents a sample. The flipping is applied across specific rows.
Returns
x
: The modified input matrix, flipped in-place.
NNFoil.airfoil_cst
— Methodairfoil_cst(x, parameters, x_split_id; N1=0.5, N2=1.0)
Reconstruct an airfoil surface from Kulfan (CST) parameters.
Arguments
x
: Vector of nondimensional chordwise coordinates (0–1)parameters::AbstractVector
: upper and lower weights, leading-edge weight, trailing-edge thicknessx_split_id::Int
: Index separating upper and lower surface coordinatesN1::Real
: Leading-edge exponent (default: 0.5)N2::Real
: Trailing-edge exponent (default: 1.0)
Returns
Vector
: Airfoil surface y-coordinates corresponding tox
NNFoil.airfoil_cst_zero_trailing_edge
— Methodairfoil_cst_zero_trailing_edge(x, parameters, x_split_id; N1=0.5, N2=1.0)
Reconstruct an airfoil surface from Kulfan (CST) parameters with a zero trailing-edge gap.
Arguments
x
: Vector of nondimensional chordwise coordinates (0–1)parameters::AbstractVector
: upper and lower weights, and leading-edge weightx_split_id::Int
: Index separating upper and lower surface coordinatesN1::Real
: Leading-edge exponent (default: 0.5)N2::Real
: Trailing-edge exponent (default: 1.0)
Returns
Vector
: Airfoil surface y-coordinates corresponding tox
NNFoil.bernstein
— Methodbernstein(x, v, n)
Evaluate the Bernstein basis polynomial of degree n
and index v
at x
.
Arguments
x
: Evaluation points (scalar, vector, or array).v::Signed
: Index of the basis polynomial.n::Signed
: Degree of the polynomial.
Returns
- Array of the same shape as
x
: Values of the Bernstein polynomial.
NNFoil.class_function
— Methodclass_function(x, N1, N2)
Evaluate the class function used in Kulfan’s parametrization.
Arguments
x
: Nondimensional chordwise coordinates [0–1].N1::Real
: Leading-edge exponent.N2::Real
: Trailing-edge exponent.
Returns
- Array of the same shape as
x
: Values of the class function
NNFoil.cst
— Methodcst(x, coefficients, leading_edge_weight, trailing_edge_thickness; N1=0.5, N2=1.0)
CST (Class–Shape Transformation) airfoil parametrization.
Arguments
x
: Nondimensional chordwise coordinates [0, 1]coefficients::AbstractVector
: Shape function weightsleading_edge_weight::Real
: Leading-edge modification termtrailing_edge_thickness::Real
: Trailing-edge thickness parameterN1::Real
: Leading-edge exponent (default: 0.5)N2::Real
: Trailing-edge exponent (default: 1.0)
Returns
- Same shape as
x
: Airfoil surface coordinates defined by the CST parametrization
NNFoil.get_aero_from_kulfan_parameters
— Methodget_aero_from_kulfan_parameters(
network_parameters,
kulfan_parameters,
alpha,
Reynolds;
n_crit=9,
xtr_upper=1,
xtr_lower=1
) -> NeuralNetworkOutput
Compute aerodynamic coefficients from Kulfan airfoil parameters using a pretrained neural network.
This function supports scalar and vector inputs for both the angle of attack alpha
and the Reynolds
number. If either is a scalar while the other is a vector, the scalar will be broadcasted to match the length of the vector.
Arguments
network_parameters::NeuralNetworkParameters
: Pretrained neural network parameters.kulfan_parameters::KulfanParameters
: Kulfan shape parameters describing the airfoil.alpha
: Angle(s) of attack in degrees (Real
orAbstractVector{<:Real}
).Reynolds
: Reynolds number(s) corresponding to eachalpha
(Real
orAbstractVector{<:Real}
).n_crit::Real=9
: Critical amplification factor for transition prediction.xtr_upper::Real=1
: Forced transition location on the upper surface.xtr_lower::Real=1
: Forced transition location on the lower surface.
Returns
NeuralNetworkOutput
: Predicted aerodynamic coefficients.
NNFoil.get_kulfan_parameters
— Methodget_kulfan_parameters(coordinates; num_coefficients=8, N1=0.5, N2=1.0)
Fit Kulfan (CST) parameters to airfoil coordinates.
Arguments
coordinates::AbstractMatrix
: Airfoil coordinates with columns[x, y]
num_coefficients::Int
: Number of Bernstein polynomial coefficients per surface (default: 8)N1::Real
: Leading-edge exponent (default: 0.5)N2::Real
: Trailing-edge exponent (default: 1.0)
Returns
KulfanParameters
: Fitted upper and lower weights, leading-edge weight, and trailing-edge thickness
NNFoil.load_network_parameters
— Methodload_network_parameters(; model_size=:xlarge, T=Float64)
Load and convert the pretrained neural network parameters used by NeuralFoil.
Arguments
model_size::Symbol=:xlarge
: Size of the pretrained model parameters to load.T::Type=Float64
: Numerical type to which all loaded arrays will be converted.
NNFoil.net
— Methodnet(network_parameters::NetworkParameters, x::AbstractMatrix{<:Real})
Evaluate the neural network using the pretrained network parameters on the given input x
.
Arguments
network_parameters::NetworkParameters
: pretrained network weights and biases.x::AbstractArray{<:Real}
: Input data.
Returns
AbstractMatrix{<:Real}
NNFoil.normalize_coordinates!
— Methodnormalize_coordinates!(coordinates)
Normalize the input coordinates in place so that the x values lie within the unit interval [0, 1].
Arguments
coordinates::AbstractMatrix
: Matrix of airfoil coordinates with columns representing the x and y values.
NNFoil.shape_function
— Methodshape_function(x, coefficients)
Kulfan shape function defined as a weighted sum of Bernstein polynomials.
Arguments
x
: Nondimensional chordwise coordinates [0, 1].coefficients::AbstractVector
: Weights for the Bernstein polynomials.
Returns
- Same shape as
x
: Values of the shape function.
NNFoil.sigmoid
— Methodsigmoid(x::T; ln_eps=log(10 / floatmax(eltype(T)))) where {T}
Apply the sigmoid activation function elementwise. Also employ input clipping for numerical stability.
Arguments
x::T
: Input value.ln_eps::Real=log(10 / floatmax(eltype(T)))
: Logarithmic bound used to clip input values for stability.
Returns
- Scalar or Array of the same shape as
x
: Values in the range (0, 1).
NNFoil.split_upper_lower_surfaces
— Methodsplit_upper_lower_surfaces(coordinates)
Split airfoil coordinates into upper and lower surfaces.
Arguments
coordinates::AbstractMatrix
: Matrix of airfoil coordinates with columns representing the x and y values.
Returns
(upper, lower)
: Two matrices containing the coordinates of the upper and lower surfaces, respectively.
NNFoil.squared_mahalanobis_distance
— Methodsquared_mahalanobis_distance(network_parameters::NetworkParameters,
x::AbstractMatrix{<:Real})
Compute the squared Mahalanobis distance between input samples and the mean of the scaled input distribution.
Arguments
network_parameters::NetworkParameters
: pretrained neural network parameters containing the mean and inverse covariance of the scaled input distribution.x::AbstractArray{<:Real}
: Input samples.
Returns
AbstractArray{<:Real}
NNFoil.swish
— Methodswish(x; beta=1)
Apply the Swish activation function elementwise.
Arguments
x
: Input value.beta::Real=1
: Slope parameter controlling smoothness.
Returns
- Scalar or Array of the same shape as
x
: Activated values.