desc.objectives.FixParameters

class desc.objectives.FixParameters(thing, params=None, target=None, bounds=None, weight=1, normalize=True, normalize_target=True, name='fixed parameters')Source

Fix specific degrees of freedom associated with a given Optimizable thing.

Parameters:
  • thing (Optimizable) – Object whose degrees of freedom are being fixed.

  • params (nested list of dicts) – Dict keys are the names of parameters to fix (str), and dict values are the indices to fix for each corresponding parameter (int array). Use True (False) instead of an int array to fix all (none) of the indices for that parameter. Must have the same pytree structure as thing.params_dict. The default is to fix all indices of all parameters.

  • target (dict of {float, ndarray}, optional) – Target value(s) of the objective. Only used if bounds is None. Should have the same tree structure as thing.params. Defaults to things.params.

  • bounds (tuple of dict {float, ndarray}, optional) – Lower and upper bounds on the objective. Overrides target. Should have the same tree structure as thing.params.

  • weight (dict of {float, ndarray}, optional) – Weighting to apply to the Objective, relative to other Objectives. Should be a scalar or have the same tree structure as thing.params.

  • 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. Has no effect for this objective.

  • 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 FixParameters

# 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))

params = [
    [
        {"current": True},  # fix "current" of the 1st TF coil
        # fix "center" and one component of "normal" for the 2nd TF coil
        {"center": True, "normal": np.array([1])},
        {"r_n": True},  # fix radius of the 3rd TF coil
        {},  # fix nothing in the 4th TF coil
    ],
    {"shift": True, "rotmat": True},  # fix "shift" & "rotmat" for all VF coils
    # fix specified indices of "X_n" and "Z_n", but not "Y_n", for other coil
    {"X_n": np.array([1, 2]), "Y_n": False, "Z_n": np.array([0])},
]
obj = FixParameters(full_coilset, params)

Methods

build([use_jit, verbose])

Build constant arrays.

compute(params[, constants])

Compute fixed degree of freedom errors.

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.