Legendre
Exposed Functions
FiniteHilbertTransform.LegendreFHT
— TypeLegendreFHT
Type representing Legendre Finite Hilbert Transform (FHT) parameters.
This struct stores the necessary parameters for performing a Finite Hilbert Transform using Legendre functions. It implements the AbstractFHT interface and provides the required data structures and computations for Legendre FHT.
Fields:
name::String
: FHT name (default: "Legendre").Ku::Int64
: Number of sample points.tabu::Array{Float64,1}
: Array containing u values (sampling points).tabw::Array{Float64,1}
: Array containing w values (weights at sampling points).tabP::Matrix{Float64}
: Matrix containing P_k(u) values (Ku x Ku).tabc::Vector{Float64}
: Vector containing prefactor values at each sampling point.tabPLeg::Array{ComplexF64,1}
: Static container for tabPLeg (Legendre functions of the first kind).tabQLeg::Array{ComplexF64,1}
: Static container for tabQLeg (Hilbert-transformed Legendre functions).tabDLeg::Array{ComplexF64,1}
: Static container for tabDLeg (Derivatives of Legendre functions).
Example:
Ku = 100
tabu = collect(-1.0:2/(Ku-1):1.0)
tabw = compute_weights(tabu) # Compute weights for Legendre functions
tabP = compute_legendre_matrix(tabu) # Compute Legendre functions matrix
tabc = compute_prefactors(tabu) # Compute Legendre prefactors
FHT = LegendreFHT("Legendre", Ku, tabu, tabw, tabP, tabc, zeros(ComplexF64, Ku), zeros(ComplexF64, Ku), zeros(ComplexF64, Ku))
FiniteHilbertTransform.GettabD!
— Functionfill fht
at a given complex frequency for the integration being considered
integration style selection is automatic: if you want to specify a type, call out to the specific integration method.
GettabD!(omg::Complex{Float64}, struct_tabLeg::LegendreFHT; verbose::Int64=0)
Fill struct_tabLeg
at a given complex frequency omg
for the integration being considered.
Arguments
omg::Complex{Float64}
: Complex frequency value.struct_tabLeg::LegendreFHT
: LegendreFHT structure to fill.verbose::Int64
: Verbosity level (default: 0).
Description
GettabD!
automatically selects the integration style based on the imaginary part of omg
. If the imaginary part is negative, it uses damped integration. If it's exactly zero, it uses neutral mode calculation. Otherwise, it uses unstable integration.
Example
tabLeg = LegendreFHT(10)
omg = 1.0 + 0.5im
GettabD!(omg, tabLeg, verbose=2)
@IMPROVE, verbose creates allocations.
FiniteHilbertTransform.Heaviside
— FunctionHeaviside(x::Float64)
Calculate the Heaviside function on the interval [-1, 1].
The Heaviside function, denoted as H(x), is defined as follows:
- H(x) = 0 for x < -1
- H(x) = 0.5 for x = -1
- H(x) = 1 for -1 < x < 1
- H(x) = 0.5 for x = 1
- H(x) = 0 for x > 1
ATTENTION: The equality tests on Float64 might not be very robust.
Arguments:
x::Float64
: Real number for which Heaviside function is calculated.
Returns:
Float64
: Value of the Heaviside function at x.
Example:
x = 0.5
result = Heaviside(x) # Returns 1.0
FiniteHilbertTransform.tabQLeg!
— FunctiontabQLeg!(omg::ComplexF64, val_0::ComplexF64, val_1::ComplexF64, tabQLeg::Array{ComplexF64,1})
Precompute Hilbert-transformed Legendre functions for a given complex frequency.
Arguments:
omg::ComplexF64
: Complex frequency. (ATTENTION: Must be complex.)val_0::ComplexF64
: Initial value fork = 0
. (ATTENTION: Must be complex.)val_1::ComplexF64
: Initial value fork = 1
. (ATTENTION: Must be complex.)tabQLeg::Array{ComplexF64,1}
: Container to store the precomputed Hilbert-transformed Legendre functions.
Details:
This function uses different recurrence relations based on the location of the complex frequency omg
. If omg
is sufficiently close to the real line [-1,1]
, it employs an upward recurrence. Otherwise, if omg
is far away from the real line [-1,1]
, it uses a backward recurrence. The transition between these regimes is determined dynamically.
The transition from the two regimes follows from the thesis Stable Implementation of Three-Term Recurrence Relations, Pascal Frederik Heiter, June, 2010 https://www.uni-ulm.de/fileadmin/websiteuniulm/mawi.inst.070/funken/bachelorarbeiten/bachelorthesis_pfh.pdf
Example:
omg = 1.0 + 2.0im
val_0 = 1.0 + 1.0im
val_1 = 2.0 + 2.0im
Ku = 10
tabQLeg = zeros(ComplexF64, Ku)
tabQLeg!(omg, val_0, val_1, tabQLeg)
FiniteHilbertTransform.tabPLeg!
— FunctiontabPLeg!(omg::Complex{Float64}, val_0::Complex{Float64}, val_1::Complex{Float64}, Ku::Int64, tabPLeg::Vector{Complex{Float64}})
Pre-computes the Legendre functions, P_k(w), for a given complex frequency.
Arguments
omg::Complex{Float64}
: Complex frequency.val_0::Complex{Float64}
: Initial value in k=0.val_1::Complex{Float64}
: Initial value in k=1.Ku::Int64
: Number of Legendre modes.tabPLeg::Vector{Complex{Float64}}
: Container to store the results.
Description
tabPLeg!
pre-computes the Legendre functions, P_k(w), for a given complex frequency omg
using the upward recurrence method. The results are stored in tabPLeg
.
Example
omg = 1.0 + 0.5im
val_0 = 0.0 + 0.0im
val_1 = 1.0 + 0.0im
Ku = 10
tabPLeg = zeros(Complex{Float64}, Ku)
tabPLeg!(omg, val_0, val_1, Ku, tabPLeg)
FiniteHilbertTransform.GetaXi!
— FunctionGetaXi!(FHT::LegendreFHT, tabGXi::AbstractVector{Float64}, res::Vector{Float64}, warnflag::Vector{Float64})
Compute the Finite Hilbert Transform for Legendre functions.
Arguments
FHT::LegendreFHT
: An object representing the Legendre Finite Hilbert Transform.tabGXi::AbstractVector{Float64}
: Vector containing precomputed values of G[u_i].res::Vector{Float64}
: Output vector to store the results of the Finite Hilbert Transform.warnflag::Vector{Float64}
: Vector to store warning flags for each Legendre function.
Details
This function computes the Finite Hilbert Transform for Legendre functions based on the provided precomputed values of G[u_i].
Output
res::Vector{Float64}
: Vector containing the results of the Finite Hilbert Transform for Legendre functions.warnflag::Vector{Float64}
: Vector containing warning flags. Each element represents the number of NaN or INF contributions for the corresponding Legendre function.
Example
FHT = LegendreFHT(parameters)
tabGXi = compute_tabGXi(...) # Precompute tabGXi values
res = zeros(Float64, FHT.Ku)
warnflag = zeros(Float64, FHT.Ku)
GetaXi!(FHT, tabGXi, res, warnflag)
compatibility signature
FiniteHilbertTransform.GetaXi
— Function@IMPROVE: how do we handle bad G values?
GetaXi(FHT::LegendreFHT, tabGXi::Array{Float64})
Calculate the Finite Hilbert Transform for Legendre functions.
This function computes the Finite Hilbert Transform for Legendre functions based on the provided Legendre Finite Hilbert Transform parameters and precomputed values of G[u_i]. The results are stored in a vector, and warning flags indicating NaN or INF contributions are also provided.
Arguments:
FHT::LegendreFHT
: An object representing the Legendre Finite Hilbert Transform.tabGXi::Array{Float64}
: Array containing precomputed values of G[u_i].
Returns:
res::Vector{Float64}
: Vector containing the results of the Finite Hilbert Transform for Legendre functions.warnflag::Vector{Float64}
: Vector containing warning flags. Each element represents the number of NaN or INF contributions for the corresponding Legendre function.
Example:
FHT = LegendreFHT(parameters)
tabGXi = compute_tabGXi(...) # Precompute tabGXi values
res, warnflag = GetaXi(FHT, tabGXi)
Implementation of Landau Integral
FiniteHilbertTransform.tabLeg!_UNSTABLE
— FunctiontabLeg!_UNSTABLE(omg::Complex{Float64}, struct_tabLeg::LegendreFHT)
Computes the Legendre coefficients D_k(w) using the UNSTABLE integration method.
Arguments
omg::Complex{Float64}
: Complex frequency.struct_tabLeg::LegendreFHT
: LegendreFHT structure.
Description
tabLeg!_UNSTABLE
computes the Legendre coefficients Dk(w) using the UNSTABLE integration method. It computes the coefficients by first calculating the coefficients Qk(w) and then setting Dk(w) equal to Qk(w) for each k.
Notes
- This function assumes that the container for Dk(w) (`structtabLeg.tabDLeg
) and Q_k(w) (
structtabLeg.tabQLeg) are already initialized in the
structtabLeg`.
FiniteHilbertTransform.tabLeg!_NEUTRAL
— FunctiontabLeg!_NEUTRAL(omg::Complex{Float64}, struct_tabLeg::LegendreFHT)
Fill in all the Legendre arrays for a NEUTRAL mode, i.e., Im[w] = 0.0.
Arguments
omg::Complex{Float64}
: Complex frequency.struct_tabLeg::LegendreFHT
: LegendreFHT structure.
Description
tabLeg!_NEUTRAL
fills in all the Legendre arrays for a NEUTRAL mode, where Im[w] = 0.0. It computes the Legendre coefficients Dk(w) and, if needed, Legendre polynomials Pk(w) for the given complex frequency.
FiniteHilbertTransform.tabLeg!_DAMPED
— FunctiontabLeg!_DAMPED(omg::Complex{Float64}, struct_tabLeg::LegendreFHT)
Fill in all the Legendre arrays for a DAMPED mode, i.e., Im[w] < 0.0.
Arguments
omg::Complex{Float64}
: Complex frequency.struct_tabLeg::LegendreFHT
: LegendreFHT structure.
Description
tabLeg!_DAMPED
fills in all the Legendre arrays for a DAMPED mode, where Im[w] < 0.0. It computes the Legendre coefficients Dk(w) and, if needed, Legendre polynomials Pk(w) for the given complex frequency.
Helper Functions
FiniteHilbertTransform.tabuwGLquad
— FunctiontabuwGLquad(K_u::Int64)
Initialize the nodes (u) and weights (w) of the Gauss-Legendre quadrature.
Arguments
K_u::Int64
: Number of Legendre modes.
Returns
tabuGLquad::Vector{Float64}
: Vector of nodes (u) of the Gauss-Legendre quadrature.tabwGLquad::Vector{Float64}
: Vector of weights (w) of the Gauss-Legendre quadrature.
Description
tabuwGLquad
initializes the nodes (u) and weights (w) of the Gauss-Legendre quadrature using the FastGaussQuadrature package. It computes the nodes and weights for the Gauss-Legendre quadrature and returns them as vectors.
Example
K_u = 5
tabuGLquad, tabwGLquad = tabuwGLquad(K_u)
Notes
- This function computes the nodes and weights for the Gauss-Legendre quadrature using external libraries.
FiniteHilbertTransform.tabINVcGLquad
— FunctiontabINVcGLquad(K_u::Int64)
Initialize the normalization constant INVc = 1/c_k = (2k+1)/2.
Arguments
K_u::Int64
: Number of Legendre modes.
Returns
tabINVcGLquad::Vector{Float64}
: Vector of normalization constants.
Description
tabINVcGLquad
initializes the normalization constant INVc for Legendre modes. It computes the normalization constant for each Legendre mode and stores the results in a vector.
Example
K_u = 5
tabINVcGLquad = tabINVcGLquad(K_u)
Notes
- This function computes the normalization constant for Legendre modes indexed from 0 to K_u-1.
- ATTENTION: It corresponds to the INVERSE normalization constant.
FiniteHilbertTransform.tabPGLquad
— FunctiontabPGLquad(K_u::Int64, tabuGLquad::Vector{Float64})
Initialize the values of the Legendre polynomials.
Arguments
K_u::Int64
: Number of Legendre modes.tabuGLquad::Vector{Float64}
: Vector of nodes (u) of the Gauss-Legendre quadrature.
Returns
tabPGLquad::Matrix{Float64}
: Matrix of Legendre polynomials.
Description
tabPGLquad
initializes the values of the Legendre polynomials and stores them in a matrix. Each column of the matrix corresponds to the Legendre polynomials evaluated at a specific node (u) from the Gauss-Legendre quadrature.
Example
K_u = 5
tabuGLquad, _ = tabuwGLquad(K_u)
tabPGLquad = tabPGLquad(K_u, tabuGLquad)
Notes
- This function computes the Legendre polynomials using Bonnet's recurrence relation.
- ATTENTION: It corresponds to the INVERSE normalization constant.
FiniteHilbertTransform.tabGLquad
— FunctiontabGLquad(K_u::Int64)
Initialize the nodes and weights of Gauss-Legendre quadrature.
Arguments
K_u::Int64
: Number of Legendre modes.
Returns
tabuGLquad::Vector{Float64}
: Vector of nodes (u) of the Gauss-Legendre quadrature.tabwGLquad::Vector{Float64}
: Vector of weights (w) of the Gauss-Legendre quadrature.INVcGLquad::Vector{Float64}
: Vector of normalization constants for Legendre polynomials.PGLquad::Matrix{Float64}
: Matrix of Legendre polynomials evaluated at quadrature points.
Description
tabGLquad
initializes the nodes (u) and weights (w) of the Gauss-Legendre quadrature, along with the normalization constants and Legendre polynomials evaluated at quadrature points.
Example
K_u = 5
tabuGLquad, tabwGLquad, INVcGLquad, PGLquad = tabGLquad(K_u)
Chebyshev
FiniteHilbertTransform.ChebyshevFHT
— TypeChebyshevFHT
FiniteHilbertTransform.GettabD!
— Methodfill fht
at a given complex frequency for the integration being considered
integration style selection is automatic: if you want to specify a type, call out to the specific integration method.