desc.objectives.GammaC
- class desc.objectives.GammaC(eq, *, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, loss_function=None, deriv_mode='auto', jac_chunk_size=None, name='Gamma_c', grid=None, X=16, Y=32, Y_B=None, alpha=Array([0.], dtype=float64), num_transit=20, num_well=None, num_quad=32, num_pitch=65, pitch_batch_size=None, surf_batch_size=1, nufft_eps=1e-07, use_bounce1d=False, Nemov=True, **kwargs)Source
Proxy for fast ion confinement.
A 3D stellarator magnetic field admits ripple wells that lead to enhanced radial drift of trapped particles. The energetic particle confinement metric γ_c quantifies whether the contours of the second adiabatic invariant close on the flux surfaces. In the limit where the poloidal drift velocity majorizes the radial drift velocity, the contours lie parallel to flux surfaces. The optimization metric Γ_c averages γ_c² over the distribution of trapped particles on each flux surface.
The radial electric field has a negligible effect, since fast particles have high energy with collisionless orbits, so it is assumed to be zero.
References
- [1] Poloidal motion of trapped particle orbits in real-space coordinates.
V. V. Nemov, S. V. Kasilov, W. Kernbichler, G. O. Leitold. Phys. Plasmas 1 May 2008; 15 (5): 052501. https://doi.org/10.1063/1.2912456. Equation 61.
- [2] A model for the fast evaluation of prompt losses of energetic ions in
stellarators. Equation 16. J.L. Velasco et al. 2021 Nucl. Fusion 61 116059. https://doi.org/10.1088/1741-4326/ac2994.
- [3] Spectrally accurate, reverse-mode differentiable bounce-averaging
algorithm and its applications. Kaya E. Unalmis et al. https://arxiv.org/abs/2412.01724.
Notes
- Performance will improve significantly by resolving these GitHub issues.
1206Upsample data above midplane to full grid assuming stellarator symmetry1034Optimizers/objectives with auxiliary output
- Parameters:
eq (Equilibrium) –
Equilibriumto be optimized.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. Default grid samples the boundary surface at ρ=1.
X (int) – Poloidal Fourier grid resolution to interpolate the angle. Preferably rounded down to power of 2.
Y (int) – Toroidal Chebyshev grid resolution over a single field period to interpolate the angle. Preferably rounded down to power of 2.
Y_B (int) –
Desired resolution for algorithm to compute bounce points. The bounce points are found with 8th order accuracy in this parameter. A reference value is 100.
An error of ε in a bounce point manifests 𝒪(ε¹ᐧ⁵) error in bounce integrals with (v_∥)¹ and 𝒪(ε⁰ᐧ⁵) error in bounce integrals with (v_∥)⁻¹.
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.
num_well (int) –
Maximum number of wells to detect for each pitch and field line. Giving
-1will detect all wells but due to current limitations in JAX this will have worse performance. Specifying a number that tightly upper bounds the number of wells will increase performance. In general, an upper bound on the number of wells per toroidal transit isAι+CwhereA,Care the poloidal and toroidal Fourier resolution of B, respectively, in straight-field line PEST coordinates, and ι is the rotational transform normalized by 2π. A tighter upper bound thannum_well=(Aι+C)*num_transitis preferable. Thecheck_pointsorplotmethods indesc.integrals.Bounce2Dare useful to select a reasonable value.This is the most important parameter to specify for performance.
num_quad (int) – Resolution for quadrature of bounce integrals. Default is 32.
num_pitch (int) – Resolution for quadrature over velocity coordinate. Default is 65.
pitch_batch_size (int) – Number of pitch values with which to compute simultaneously. If given
None, thenpitch_batch_sizeisnum_pitch. Default isnum_pitch.surf_batch_size (int) – Number of flux surfaces with which to compute simultaneously. If given
None, thensurf_batch_sizeisgrid.num_rho. Default is1. Only consider increasing ifpitch_batch_sizeisNone.nufft_eps (float) – Precision requested for interpolation with non-uniform fast Fourier transform (NUFFT). If less than
1e-14then NUFFT will not be used.use_bounce1d (bool) – Set to
Trueto useBounce1Dinstead ofBounce2D, basically replacing some pseudo-spectral methods with local splines. This can be efficient ifnum_transitandalpha.sizeare small, depending on hardware and hardware features used by the JIT compiler. IfTrue, then parametersX,Y,nufft_epsare ignored.Nemov (bool) –
Whether to use the Γ_c as defined by Nemov et al. or Velasco et al. Default is Nemov. Set to
Falseto use Velascos’s.Nemov’s Γ_c converges to a finite nonzero value in the infinity limit of the number of toroidal transits. Velasco’s expression has a secular term that drives the result to zero as the number of toroidal transits increases if the secular term is not averaged out from the singular integrals. At finite resolution, an optimization using Velasco’s metric may need to be evaluated by measuring decrease in Γ_c at a fixed number of toroidal transits.
target ({float, ndarray}, optional) – Target value(s) of the objective. Only used if
boundsisNone. Must be broadcastable toObjective.dim_f. Defaults totarget=0.bounds (tuple of {float, ndarray}, optional) – Lower and upper bounds on the objective. Overrides
target. Both bounds must be broadcastable toObjective.dim_f. Defaults totarget=0.weight ({float, ndarray}, optional) – Weighting to apply to the Objective, relative to other Objectives. Must be broadcastable to
Objective.dim_f.normalize (bool, optional) – Whether to compute the error in physical units or non-dimensionalize. Note: Has no effect for this objective.
normalize_target (bool, optional) – Whether target and bounds should be normalized before comparing to computed values. If
normalizeisTrueand the target is in physical units, this should also be set toTrue. Note: Has no effect for this objective.loss_function ({None, 'mean', 'min', 'max','sum'}, optional) – Loss function to apply to the objective values once computed. This loss function is called on the raw compute value, before any shifting, scaling, or normalization.
deriv_mode ({"auto", "fwd", "rev"}) –
Specify how to compute Jacobian matrix, either forward mode or reverse mode AD.
autoselects forward or reverse mode based on the size of the input and output of the objective. Has no effect onself.gradorself.hesswhich always use reverse mode and forward over reverse mode respectively.The default mode of
autowill likely chooserevfor this objective. Inrevmode, reducing the value of the parameterpitch_batch_sizedoes not reduce memory consumption, so it is recommended to retain the default unless you have explicitly requested to usefwdmode.name (str, optional) – Name of the objective.
jac_chunk_size (int or
auto, optional) – Will calculate the Jacobianjac_chunk_sizecolumns at a time, instead of all at once. The memory usage of the Jacobian calculation is roughlymemory usage = m0+m1*jac_chunk_size: the smaller the chunk size, the less memory the Jacobian calculation will require (with some baseline memory usage). The time it takes to compute the Jacobian is roughlyt = t0+t1/jac_chunk_sizeso the larger thejac_chunk_size, the faster the calculation takes, at the cost of requiring more memory. IfNone, it will use the largest size i.eobj.dim_x. Can also help with Hessian computation memory, as Hessian is essentiallyjacfwd(jacrev(f)), and each of these operations may be chunked. Defaults tochunk_size=None. Note: When running on a CPU (not a GPU) on a HPC cluster, DESC is unable to accurately estimate the available device memory, so theautochunk_size option will yield a larger chunk size than may be needed. It is recommended to manually choose a chunk_size if an OOM error is experienced in this case.
Methods
build([use_jit, verbose])Build constant arrays.
compute(params[, constants])Compute Γ_c.
compute_scalar(*args, **kwargs)Compute the scalar form of the objective.
compute_scaled(*args, **kwargs)Compute and apply weighting and normalization.
compute_scaled_error(*args, **kwargs)Compute and apply the target/bounds, weighting, and normalization.
compute_unscaled(*args, **kwargs)Compute the raw value of the objective.
copy([deepcopy])Return a (deep)copy of this object.
equiv(other)Compare equivalence between DESC objects.
grad(*args, **kwargs)Compute gradient vector of self.compute_scalar wrt x.
hess(*args, **kwargs)Compute Hessian matrix of self.compute_scalar wrt x.
jac_scaled(*args, **kwargs)Compute Jacobian matrix of self.compute_scaled wrt x.
jac_scaled_error(*args, **kwargs)Compute Jacobian matrix of self.compute_scaled_error wrt x.
jac_unscaled(*args, **kwargs)Compute Jacobian matrix of self.compute_unscaled wrt x.
jvp_scaled(v, x[, constants])Compute Jacobian-vector product of self.compute_scaled.
jvp_scaled_error(v, x[, constants])Compute Jacobian-vector product of self.compute_scaled_error.
jvp_unscaled(v, x[, constants])Compute Jacobian-vector product of self.compute_unscaled.
load(load_from[, file_format])Initialize from file.
print_value(args[, args0])Print the value of the objective and return a dict of values.
save(file_name[, file_format, file_mode])Save the object.
xs(*things)Return a tuple of args required by this objective from optimizable things.
Attributes
Lower and upper bounds of the objective.
Whether the transforms have been precomputed (or not).
Constant parameters such as transforms and profiles.
Number of objective equations.
Whether the objective fixes individual parameters (or linear combo).
Whether the objective is a linear function (or nonlinear).
Name of objective (str).
normalizing scale factor.
Whether default "compute" method is a scalar or vector.
Target value(s) of the objective.
Optimizable things that this objective is tied to.
Weighting to apply to the Objective, relative to other Objectives.