desc.objectives.FixSumCoilCurrent

class desc.objectives.FixSumCoilCurrent(coil, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, indices=True, name='summed coil current')Source

Fixes the sum of coil current(s) in a Coil or CoilSet.

NOTE: When using this objective, take care in knowing the signs of the current in the coils and the orientations of the coils. It is possible for coils with the same signs of their current to have currents flowing in differing directions in physical space due to the orientation of the coils.

Parameters:
  • coil (Coil) – Coil(s) that will be optimized to satisfy the Objective.

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

  • bounds (tuple of {float, ndarray}, optional) – Lower and upper bounds on the objective. Overrides target. Both bounds must be broadcastable to Objective.dim_f. Default is to use the target instead.

  • 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 normalize is True and the target is in physical units, this should also be set to True.

  • indices (nested list of bool, optional) – Pytree of bool specifying which coil currents to sum together. See the example for how to use this on a mixed coil set. If True/False sums all/none of the coil currents.

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

Examples

import numpy as np
from desc.coils import (
    CoilSet, FourierPlanarCoil, FourierRZCoil, FourierXYZCoil, MixedCoilSet
)
from desc.objectives import FixSumCoilCurrent

# toroidal field coil set with 4 coils
tf_coil = FourierPlanarCoil(
    current=3, center=[2, 0, 0], normal=[0, 1, 0], r_n=[1]
)
tf_coilset = CoilSet.linspaced_angular(tf_coil, n=4)
# vertical field coil set with 3 coils
vf_coil = FourierRZCoil(current=-1, R_n=3, Z_n=-1)
vf_coilset = CoilSet.linspaced_linear(
    vf_coil, displacement=[0, 0, 2], n=3, endpoint=True
)
# another single coil
xyz_coil = FourierXYZCoil(current=2)
# full coil set with TF coils, VF coils, and other single coil
full_coilset = MixedCoilSet((tf_coilset, vf_coilset, xyz_coil))

# equilibrium G(rho=1) determines the necessary net poloidal current through
# the coils (as dictated by Ampere's law)
# the sign convention is positive poloidal current flows up through the torus
# hole
grid_at_surf = LinearGrid(rho=1.0, M=eq.M_grid, N=eq.N_grid)
G_tot = 2*jnp.pi*eq.compute("G", grid=grid_at_surf)["G"][0] / mu_0

# to use this objective to satisfy Ampere's law for the targeted equilibrium,
# only coils that link the equilibrium poloidally should be included in the sum,
# which is the TF coil set and the FourierXYZ coil, but not the VF coil set
obj = FixSumCoilCurrent(full_coilset, indices=[True, False, True], target=G_tot)

Methods

build([use_jit, verbose])

Build constant arrays.

compute(params[, constants])

Compute sum of coil currents.

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.

update_target(thing)

Update target values using an Optimizable 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.