desc.integrals.Bounce2D
- class desc.integrals.Bounce2D(grid, data, angle, Y_B=None, alpha=Array([0.], dtype=float64), num_transit=20, quad=None, *, automorphism=None, nufft_eps=1e-06, is_reshaped=False, is_fourier=False, Bref=1.0, Lref=1.0, spline=True, check=False, vander=None, **kwargs)Source
Computes bounce integrals using pseudo-spectral methods.
The bounce integral is defined as ∫ f(ρ,α,λ,ℓ) dℓ where
dℓ parametrizes the distance along the field line in meters.
f(ρ,α,λ,ℓ) is the quantity to integrate along the field line.
The boundaries of the integral are bounce points ℓ₁, ℓ₂ s.t. λB(ρ,α,ℓᵢ) = 1.
λ is a constant defining the integral proportional to the magnetic moment over energy.
B is the norm of the magnetic field.
For a particle with fixed λ, bounce points are defined to be the location on the field line such that the particle’s velocity parallel to the magnetic field is zero. The bounce integral is defined up to a sign. We choose the sign that corresponds to the particle’s guiding center trajectory traveling in the direction of increasing field-line-following coordinate ζ.
Refrences
Spectrally accurate, reverse-mode differentiable bounce-averaging algorithm and its applications. Kaya E. Unalmis et al. https://arxiv.org/abs/2412.01724.
Examples
tests/test_integrals.py::TestBounce2D::test_bounce2d_checksdesc/compute/_fast_ion.py::_little_gamma_c_Nemovdesc/compute/_neoclassical.py::_epsilon_32desc/objectives/_fast_ion.py::GammaCdesc/objectives/_neoclassical.py::EffectiveRipple
See also
Bounce1DSome comments comparing
Bounce1DtoBounce2Dare given below.Bounce1Duses lower order accurate, one-dimensional splines.Bounce2Dis superior for optimization objectives in DESC as it solves the moving grid interpolation problem and avoids recomputing 3D Fourier-Zernike series on a time-dependent grid.
- param grid:
Tensor-product grid in (ρ, θ, ζ) with uniformly spaced nodes (θ, ζ) ∈ [0, 2π) × [0, 2π/NFP). Number of poloidal and toroidal nodes preferably rounded down to powers of two. Determines the flux surfaces to compute on and resolution of FFTs. The ζ coordinates (the unique values prior to taking the tensor-product) must be strictly increasing.
- type grid:
Grid
- param data:
Data evaluated on
grid. Must include names inBounce2D.required_names.- type data:
dict[str, jnp.ndarray]
- param angle:
Shape (num ρ, X, Y). Angle returned by
Bounce2D.angle.- type angle:
jnp.ndarray
- param Y_B:
Desired resolution for algorithm to compute bounce points. If the option
splineisTrue, the bounce points are found with 8th order accuracy in this parameter. If the optionsplineisFalse, then the bounce points are found with spectral accuracy in this parameter. A reference value for thesplineoption is 100.An error of ε in a bounce point manifests 𝒪(ε¹ᐧ⁵) error in bounce integrals with (v_∥)¹ and 𝒪(ε⁰ᐧ⁵) error in bounce integrals with (v_∥)⁻¹.
- type Y_B:
int
- param alpha:
Shape (num α, ). Starting field line poloidal labels. Default is single field line. To compute a surface average on a rational surface, it is necessary to average over multiple field lines until the surface is covered sufficiently.
- type alpha:
jnp.ndarray
- param num_transit:
Number of toroidal transits to follow field line. In an axisymmetric device, field line integration over a single poloidal transit is sufficient to capture a surface average. For a 3D configuration, more transits will approximate surface averages on an irrational magnetic surface better, with diminishing returns.
- type num_transit:
int
- param quad:
Quadrature points xₖ and weights wₖ for the approximate evaluation of an integral ∫₋₁¹ g(x) dx = ∑ₖ wₖ g(xₖ). Default is 32 points.
- type quad:
tuple[jnp.ndarray]
- param automorphism:
The first callable should be an automorphism of the real interval [-1, 1]. The second callable should be the derivative of the first. This map defines a change of variable for the bounce integral. The choice made for the automorphism will affect the performance of the quadrature.
- type automorphism:
tuple[Callable] or None
- param nufft_eps:
Precision requested for interpolation with non-uniform fast Fourier transform (NUFFT). If less than
1e-14then NUFFT will not be used.- type nufft_eps:
float
- param is_reshaped:
Whether the arrays in
dataare already reshaped to the expected form of shape (…, num ζ, num θ) or (num ρ, num ζ, num θ). This option can be used to iteratively compute bounce integrals one flux surface at a time, reducing memory usage. To do so, set toTrueand provide only those chunks of the reshaped data. If set toTrue, then it is assumed thatdata["iota"]has shape(grid.num_rho,)or is a scalar.- type is_reshaped:
bool
- param is_fourier:
If true, then it is assumed that
dataholds Fourier transforms as returned byBounce2D.fourieranddata["iota"]has shape(grid.num_rho,)or is a scalar. Default is false.- type is_fourier:
bool
- param Bref:
Optional. Reference magnetic field strength for normalization.
- type Bref:
float
- param Lref:
Optional. Reference length scale for normalization.
- type Lref:
float
- param spline:
Whether to use cubic splines to compute bounce points instead of Chebyshev series. Default is
True.- type spline:
bool
- param check:
Flag for debugging. Must be false for JAX transformations.
- type check:
bool
Methods
angle(eq[, X, Y, rho, iota, params, ...])Return the angle for mapping boundary coordinates to field line coordinates.
batch(fun, fun_data, desc_data, angle, grid, ...)Compute function
funover phase space in batches.check_points(points, pitch_inv, *[, plot])Check that bounce points are computed correctly.
compute_fieldline_length([quad])Compute the (mean) proper length of the field line ∫ dℓ / B.
compute_theta(eq[, X, Y, rho, iota, params, ...])Method has been deprecated in favor of Bounce2D.angle.
fourier(f)Transform to DESC spectral domain.
get_pitch_inv_quad(min_B, max_B, num_pitch)Return 1/λ values and weights for quadrature between
min_Bandmax_B.integrate(integrand, pitch_inv[, data, ...])Bounce integrate ∫ f(ρ,α,λ,ℓ) dℓ.
interp_to_argmin(f, points, *[, nufft_eps, ...])Interpolate
fto the deepest point pⱼ in magnetic well j.plot(l, m[, pitch_inv])Plot B and bounce points on the specified field line.
plot_angle_spectrum(angle, l, *[, truncate, ...])Plot frequency spectrum of the given stream map.
plot_theta(l, m, **kwargs)Plot θ on the specified field line.
points(pitch_inv[, num_well])Compute bounce points.
reshape(grid, f)Reshape arrays for acceptable input to
integrate.Attributes