desc.objectives.BallooningStability
- class desc.objectives.BallooningStability(eq, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, loss_function=None, deriv_mode='auto', rho=array([0.5]), alpha=None, nturns=3, nzetaperturn=200, zeta0=None, Neigvals=1, lambda0=0.0, w0=1.0, w1=10.0, name='ideal ballooning lambda', jac_chunk_size=None)Source
A type of ideal MHD instability.
Infinite-n ideal MHD ballooning modes are of significant interest. These instabilities are also related to smaller-scale kinetic instabilities. With this class, we optimize MHD equilibria against the ideal ballooning mode.
Targets the following metric:
f = w₀ sum(ReLU(λ-λ₀)) + w₁ max(ReLU(λ-λ₀))
where λ is the negative squared growth rate for each field line (such that λ>0 is unstable), λ₀ is a cutoff, and w₀ and w₁ are weights.
- Parameters:
eq (Equilibrium) –
Equilibriumto be optimized.rho (float) – Flux surface to optimize on. Instabilities often peak near the middle.
alpha (float, ndarray) – Field line labels to optimize. Values should be in [0, 2π). Default is
alpha=0for axisymmetric equilibria, or 8 field lines linearly spaced in [0, π] for non-axisymmetric cases.nturns (int) – Number of toroidal transits of a field line to consider. Field line will run from -π*``nturns`` to π*``nturns``. Default 3.
nzetaperturn (int) – Number of points along the field line per toroidal transit. Total number of points is
nturns*nzetaperturn. Default 100.zeta0 (array-like) – Points of vanishing integrated local shear to scan over. Default 15 points in [-π/2,π/2]. The values
zeta0correspond to values of ι ζ₀ and not ζ₀.Neigvals (int) – Number of top eigenvalues to select. Default is 1.
lambda0 (float) – Threshold for penalizing growth rates in metric above.
w0 (float) – Weights for sum and max terms in metric above.
w1 (float) – Weights for sum and max terms in metric above.
name (str, optional) – Name of the objective function.
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.name – 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 the ballooning stability growth rate.
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.