desc.objectives.EffectiveRipple

class desc.objectives.EffectiveRipple(eq, *, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, loss_function=None, deriv_mode='auto', jac_chunk_size=None, name='Effective ripple', grid=None, X=16, Y=32, Y_B=None, alpha=array([0.]), num_transit=20, num_well=None, num_quad=32, num_pitch=51, pitch_batch_size=None, surf_batch_size=1, spline=False)Source

Proxy for neoclassical transport in the banana regime.

A 3D stellarator magnetic field admits ripple wells that lead to enhanced radial drift of trapped particles. In the banana regime, neoclassical (thermal) transport from ripple wells can become the dominant transport channel. The effective ripple (ε) proxy estimates the neoclassical transport coefficients in the banana regime. To ensure low neoclassical transport, a stellarator is typically optimized so that ε < 0.02.

References

https://doi.org/10.1063/1.873749. Evaluation of 1/ν neoclassical transport in stellarators. V. V. Nemov, S. V. Kasilov, W. Kernbichler, M. F. Heyn. Phys. Plasmas 1 December 1999; 6 (12): 4622–4632.

Notes

Performance will improve significantly by resolving these GitHub issues.
  • 1154 Improve coordinate mapping performance

  • 1294 Nonuniform fast transforms

  • 1303 Patch for differentiable code with dynamic shapes

  • 1206 Upsample data above midplane to full grid assuming stellarator symmetry

  • 1034 Optimizers/objectives with auxiliary output

Parameters:
  • eq (Equilibrium) – Equilibrium to 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 poloidal coordinate. Preferably rounded down to power of 2.

  • Y (int) – Toroidal Chebyshev grid resolution to interpolate the poloidal coordinate. Preferably rounded down to power of 2.

  • Y_B (int) – Desired resolution for algorithm to compute bounce points. Default is double Y. Something like 100 is usually sufficient. Currently, this is the number of knots per toroidal transit over to approximate B with cubic splines.

  • alpha (np.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 None will 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 is Aι+B where A, B are 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 than num_well=(Aι+B)*num_transit is preferable. The check_points or plot methods in desc.integrals.Bounce2D are useful to select a reasonable value.

  • num_quad (int) – Resolution for quadrature of bounce integrals. Default is 32.

  • num_pitch (int) – Resolution for quadrature over velocity coordinate. Default is 51.

  • pitch_batch_size (int) – Number of pitch values with which to compute simultaneously. If given None, then pitch_batch_size is num_pitch. Default is num_pitch.

  • surf_batch_size (int) – Number of flux surfaces with which to compute simultaneously. If given None, then surf_batch_size is grid.num_rho. Default is 1. Only consider increasing if pitch_batch_size is None.

  • spline (bool) – Set to True to replace pseudo-spectral methods with local splines. This can be efficient if num_transit and alpha.size are small, depending on hardware and hardware features used by the JIT compiler. If True, then parameters X and Y are ignored.

  • target ({float, ndarray}, optional) – Target value(s) of the objective. Only used if bounds is None. Must be broadcastable to Objective.dim_f. Defaults to target=0.

  • bounds (tuple of {float, ndarray}, optional) – Lower and upper bounds on the objective. Overrides target. Both bounds must be broadcastable to Objective.dim_f. Defaults to target=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 normalize is True and the target is in physical units, this should also be set to True. Note: Has no effect for this objective.

  • loss_function ({None, 'mean', 'min', 'max'}, 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. auto selects forward or reverse mode based on the size of the input and output of the objective. Has no effect on self.grad or self.hess which always use reverse mode and forward over reverse mode respectively.

    Unless fwd is specified, jac_chunk_size=1 is recommended to reduce memory consumption. In rev mode, reducing the pitch angle parameter pitch_batch_size does not reduce memory consumption, so it is recommended to retain the default for that.

  • name (str, optional) – Name of the objective.

  • jac_chunk_size (int or auto, optional) – Will calculate the Jacobian jac_chunk_size columns at a time, instead of all at once. The memory usage of the Jacobian calculation is roughly memory 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 roughly t = t0+t1/jac_chunk_size so the larger the jac_chunk_size, the faster the calculation takes, at the cost of requiring more memory. If None, it will use the largest size i.e obj.dim_x. Can also help with Hessian computation memory, as Hessian is essentially jacfwd(jacrev(f)), and each of these operations may be chunked. Defaults to chunk_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 the auto chunk_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 the effective ripple.

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

bounds

Lower and upper bounds of the objective.

built

Whether the transforms have been precomputed (or not).

constants

Constant parameters such as transforms and profiles.

dim_f

Number of objective equations.

fixed

Whether the objective fixes individual parameters (or linear combo).

linear

Whether the objective is a linear function (or nonlinear).

name

Name of objective (str).

normalization

normalizing scale factor.

scalar

Whether default "compute" method is a scalar or vector.

target

Target value(s) of the objective.

things

Optimizable things that this objective is tied to.

weight

Weighting to apply to the Objective, relative to other Objectives.