QPGreen
Documentation for QPGreen.
QPGreen.FFTCacheQPGreen.IntegrationCacheQPGreen.IntegrationParametersQPGreen.S_evenQPGreen.S_oddQPGreen.S₀QPGreen.YεQPGreen.Yε_1st_derQPGreen.Yε_2nd_derQPGreen.asin_lsQPGreen.bernoulliQPGreen.check_compatibilityQPGreen.eigfunc_expansionQPGreen.eigfunc_expansion_gradQPGreen.eigfunc_expansion_hessQPGreen.eval_qp_greenQPGreen.eval_qp_green_asymptoticQPGreen.eval_smooth_qp_greenQPGreen.eval_smooth_qp_green_asymptoticQPGreen.ewaldQPGreen.f_hankelQPGreen.f₁QPGreen.f₂QPGreen.g_singQPGreen.get_F̂ⱼQPGreen.get_K̂ⱼ!QPGreen.get_tQPGreen.get_ĤⱼQPGreen.grad_f_hankelQPGreen.grad_qp_greenQPGreen.grad_qp_green_asymptoticQPGreen.grad_smooth_qp_greenQPGreen.grad_smooth_qp_green_asymptoticQPGreen.h_singQPGreen.hess_f_hankelQPGreen.hess_qp_greenQPGreen.hess_smooth_qp_greenQPGreen.h₁QPGreen.h₁_reducedQPGreen.h₂QPGreen.h₂_reducedQPGreen.image_expansionQPGreen.image_expansion_gradQPGreen.image_expansion_grad_smoothQPGreen.image_expansion_hessQPGreen.image_expansion_hess_smoothQPGreen.image_expansion_smoothQPGreen.init_qp_green_fftQPGreen.init_qp_green_fft_asymptoticQPGreen.lattice_sums_calculationQPGreen.lattice_sums_preparationQPGreen.process_frequency_component!QPGreen.rfftshift_normalization!QPGreen.singularity_hessianQPGreen.Φ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.AbstractIntegrationCacheStructure 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_grad — Methodeigfunc_expansion_grad(z, params::NamedTuple; period=2π, nb_terms=50)Compute the gradient 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 gradient of the α-quasi-periodic Green's function for 2D Helmholtz equation at the point z, computed via the basic eigenfunction expansion.
QPGreen.eigfunc_expansion_hess — Methodeigfunc_expansion_hess(z, params::NamedTuple; period=2π, nb_terms=50)Compute the Hessian 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 hessian 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_qp_green_asymptotic — Methodeval_qp_green_asymptotic(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. Here to remove the singularity we use its asymptotic expansion.
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; 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.
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.eval_smooth_qp_green_asymptotic — Methodeval_smooth_qp_green_asymptotic(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. Here to remove the singularity we use its asymptotic expansion.
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_hankel — Methodf_hankel(x_norm, k, cache::IntegrationCache)Calculate the function f_hankel.
Input arguments
- `x_norm`: norm of the point at which the function is evaluated
- `k`: wavenumber
- `cache`: cache for the cut-off function `Yε`Returns
- The value of the function `f_hankel` at the point `x`.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₁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.g_sing — Methodg_sing(x_norm, k, cache::IntegrationCache)Calculate the function g_sing (removal of the singularity - gradient case - in the Fourier space in the case where you use the Hankel function directly and not its asymptotic form).
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| = 0T: 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_f_hankel — Methodgrad_f_hankel(x_norm, k, cache::IntegrationCache)Calculate the gradient of the function f_hankel.
Input arguments
- `x_norm`: norm of the point at which the function is evaluated
- `k`: wavenumber
- `cache`: cache for the cut-off function `Yε`Returns
- The value of the gradient of `f_hankel` at the point `x`.QPGreen.grad_qp_green — Methodgrad_qp_green(x, params::NamedTuple, grad::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)Compute the gradient 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 gradient of the quasiperiodic Green's function at pointx
QPGreen.grad_qp_green_asymptotic — Methodgrad_qp_green_asymptotic(x, params::NamedTuple, grad::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)Compute the gradient of the α-quasi-periodic Green's function $G(x)$ using the FFT-based method [1] with series expansion fallback. Here to remove the singularity we use its asymptotic expansion.
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 gradient of the quasiperiodic Green's function at pointx
QPGreen.grad_smooth_qp_green — Methodgrad_smooth_qp_green(x, params::NamedTuple, grad::NamedTuple; nb_terms=50)Compute the gradient 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.
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 gradient of the smooth quasiperiodic Green's function at pointx
QPGreen.grad_smooth_qp_green_asymptotic — Methodgrad_smooth_qp_green_asymptotic(x, params::NamedTuple, grad::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)Compute the gradient of the smooth α-quasi-periodic Green's function using the FFT-based method [1] with series expansion fallback. Here to remove the singularity we use its asymptotic expansion.
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 gradient of the smooth quasiperiodic Green's function at pointx
QPGreen.h_sing — Methodh_sing(x_norm, k, cache::IntegrationCache)Calculate the function h_sing (removal of the singularity - hessian case - in the Fourier space in the case where you use the Hankel function directly and not its asymptotic form).
QPGreen.hess_f_hankel — Methodhess_f_hankel(x, r, k, α, cache::IntegrationCache)Calculate the Hessian of the function f_hankel.
Input arguments
- `x`: point at which the function is evaluated
- `r`: norm of the point at which the function is evaluated
- `k`: wavenumber
- `α`: quasi-periodicity parameter
- `cache`: cache for the cut-off function `Yε`Returns
- The value of the Hessian of `f_hankel` at the point `x`.QPGreen.hess_qp_green — Methodhess_qp_green(x, params::NamedTuple, hess::NamedTuple, Yε_cache::IntegrationCache; nb_terms=50)Compute the Hessian 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
hess: Bicubic spline interpolator for the Hessian 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
HG: The approximate value of the Hessian of the quasiperiodic Green's function at pointx
QPGreen.hess_smooth_qp_green — Methodhess_smooth_qp_green(x, params::NamedTuple, hess::NamedTuple; nb_terms=50)Compute the Hessian of the smooth α-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
hess: Bicubic spline interpolator for the Hessian ofLn.
Keyword Arguments
nb_terms: Number of terms to use in the eigenfunction expansion fallback (whenx₂ ∉ [-c, c])
Returns
HG: The approximate value of the Hessian 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_grad — Methodimage_expansion_grad(z, params::NamedTuple; period=2π, nb_terms=50)Compute the gradient 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 gradient 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_grad_smooth — Methodimage_expansion_grad_smooth(z, params::NamedTuple; period=2π, nb_terms=50)Compute the gradient 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 gradient 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_hess — Methodimage_expansion_hess(z, params::NamedTuple; period=2π, nb_terms=50)Compute the Hessian 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 Hessian 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_hess_smooth — Methodimage_expansion_hess_smooth(z, params::NamedTuple; period=2π, nb_terms=50)Compute the Hessian 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 Hessian of the smooth α-quasi-periodic Green's function for the 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; grad=false, hess=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
grad: iftrue, computes additionally the gradient∇Gof the quasi periodic Green's function.hess: iftrue, computes additionally the HessianHGof the quasi periodic Green's function.
Returns
- A NamedTuple with fields
+ `value`: Spline interpolator for the function `Ln`.
+ `grad`: Tuple of spline interpolators for the first derivatives of `Ln` (`∂/∂x₁`, `∂/∂x₂`), if `grad=true`.
+ `hess`: Tuple of spline interpolators for the second derivatives of `Ln` (`∂²/∂x₁²`, `∂²/∂x₁∂x₂`, `∂²/∂x₂²`), if `hess=true`.
+ `cache`: Precomputed integration cache for reuse in later computations.QPGreen.init_qp_green_fft_asymptotic — Methodinit_qp_green_fft_asymptotic(params::NamedTuple, grid_size::Integer; grad=false)Preparation step of the FFT-based algorithm where we remove the singularity by using its asymptotic expansion.
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
grad: iftrue, computes additionally the gradient∇Gof the quasi periodic Green's function.
Returns
- If `grad=false`: a NamedTuple with fields
+ `value`: Spline interpolator for the function `Ln`.
+ `cache`: Precomputed integration cache for reuse in later computations.
- If `grad=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.singularity_hessian — Methodsingularity_hessian(Z, r, k)Calculate the Hessian of the singularity of the QP Green function.
QPGreen.Φ — MethodΦ(x, k, cache::IntegrationCache)Calculate the function Φ (removal of the singularity in the Fourier space in the case where you use the Hankel function directly and not its asymptotic form).
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.