QPGreen

Documentation for QPGreen.

QPGreen.FFTCacheMethod
FFTCache(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 constants
  • T: Floating-point type for allocations. Defaults to Float64.

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.
source
QPGreen.IntegrationCacheType
struct 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
source
QPGreen.IntegrationParametersType
struct 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
source
QPGreen.S_evenMethod
S_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.

source
QPGreen.S_oddMethod
S_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.

source
QPGreen.S₀Method
S₀(β, 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₀.

source
QPGreen.YεMethod
Yε(x, cache::IntegrationCache)

Evaluate the cutoff function at the point x (C^∞ function).

Input arguments

Returns

  • The value of the cutoff function at x.
source
QPGreen.Yε_1st_derMethod
Yε_1st_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 at x.
source
QPGreen.Yε_2nd_derMethod
Yε_2nd_der(x, cache::IntegrationCache)

Evaluate the 2nd order derivative of the cutoff function 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 at x.
source
QPGreen.asin_lsMethod
asin_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.

source
QPGreen.bernoulliMethod
bernoulli(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.

source
QPGreen.check_compatibilityMethod
check_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.

source
QPGreen.eigfunc_expansionMethod
eigfunc_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.
source
QPGreen.eigfunc_expansion_derivativeMethod
eigfunc_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.
source
QPGreen.eval_qp_greenMethod
eval_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 coefficient
    • k: Wavenumber
    • c: Lower cutoff parameter for function χ.
    • c_tilde: Upper cutoff parameter for function χ.
    • epsilon: cutoff parameter for function (recommended: 0.4341).
    • order: Quadrature order.
  • value_interpolator: Bicubic spline interpolator of Ln.

  • Yε_cache: Precomputed cache for cutoff function evaluations.

Keyword Arguments

  • nb_terms: Number of terms to use in the eigenfunction expansion fallback (when x₂ ∉ [-c, c])

Returns

  • G: The approximate value of the quasiperiodic Green's function at point x
source
QPGreen.eval_smooth_qp_greenMethod
eval_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 coefficient
    • k: Wavenumber
    • c: Lower cutoff parameter for function χ.
    • c_tilde: Upper cutoff parameter for function χ.
    • epsilon: cutoff parameter for function (recommended: 0.4341).
    • order: Quadrature order
  • value_interpolator: Bicubic spline interpolator for Ln.

  • Yε_cache: Precomputed cache for cutoff function evaluations.

Keyword Arguments

  • nb_terms: Number of terms to use in the eigenfunction expansion fallback (when x₂ ∉ [-c, c])

Returns

  • G_0: The approximate value of the quasiperiodic Green's function at point x
source
QPGreen.ewaldMethod
ewald(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.

source
QPGreen.f₁Method
f₁(x, cache::IntegrationCache)

Calculate the function f₁.

Input arguments

  • x: point at which the function is evaluated
  • cache: cache for the cut-off functionYε`

Returns

  • The value of the function f₁ at the point x.
source
QPGreen.f₂Method
f₂(x, cache::IntegrationCache)

Calculate the function f₂.

Input arguments

  • x: point at which the function is evaluated
  • cache: cache for the cut-off function

Returns

  • The value of the function f₂.
source
QPGreen.get_F̂ⱼMethod
get_F̂ⱼ(j₁, j₂, c̃, Φ̂₁ⱼ, Φ̂₂ⱼ, F̂₁ⱼ₀::Complex{T}, ::Type{T}) where {T}

Calculate the Fourier coefficients F̂ⱼ.

Input arguments

  • j₁: first index
  • j₂: second index
  • : parameter of the periodicity
  • Φ̂₁ⱼ: Fourier coefficients of the function Φ₁
  • Φ̂₂ⱼ: Fourier coefficients of the function Φ₂
  • F̂₁ⱼ₀: Fourier coefficient of the function F₁ at |j| = 0
  • T: type of the parameter α

Returns

  • The Fourier coefficients F̂ⱼ.
source
QPGreen.get_K̂ⱼ!Method
get_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 coefficients
  • params: Named tuple containing physical and numerical constants
  • N: size of the grid
  • i: index of the grid points
  • fft_cache: cache for the FFT
  • cache: cache for the cut-off function
  • p: plan for the FFT

Returns

  • The Fourier coefficients K̂ⱼ.
source
QPGreen.get_tMethod
get_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.
source
QPGreen.get_ĤⱼMethod
get_Ĥⱼ(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 Ĥⱼ.
source
QPGreen.grad_qp_greenMethod
grad_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 coefficient
    • k: Wavenumber
    • c: Lower cutoff parameter for function χ.
    • c_tilde: Upper cutoff parameter for function χ.
    • epsilon: cutoff parameter for function (recommended: 0.4341).
    • order: Quadrature order
  • grad: Bicubic spline interpolator for the gradient of Ln.

  • Yε_cache: Precomputed cache for cutoff function evaluations.

Keyword Arguments

  • nb_terms: Number of terms to use in the eigenfunction expansion fallback (when x₂ ∉ [-c, c])

Returns

  • ∇G: The approximate value of the fist order derivative of the quasiperiodic Green's function at point x
source
QPGreen.grad_smooth_qp_greenMethod
grad_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 coefficient
    • k: Wavenumber
    • c: Lower cutoff parameter for function χ.
    • c_tilde: Upper cutoff parameter for function χ.
    • epsilon: cutoff parameter for function (recommended: 0.4341).
    • order: Quadrature order
  • grad: Bicubic spline interpolator for the gradient of Ln.

  • Yε_cache: Precomputed cache for cutoff function evaluations.

Keyword Arguments

  • nb_terms: Number of terms to use in the eigenfunction expansion fallback (when x₂ ∉ [-c, c])

Returns

  • ∇G_0: The approximate value of the fist order derivative of the smooth quasiperiodic Green's function at point x
source
QPGreen.h₁Method
h₁(x, params::NamedTuple, cache::IntegrationCache)

Calculate the function h₁.

Input arguments

  • x: point at which the function is evaluated
  • params: Named tuple containing physical and numerical constants

Returns

  • The value of the function h₁.
source
QPGreen.h₁_reducedMethod
h₁_reduced(x, x_norm, α, cache::IntegrationCache)

Calculate the function h₁_reduced (removal of the singularity in the Fourier space for the 1st order derivative).

source
QPGreen.h₂Method
h₂(x, params::NamedTuple, cache)

Calculate the function h₂.

Input arguments

  • x: point at which the function is evaluated
  • params: Named tuple containing physical and numerical constants

Returns

  • The value of the function h₂.
source
QPGreen.h₂_reducedMethod
h₂_reduced(x, x_norm, α, cache::IntegrationCache)

Calculate the function h₂_reduced (removal of the singularity in the Fourier space for the 1st order derivative).

source
QPGreen.image_expansionMethod
image_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.
source
QPGreen.image_expansion_derivativeMethod
image_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.
source
QPGreen.image_expansion_derivative_smoothMethod
image_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.
source
QPGreen.image_expansion_smoothMethod
image_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.
source
QPGreen.init_qp_green_fftMethod
init_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 is 2 * grid_size × 2 * grid_size).

Keyword Arguments

  • derivative: if true, 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.
source
QPGreen.lattice_sums_calculationMethod
lattice_sums_calculation(x, csts, Sₗ; c=0.6, nb_terms=100)

Compute the Green's function using lattice sums. Input arguments:

  • x: evaluation point
  • csts: tuple of constants (β, k, d, M, L)
  • Sₗ: lattice sum coefficients
  • c: cutoff value
  • nb_terms: number of terms in the sum

Keyword arguments:

  • c: cutoff value
  • nb_terms: number of terms in the sum

Returns the value of the Green's function at the point x.

source
QPGreen.lattice_sums_preparationMethod
lattice_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.

source
QPGreen.process_frequency_component!Method
process_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.

source
QPGreen.Φ₁Method
Φ₁(x, cache::IntegrationCache)

Calculate the function Φ₁ (removal of the singularity in the Fourier space).

source
QPGreen.ρₓMethod
ρₓ(x, x_norm, cache::IntegrationCache)

Calculate the function ρₓ (removal of the singularity in the Fourier space for the 1st order derivative).

source
QPGreen.χMethod
χ(x, cache::IntegrationCache)

Evaluate the cutoff function χ at the point x (C^∞ function).

Input arguments

Returns

  • The value of the cutoff function χ at x.
source
QPGreen.χ_derMethod
χ_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 χ at x.
source