freegsnke.inverse module

Implements the optimiser for the inverse Grad-Shafranov problem.

Copyright 2025 UKAEA, UKRI-STFC, and The Authors, as per the COPYRIGHT and README files.

This file is part of FreeGSNKE.

FreeGSNKE is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

FreeGSNKE is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

You should have received a copy of the GNU Lesser General Public License along with FreeGSNKE. If not, see <http://www.gnu.org/licenses/>.

class freegsnke.inverse.Inverse_optimizer(isoflux_set=None, null_points=None, psi_vals=None, curr_vals=None)[source]

Bases: object

This class implements a gradient based optimiser for the coil currents, used to perform (static) inverse GS solves.

__init__(isoflux_set=None, null_points=None, psi_vals=None, curr_vals=None)[source]

Instantiates the object and sets all magnetic constraints to be used.

Parameters:
  • isoflux_set (list or np.array, optional) – list of isoflux objects, each with structure [Rcoords, Zcoords] with Rcoords and Zcoords being 1D lists of the coords of all points that are requested to have the same flux value

  • null_points (list or np.array, optional) – structure [Rcoords, Zcoords], with Rcoords and Zcoords being 1D lists Sets the coordinates of the desired null points, including both Xpoints and Opoints

  • psi_vals (list or np.array, optional) – structure [Rcoords, Zcoords, psi_values] with Rcoords, Zcoords and psi_values having the same shape Sets the desired values of psi for a set of coordinates, possibly an entire map

  • curr_vals (list, optional) – structure [[coil indexes in the array of coils available for control], [coil current values]]

build_control_coils(eq)[source]

Records what coils are available for control

Parameters:

eq (freegsnke equilibrium object)

build_control_currents(eq)[source]

Builds vector of coil current values, including only those coils that are available for control. Values are extracted from the equilibrium itself.

Parameters:

eq (freegsnke equilibrium object)

build_control_currents_Vec(full_currents_vec)[source]

Builds vector of coil current values, including only those coils that are available for control. Values are extracted from the full current vector.

Parameters:

full_currents_vec (np.array) – Vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

build_curr_vals_lsq(full_currents_vec)[source]

Builds for the ordinary least sq problem associated to the psi values

Parameters:

full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

build_full_current_vec(eq)[source]

Builds full vector of coil current values.

Parameters:

eq (freegsnke equilibrium object)

build_greens(eq)[source]

Calculates and stores all of the needed green function values.

Parameters:

eq (freegsnke equilibrium object)

build_isoflux_lsq(full_currents_vec)[source]

Builds for the ordinary least sq problem associated to the isoflux constraints

Parameters:

full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

build_lsq(full_currents_vec)[source]

Fetches all terms for the least sq problem, combining all types of magnetic constraints

Parameters:

full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

build_null_points_lsq(full_currents_vec)[source]

Builds for the ordinary least sq problem associated to the null points

Parameters:

full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

build_plasma_isoflux_lsq(full_currents_vec, trial_plasma_psi)[source]
build_plasma_vals(trial_plasma_psi)[source]

Builds and stores all the values relative to the plasma, based on the provided plasma_psi

Parameters:

trial_plasma_psi (np.array) – Flux due to the plasma. Same shape as eq.R

build_psi_vals_lsq(full_currents_vec)[source]

Builds for the ordinary least sq problem associated to the psi values

Parameters:

full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

optimize_currents(full_currents_vec, trial_plasma_psi, l2_reg)[source]

Solves the least square problem. Tikhonov regularization is applied.

Parameters:
  • full_currents_vec (np.array) – Full vector of all coil current values. For example as returned by eq.tokamak.getCurrentsVec()

  • trial_plasma_psi (np.array) – Flux due to the plasma. Same shape as eq.R

  • l2_reg (either float or 1d np.array with len=self.n_control_coils) – The regularization factor

optimize_plasma_psi(full_currents_vec, trial_plasma_psi, l2_reg)[source]
plot(axis=None, show=True)[source]

Plots constraints used for coil current control

axis - Specify the axis on which to plot show - Call matplotlib.pyplot.show() before returning

prepare_for_plasma_optimization(eq)[source]
prepare_for_solve(eq)[source]

To be called after object is instantiated. Prepares necessary quantities for loss/gradient calculations.

Parameters:

eq (freegsnke equilibrium object) – Sources information on: - coils available for control - coil current values - green functions

prepare_plasma_psi(trial_plasma_psi)[source]
prepare_plasma_vals_for_plasma(trial_plasma_psi)[source]
rebuild_full_current_vec(control_currents, filling=0)[source]

Builds a full_current vector using the input values. Only the coil currents of the coils available for control are filled in.

Parameters:

control_currents (np.array) – Vector of coil currents for those coils available for control.

source_domain_properties(eq)[source]