QPGreen
Documentation for QPGreen.
QPGreen.FFTCache
QPGreen.IntegrationCache
QPGreen.IntegrationParameters
QPGreen.S_even
QPGreen.S_odd
QPGreen.S₀
QPGreen.Yε
QPGreen.Yε_1st_der
QPGreen.Yε_2nd_der
QPGreen.asin_ls
QPGreen.bernoulli
QPGreen.check_compatibility
QPGreen.eigfunc_expansion
QPGreen.eigfunc_expansion_derivative
QPGreen.eval_qp_green
QPGreen.eval_smooth_qp_green
QPGreen.ewald
QPGreen.f₁
QPGreen.f₂
QPGreen.get_F̂ⱼ
QPGreen.get_K̂ⱼ!
QPGreen.get_t
QPGreen.get_Ĥⱼ
QPGreen.grad_qp_green
QPGreen.grad_smooth_qp_green
QPGreen.h₁
QPGreen.h₁_reduced
QPGreen.h₂
QPGreen.h₂_reduced
QPGreen.image_expansion
QPGreen.image_expansion_derivative
QPGreen.image_expansion_derivative_smooth
QPGreen.image_expansion_smooth
QPGreen.init_qp_green_fft
QPGreen.lattice_sums_calculation
QPGreen.lattice_sums_preparation
QPGreen.process_frequency_component!
QPGreen.rfftshift_normalization!
QPGreen.Φ₁
QPGreen.ρₓ
QPGreen.χ
QPGreen.χ_der
QPGreen.FFTCache
— MethodFFTCache(N::Integers, grid_size::Integer, params::NamedTuple, T=Float64)
Construct a containers for FFT operations.
Arguments
N
: total number of grid points in one dimension.grid_size
: Number of half of the grid points in one dimension.params
: Physical and numerical constantsT
: Floating-point type for allocations. Defaults toFloat64
.
Returns
An FFTCache
object containing:
j_idx
: Vector of integers for the computations of Fourier coefficients[-grid_size, grid_size-1]
.t_j_fft
: Spatial grid points in[-c̃, c̃]
.Preallocated complex vectors for FFT operations:
eval_int_fft_1D
: 1D integration using FFT.shift_sample_eval_int
: Shifted samples for FFT.fft_eval
: FFT evaluation.shift_fft_1d
: Shifted FFT result.fft_eval_flipped
: Transposed result of FFT.
QPGreen.IntegrationCache
— Typestruct IntegrationCache{T1<:Real, T2<:Signed} <: QPGreen.AbstractIntegrationCache
Structure storing the normalization factor and the parameters of integration for the cutoff functions.
normalization::Real
: Normalization factor
params::QPGreen.IntegrationParameters
: Parameters of integration
QPGreen.IntegrationParameters
— Typestruct IntegrationParameters{T1<:Real, T2<:Signed}
Structure storing the parameters of integration for the cutoff functions.
a::Real
: Lower bound
b::Real
: Upper bound
order::Signed
: Order of the cutoff function
QPGreen.S_even
— MethodS_even(l, β, k, d, M)
Compute the lattice sum Sₗ for even values of l
. Input arguments:
l
: integer value.β
: parameter β.k
: parameter k.d
: parameter d.M
: number of terms in the sum.
Returns the value of the lattice sum Sₗ for even values of l
.
QPGreen.S_odd
— MethodS_odd(l, β, k, d, M)
Compute the lattice sum Sₗ for odd values of l
.
Input arguments:
l
: integer value.β
: parameter β.k
: parameter k.d
: parameter d.M
: number of terms in the sum.
Returns the value of the lattice sum Sₗ for odd values of l
.
QPGreen.S₀
— MethodS₀(β, k, d, M)
Compute the lattice sum S₀.
Input arguments:
β
: parameter β.k
: parameter k.d
: parameter d.M
: number of terms in the sum.
Returns the value of the lattice sum S₀.
QPGreen.Yε
— MethodYε(x, cache::IntegrationCache)
Evaluate the cutoff function Yε
at the point x
(C^∞
function).
Input arguments
- x: point at which the cutoff function is evaluated.
- cache: see
IntegrationCache
.
Returns
- The value of the cutoff function
Yε
atx
.
QPGreen.Yε_1st_der
— MethodYε_1st_der(x, cache::IntegrationCache)
Evaluate the derivative of the cutoff function Yε
at the point x
.
Input arguments
- x: point at which the derivative of the cutoff function is evaluated.
- cache: see
IntegrationCache
.
Returns
- The value of the derivative of the cutoff function
Yε
atx
.
QPGreen.Yε_2nd_der
— MethodYε_2nd_der(x, cache::IntegrationCache)
Evaluate the 2nd order derivative of the cutoff function Yε
at the point x
.
Input arguments
- x: point at which the 2nd order derivative of the cutoff function is evaluated.
- cache: see
IntegrationCache
.
Returns
- The value of the 2nd order derivative of the cutoff function
Yε
atx
.
QPGreen.asin_ls
— Methodasin_ls(x)
Compute the inverse sine function (arcsine function) for values outside of the interval [-1,1]. Input arguments:
x
: real values outside the interval [-1,1].
Returns the inverse sine function (arcsine function) evaluated at x
. For values inside the interval [-1,1], use the standard asin
function.
QPGreen.bernoulli
— Methodbernoulli(n, x)
Compute the Bernoulli polynomial of degree n
evaluated at x
. Input arguments:
n
: degree of the Bernoulli polynomial.x
: value at which the polynomial is evaluated.
Returns the nth Bernoulli polynomial evaluated at x
.
QPGreen.check_compatibility
— Methodcheck_compatibility(alpha, k)
Check if the parameters alpha
and k
are compatible, i.e., if βₙ
is different from zero. For βₙ
equal to zero, the algorithm fails.
QPGreen.eigfunc_expansion
— Methodeigfunc_expansion(z, params::NamedTuple; period=2π, nb_terms=50)
Compute the α-quasi-periodic Green's function using the eigenfunction expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the α-quasi-periodic Green's function for 2D Helmholtz equation at the point z, computed by the basic eigenfunction expansion.
QPGreen.eigfunc_expansion_derivative
— Methodeigfunc_expansion_derivative(z, params::NamedTuple; period=2π, nb_terms=50)
Compute the derivative of the α-quasi-periodic Green's function using the eigenfunction expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the derivative of the α-quasi-periodic Green's function for 2D Helmholtz equation at the point z, computed via the basic eigenfunction expansion.
QPGreen.eval_qp_green
— Methodeval_qp_green(x, params::NamedTuple, interpolator, Yε_cache::IntegrationCache; nb_terms=50)
Compute the quasiperiodic Green's function $G(x)$ using the FFT-based method [1] with series expansion fallback.
Input arguments
x
: 2D point at which to evaluate the Green's function.params
: Physical and numerical parameters containing:alpha
: Quasiperiodicity coefficientk
: Wavenumberc
: Lower cutoff parameter for functionχ
.c_tilde
: Upper cutoff parameter for functionχ
.epsilon
: cutoff parameter for functionYε
(recommended:0.4341
).order
: Quadrature order.
value_interpolator
: Bicubic spline interpolator ofLn
.Yε_cache
: Precomputed cache for cutoff functionYε
evaluations.
Keyword Arguments
nb_terms
: Number of terms to use in the eigenfunction expansion fallback (whenx₂ ∉ [-c, c]
)
Returns
G
: The approximate value of the quasiperiodic Green's function at pointx
QPGreen.eval_smooth_qp_green
— Methodeval_smooth_qp_green(x, params::NamedTuple, value_interpolator, Yε_cache::IntegrationCache; nb_terms=50)
Compute the smooth α-quasi-periodic Green's function $G_0(x)$ (i.e. without the term $H_0^{(1)(k|x|)}$ using the FFT-based method [1] with series expansion fallback.
Input arguments
x
: 2D point at which to evaluate the Green's function.params
: Physical and numerical parameters containing:alpha
: Quasiperiodicity coefficientk
: Wavenumberc
: Lower cutoff parameter for functionχ
.c_tilde
: Upper cutoff parameter for functionχ
.epsilon
: cutoff parameter for functionYε
(recommended:0.4341
).order
: Quadrature order
value_interpolator
: Bicubic spline interpolator forLn
.Yε_cache
: Precomputed cache for cutoff functionYε
evaluations.
Keyword Arguments
nb_terms
: Number of terms to use in the eigenfunction expansion fallback (whenx₂ ∉ [-c, c]
)
Returns
G_0
: The approximate value of the quasiperiodic Green's function at pointx
QPGreen.ewald
— Methodewald(z, csts)
Evaluation of the Green's function using Ewald's method. Input arguments:
z
: coordinates.csts
: tuple of constants(a, M₁, M₂, N, β, k, d)
.
Returns the value of the Green's function at the point z
.
QPGreen.f₁
— Methodf₁(x, cache::IntegrationCache)
Calculate the function f₁
.
Input arguments
x
: point at which the function is evaluatedcache
: cache for the cut-off function
Yε`
Returns
- The value of the function
f₁
at the pointx
.
QPGreen.f₂
— Methodf₂(x, cache::IntegrationCache)
Calculate the function f₂
.
Input arguments
x
: point at which the function is evaluatedcache
: cache for the cut-off functionYε
Returns
- The value of the function
f₂
.
QPGreen.get_F̂ⱼ
— Methodget_F̂ⱼ(j₁, j₂, c̃, Φ̂₁ⱼ, Φ̂₂ⱼ, F̂₁ⱼ₀::Complex{T}, ::Type{T}) where {T}
Calculate the Fourier coefficients F̂ⱼ
.
Input arguments
j₁
: first indexj₂
: second indexc̃
: parameter of the periodicityΦ̂₁ⱼ
: Fourier coefficients of the functionΦ₁
Φ̂₂ⱼ
: Fourier coefficients of the functionΦ₂
F̂₁ⱼ₀
: Fourier coefficient of the functionF₁
at|j| = 0
T
: type of the parameterα
Returns
- The Fourier coefficients
F̂ⱼ
.
QPGreen.get_K̂ⱼ!
— Methodget_K̂ⱼ!(K̂ⱼ, params::NamedTuple, N, i, fft_cache::FFTCache{T}, cache::IntegrationCache, p) where {T}
Mutating function that computes the Fourier coefficients K̂ⱼ
.
Input arguments
K̂ⱼ
: matrix to store the Fourier coefficientsparams
: Named tuple containing physical and numerical constantsN
: size of the gridi
: index of the grid pointsfft_cache
: cache for the FFTcache
: cache for the cut-off functionYε
p
: plan for the FFT
Returns
- The Fourier coefficients
K̂ⱼ
.
QPGreen.get_t
— Methodget_t(x)
Find the value of t
in the interval [-π, π[
such that x = 2nπ + t
.
Input arguments
x
: point at which the function is evaluated
Returns
- The value of
t
.
QPGreen.get_Ĥⱼ
— Methodget_Ĥⱼ(j₁, j₂, params::NamedTuple, F̂₁ⱼ, F̂₂ⱼ, ρ̂ₓⱼ, Φ̂₃ⱼ, ĥ₁ⱼ, ĥ₂ⱼ, Ĥ₁ⱼ₀, T)
Calculate the Fourier coefficients Ĥⱼ
.
Input arguments
- `j₁`: first index
- `j₂`: second index
- `params`: named tuple containing physical and numerical constants
- `F̂₁ⱼ`: Fourier coefficients of the function `F₁`
- `F̂₂ⱼ`: Fourier coefficients of the function `F₂`
- `ρ̂ₓⱼ`: Fourier coefficients of the function `ρₓ`
- `Φ̂₃ⱼ`: Fourier coefficients of the function `Φ₃`
- `ĥ₁ⱼ`: Fourier coefficients of the function `h₁_reduced`
- `ĥ₂ⱼ`: Fourier coefficients of the function `h₂_reduced`
- `Ĥ₁ⱼ₀`: Fourier coefficient of the function `H₁` at `|j| = 0`
- `T`: type of the parameter `α`
Returns
- The Fourier coefficients
Ĥⱼ
.
QPGreen.grad_qp_green
— Methodgrad_qp_green(x, params::NamedTuple, grad::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)
Compute the first order derivative of the α-quasi-periodic Green's function $G(x)$ using the FFT-based method [1] with series expansion fallback.
Input arguments
x
: 2D point at which to evaluate the Green's function.params
: Physical and numerical parameters containing:alpha
: Quasiperiodicity coefficientk
: Wavenumberc
: Lower cutoff parameter for functionχ
.c_tilde
: Upper cutoff parameter for functionχ
.epsilon
: cutoff parameter for functionYε
(recommended:0.4341
).order
: Quadrature order
grad
: Bicubic spline interpolator for the gradient ofLn
.Yε_cache
: Precomputed cache for cutoff functionYε
evaluations.
Keyword Arguments
nb_terms
: Number of terms to use in the eigenfunction expansion fallback (whenx₂ ∉ [-c, c]
)
Returns
∇G
: The approximate value of the fist order derivative of the quasiperiodic Green's function at pointx
QPGreen.grad_smooth_qp_green
— Methodgrad_smooth_qp_green(x, params::NamedTuple, grad::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)
Compute the first order derivative of the smooth α-quasi-periodic Green's function using the FFT-based method [1] with series expansion fallback.
Input arguments
x
: 2D point at which to evaluate the Green's function.params
: Physical and numerical parameters containing:alpha
: Quasiperiodicity coefficientk
: Wavenumberc
: Lower cutoff parameter for functionχ
.c_tilde
: Upper cutoff parameter for functionχ
.epsilon
: cutoff parameter for functionYε
(recommended:0.4341
).order
: Quadrature order
grad
: Bicubic spline interpolator for the gradient ofLn
.Yε_cache
: Precomputed cache for cutoff functionYε
evaluations.
Keyword Arguments
nb_terms
: Number of terms to use in the eigenfunction expansion fallback (whenx₂ ∉ [-c, c]
)
Returns
∇G_0
: The approximate value of the fist order derivative of the smooth quasiperiodic Green's function at pointx
QPGreen.h₁
— Methodh₁(x, params::NamedTuple, cache::IntegrationCache)
Calculate the function h₁
.
Input arguments
x
: point at which the function is evaluatedparams
: Named tuple containing physical and numerical constants
Returns
- The value of the function
h₁
.
QPGreen.h₁_reduced
— Methodh₁_reduced(x, x_norm, α, cache::IntegrationCache)
Calculate the function h₁_reduced (removal of the singularity in the Fourier space for the 1st order derivative).
QPGreen.h₂
— Methodh₂(x, params::NamedTuple, cache)
Calculate the function h₂
.
Input arguments
x
: point at which the function is evaluatedparams
: Named tuple containing physical and numerical constants
Returns
- The value of the function
h₂
.
QPGreen.h₂_reduced
— Methodh₂_reduced(x, x_norm, α, cache::IntegrationCache)
Calculate the function h₂_reduced (removal of the singularity in the Fourier space for the 1st order derivative).
QPGreen.image_expansion
— Methodimage_expansion(z, params::NamedTuple; period=2π, nb_terms=50)
Compute the α-quasi-periodic Green's function using the image expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the α-quasi-periodic Green's function for the 2D Helmholtz equation at the point z, computed via the basic image expansion.
QPGreen.image_expansion_derivative
— Methodimage_expansion_derivative(z, params::NamedTuple; period=2π, nb_terms=50)
Compute the derivative of the α-quasi-periodic Green's function using the image expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the derivative of the α-quasi-periodic Green's function for the 2D Helmholtz equation at the point z, computed via the basic image expansion.
QPGreen.image_expansion_derivative_smooth
— Methodimage_expansion_derivative_smooth(z, params:NamedTuple; period=2π, nb_terms=50)
Compute the derivative of the smooth α-quasi-periodic Green's function using the image expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the derivative of the smooth α-quasi-periodic Green's function for 2D Helmholtz equation at the point z, computed via the basic image expansion.
QPGreen.image_expansion_smooth
— Methodimage_expansion_smooth(z, params::NamedTuple; period=2π, nb_terms=50)
Compute the smooth α-quasi-periodic Green's function (i.e. without the term $H_0^{(1)(k|x|)}$ using the image expansion.
Input arguments
- z: coordinates of the difference between the target point and source point.
- params: named tuple of the physical and numerical parameters for the problem definition.
Returns
- The value of the smooth α-quasi-periodic Green's function for the 2D Helmholtz equation at the point z, computed via the basic image expansion.
QPGreen.init_qp_green_fft
— Methodinit_qp_green_fft(params::NamedTuple, grid_size::Integer; derivative=false)
Preparation step of the FFT-based algorithm.
Input arguments
params
: Physical and numerical parameters, containing:alpha
: Quasiperiodicity coefficient.k
: Wave number.c
: Lower cutoff parameter for function χ.c_tilde
: Upper cutoff parameter for function χ.epsilon
: cutoff parameter for function Yε (recommended:0.4341
).order
: Quadrature order for integration.
grid_size
: Number of grid points per dimension (grid is2 * grid_size × 2 * grid_size
).
Keyword Arguments
derivative
: iftrue
, computes additionally the first-order derivative (∇G
) of the quasi periodic Green's function.
Returns
- If `derivative=false`: a NamedTuple with fields
+ `value`: Spline interpolator for the function `Ln`.
+ `cache`: Precomputed integration cache for reuse in later computations.
- If `derivative=true`: a NamedTuple with fields
+ `value`: Spline interpolator for the function `Ln`.
+ `grad`: Tuple of spline interpolators for the first derivatives of `Ln` (`∂/∂x₁`, `∂/∂x₂`).
+ `cache`: Precomputed integration cache for reuse in later computations.
QPGreen.lattice_sums_calculation
— Methodlattice_sums_calculation(x, csts, Sₗ; c=0.6, nb_terms=100)
Compute the Green's function using lattice sums. Input arguments:
x
: evaluation pointcsts
: tuple of constants(β, k, d, M, L)
Sₗ
: lattice sum coefficientsc
: cutoff valuenb_terms
: number of terms in the sum
Keyword arguments:
c
: cutoff valuenb_terms
: number of terms in the sum
Returns the value of the Green's function at the point x
.
QPGreen.lattice_sums_preparation
— Methodlattice_sums_preparation(csts)
Perform the preparation step for the lattice sums algorithm. Input arguments:
csts
: tuple of constants(β, k, d, M, L)
representing quasi-periodicty parameter, wave number , period, number of terms, and number of terms in the sum respectively.
Returns the lattice sum coefficients.
QPGreen.process_frequency_component!
— Methodprocess_frequency_component!(i, N, params, fft_cache, χ_cache, fft_plan, K̂ⱼ)
Internal helper for processing a single frequency component during quasi-periodic Green's function computation. Handles coefficient calculations and frequency index mapping for position i
in the FFT grid.
QPGreen.rfftshift_normalization!
— Methodrfftshift_normalization!(Φ̂₁ⱼ, fft_Φ₁_eval, N, c̃)
Shift the Fourier coefficients obtained by rfft and normalize them.
QPGreen.Φ₁
— MethodΦ₁(x, cache::IntegrationCache)
Calculate the function Φ₁ (removal of the singularity in the Fourier space).
QPGreen.ρₓ
— Methodρₓ(x, x_norm, cache::IntegrationCache)
Calculate the function ρₓ (removal of the singularity in the Fourier space for the 1st order derivative).
QPGreen.χ
— Methodχ(x, cache::IntegrationCache)
Evaluate the cutoff function χ
at the point x
(C^∞
function).
Input arguments
x
: point at which the cutoff function is evaluated.cache
: seeIntegrationCache
.
Returns
- The value of the cutoff function
χ
atx
.
QPGreen.χ_der
— Methodχ_der(x, cache::IntegrationCache)
Evaluate the derivative of the cutoff function χ
at the point x
.
Input arguments
- x: point at which the derivative of the cutoff function is evaluated.
- cache: see
IntegrationCache
.
Returns
- The value of the derivative of the cutoff function
χ
atx
.