desc.integrals.Bounce2D
- class desc.integrals.Bounce2D(grid, data, theta, 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)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 ζ.
Notes
- Magnetic field line with label α, defined by B = ∇ψ × ∇α, is determined from
α : ρ, θ, ζ ↦ θ + Λ(ρ,θ,ζ) − ι(ρ) [ζ + ω(ρ,θ,ζ)]
- Interpolate Fourier-Chebyshev series to DESC poloidal coordinate.
θ : ρ, α, ζ ↦ tₘₙ(ρ) exp(jmα) Tₙ(ζ)
- Compute bounce points.
λ B(ζₖ) = 1
- Interpolate smooth periodic parts of integrand with FFTs.
G : ρ, α, ζ ↦ gₘₙ(ρ) exp(j [m θ(ρ,α,ζ) + n ζ])
- Perform Gaussian quadrature after removing singularities.
Fᵢ : ρ, α, λ, ζ₁, ζ₂ ↦ ∫ᵢ f(ρ,α,λ,ζ,{Gⱼ}) dζ
If the map G is multivalued at a physical location, then it is still permissible if separable into periodic and secular parts. In that case, supply the periodic part, which will be interpolated with FFTs, and use the provided coordinates θ,ζ ∈ ℝ to compose G.
Examples
See
tests/test_integrals.py::TestBounce2D::test_bounce2d_checks.See also
Bounce1DBounce1Duses one-dimensional splines for the same task.Bounce2Dsolves the dominant cost of optimization objectives in DESC relying onBounce1D: Computing a dense optimization-step dependent grid along field lines and interpolating 3D FourierZernike series to this grid. The function approximation done here requires FourierZernike series on a smaller fixed grid and uses FFTs to compute the map between coordinate systems. 2D interpolation enables tracing the field line for more toroidal transits. Performance will improve significantly by resolving GitHub issue1303: Patch for differentiable code with dynamic shapes.
- Parameters:
grid (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.
data (dict[str, jnp.ndarray]) – Data evaluated on
grid. Must include names inBounce2D.required_names.theta (jnp.ndarray) – Shape (num rho, X, Y). DESC coordinates θ from
Bounce2D.compute_theta.XandYare preferably rounded down to powers of two.Y_B (int) – Desired resolution for algorithm to compute bounce points. A reference value is 100. Default is double
Y.alpha (jnp.ndarray) – Shape (num alpha, ). 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.
num_transit (int) – 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.
quad (tuple[jnp.ndarray]) – Quadrature points xₖ and weights wₖ for the approximate evaluation of an integral ∫₋₁¹ g(x) dx = ∑ₖ wₖ g(xₖ). Default is 32 points.
automorphism (tuple[Callable] or None) – 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.
nufft_eps (float) – Precision requested for interpolation with non-uniform fast Fourier transform (NUFFT). If less than
1e-14then NUFFT will not be used.is_reshaped (bool) – Whether the arrays in
dataare already reshaped to the expected form of shape (…, num zeta, num theta) or (num rho, num zeta, num theta). 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.is_fourier (bool) – 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.Bref (float) – Optional. Reference magnetic field strength for normalization.
Lref (float) – Optional. Reference length scale for normalization.
spline (bool) – Whether to use cubic splines to compute bounce points instead of Chebyshev series. Note the algorithm for efficient root-finding on Chebyshev series has not been implemented.
check (bool) – Flag for debugging. Must be false for JAX transformations.
Methods
check_points(points, pitch_inv, *[, plot])Check that bounce points are computed correctly.
compute_fieldline_length([quad, vander])Compute the (mean) proper length of the field line ∫ dℓ / B.
compute_theta(eq[, X, Y, rho, iota, params, ...])Return DESC coordinates θ of (α,ζ) Fourier Chebyshev basis nodes.
copy([deepcopy])Return a (deep)copy of this object.
equiv(other)Compare equivalence between DESC objects.
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.load(load_from[, file_format])Initialize from file.
plot(l, m[, pitch_inv])Plot B and bounce points on the specified field line.
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.save(file_name[, file_format, file_mode])Save the object.
Attributes