desc.integrals.Bounce1D
- class desc.integrals.Bounce1D(grid, data, quad=None, *, automorphism=None, Bref=1.0, Lref=1.0, is_reshaped=False, check=False)Source
Computes bounce integrals using one-dimensional local spline 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 ζ.
Examples
tests/test_integrals.py::TestBounce::test_bounce1d_checks
See also
Bounce2DBounce2Duses 2D pseudo-spectral methods for the same task. The function approximation inBounce1Dis ignorant that the objects to approximate are defined on a bounded subset of ℝ². The domain is projected to ℝ, where information sampled about the function at infinity cannot support reconstruction of the function near the origin. As the functions of interest do not vanish at infinity, pseudo-spectral techniques are not used. Instead, function approximation is done with local splines. This is useful if one can efficiently obtain data along field lines and the number of toroidal transits to follow a field line is not large.
- Parameters:
grid (Grid) – Tensor-product grid in (ρ, α, ζ) Clebsch coordinates. The ζ coordinates (the unique values prior to taking the tensor-product) must be strictly increasing and preferably uniformly spaced. These are used as knots to construct splines. A reference knot density is 100 knots per toroidal transit. Also, the minimum value of the zeta coordinate must be greater than the sentinel value of
-1e5. If this requirement is limiting make a GitHub issue requesting to lower this value.data (dict[str, jnp.ndarray]) – Data evaluated on
grid. Must include names inBounce1D.required_names.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.
Bref (float) – Optional. Reference magnetic field strength for normalization.
Lref (float) – Optional. Reference length scale for normalization.
is_reshaped (bool) – Whether the arrays in
dataare already reshaped to the expected form of shape (…, num ζ) or (…, num α, num ζ) or (num ρ, num α, num ζ). This option can be used to iteratively compute bounce integrals one flux surface or one field line at a time, respectively, reducing memory usage. To do so, set toTrueand provide only those chunks of the reshaped data.check (bool) – Flag for debugging. Must be false for JAX transformations.
Methods
batch(fun, fun_data, desc_data, grid, num_pitch)Compute function
funover phase space in batches.check_points(points, pitch_inv, *[, plot])Check that bounce points are computed correctly.
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, *[, method])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.
points(pitch_inv[, num_well])Compute bounce points.
reshape(grid, f)Reshape arrays for acceptable input to
integrate.Attributes