desc.objectives.ToroidalFlux
- class desc.objectives.ToroidalFlux(eq, field, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, loss_function=None, deriv_mode='auto', field_grid=None, eval_grid=None, name='toroidal-flux', field_fixed=False, eq_fixed=False, jac_chunk_size=None, *, bs_chunk_size=None, **kwargs)Source
Target the toroidal flux in an equilibrium from a magnetic field.
This objective is needed when performing stage-two coil optimization on a vacuum equilibrium, to avoid the trivial solution of minimizing Bn by making the coil currents zero. Instead, this objective ensures the coils create the necessary toroidal flux for the equilibrium field.
Will try to use the vector potential method to calculate the toroidal flux (Φ = ∮ 𝐀 ⋅ 𝐝𝐥 over the perimeter of a constant zeta plane) instead of the brute force method using the magnetic field (Φ = ∯ 𝐁 ⋅ 𝐝𝐒 over a constant zeta XS). The vector potential method is much more efficient, however not every
MagneticFieldobject has a vector potential available to compute, so in those cases the magnetic field method is used.- Parameters:
eq (Equilibrium or FourierRZToroidalSurface) – Equilibrium (or QFM surface) for which the toroidal flux will be calculated.
field (MagneticField) – MagneticField object, the parameters of this will be optimized to minimize the objective.
field_grid (Grid, optional) – Grid containing the nodes to evaluate field source at on the winding surface. (used if e.g. field is a CoilSet or FourierCurrentPotentialField). Defaults to the default for the given field, see the docstring of the field object for the specific default.
eval_grid (Grid, optional) – Collocation grid containing the nodes to evaluate the normal magnetic field at plasma geometry at. Defaults to a LinearGrid(L=eq.L_grid, M=eq.M_grid, zeta=jnp.array(0.0), NFP=eq.NFP).
field_fixed (bool) – Whether to fix the field’s DOFs during the optimization.
eq_fixed (bool) – Whether to fix the equilibrium (or QFM surface) DOFs during the optimization.
bs_chunk_size (int or None) – Size to split Biot-Savart computation into chunks of evaluation points. If no chunking should be done or the chunk size is the full input then supply
None.target ({float, ndarray}, optional) – Target value(s) of the objective. Only used if
boundsisNone. Must be broadcastable toObjective.dim_f. Defaults totarget=eq.Psiif an Equilibrium is passed, ortarget=1.0if a surface.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=eq.Psiif an Equilibrium is passed, ortarget=1.0if a surface.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.
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.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. Note: has no effect for this objective.
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 (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_1[, params_2, constants])Compute toroidal flux.
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.