freegsnke.jtor_update module

Defines the FreeGSNKE profile Object, which inherits from the FreeGS4E profile object.

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.jtor_update.ConstrainBetapIp(eq, *args, **kwargs)[source]

Bases: ConstrainBetapIp, Jtor_universal

FreeGSNKE profile class adapting the original FreeGS object with the same name, with a few modifications, to: - retain memory of critical point calculation; - deal with limiter plasma configurations

Lao_parameters(n_alpha, n_beta, alpha_logic=True, beta_logic=True, Ip_logic=True, nn=100)[source]

Finds best fitting alpha, beta parameters for a Lao85 profile, to reproduce the input pprime_ and ffprime_ n_alpha and n_beta represent the number of free parameters

See Lao_parameters_finder.

__init__(eq, *args, **kwargs)[source]

Instantiates the object.

Parameters:

eq (FreeGSNKE Equilibrium object) – Specifies the domain properties

class freegsnke.jtor_update.ConstrainPaxisIp(eq, *args, **kwargs)[source]

Bases: ConstrainPaxisIp, Jtor_universal

FreeGSNKE profile class adapting the original FreeGS object with the same name, with a few modifications, to: - retain memory of critical point calculation; - deal with limiter plasma configurations

Lao_parameters(n_alpha, n_beta, alpha_logic=True, beta_logic=True, Ip_logic=True, nn=100)[source]

Finds best fitting alpha, beta parameters for a Lao85 profile, to reproduce the input pprime_ and ffprime_ n_alpha and n_beta represent the number of free parameters

See Lao_parameters_finder.

__init__(eq, *args, **kwargs)[source]

Instantiates the object.

Parameters:

eq (FreeGSNKE Equilibrium object) – Specifies the domain properties

class freegsnke.jtor_update.Fiesta_Topeol(eq, *args, **kwargs)[source]

Bases: Fiesta_Topeol, Jtor_universal

FreeGSNKE profile class adapting the FreeGS4E object with the same name, with a few modifications, to: - retain memory of critical point calculation; - deal with limiter plasma configurations

Lao_parameters(n_alpha, n_beta, alpha_logic=True, beta_logic=True, Ip_logic=True, nn=100)[source]

Finds best fitting alpha, beta parameters for a Lao85 profile, to reproduce the input pprime_ and ffprime_ n_alpha and n_beta represent the number of free parameters

See Lao_parameters_finder.

__init__(eq, *args, **kwargs)[source]

Instantiates the object.

Parameters:

eq (FreeGSNKE Equilibrium object) – Specifies the domain properties

class freegsnke.jtor_update.Jtor_universal(refine_jtor=False)[source]

Bases: object

Jtor_build(Jtor_part1, Jtor_part2, core_mask_limiter, R, Z, psi, psi_bndry, mask_outside_limiter, limiter_mask_out)[source]

Universal function that calculates the plasma current distribution, common to all of the different types of profile parametrizations used in FreeGSNKE.

Parameters:
  • Jtor_part1 (method) – method from the freegs4e Profile class returns opt, xpt, diverted_core_mask

  • Jtor_part2 (method) – method from each individual profile class returns jtor itself

  • core_mask_limiter (method) – method of the limiter_handler class returns the refined core_mask where jtor>0 accounting for the limiter

  • R (np.ndarray) – R coordinates of the domain grid points

  • Z (np.ndarray) – Z coordinates of the domain grid points

  • psi (np.ndarray) – Poloidal field flux / 2*pi at each grid points (for example as returned by Equilibrium.psi())

  • psi_bndry (float, optional) – Value of the poloidal field flux at the boundary of the plasma (last closed flux surface), by default None

  • mask_outside_limiter (np.ndarray) – Mask of points outside the limiter, if any, optional

  • limiter_mask_out (np.ndarray) – The mask identifying the border of the limiter, including points just inside it, the ‘last’ accessible to the plasma. Same size as psi.

Jtor_refined(R, Z, psi, psi_bndry=None, thresholds=None)[source]

Implements the call to the Jtor method for the case in which the subgrid refinement is used.

Parameters

Znp.ndarray

Z coordinates of the domain grid points

psinp.ndarray

Poloidal field flux / 2*pi at each grid points (for example as returned by Equilibrium.psi())

psi_bndryfloat, optional

Value of the poloidal field flux at the boundary of the plasma (last closed flux surface), by default None

thresholdstuple (threshold for jtor criterion, threshold for gradient criterion)

tuple of values used to identify where to apply refinement when None, the default refinement_thresholds are used

Returns:

2d map of toroidal current values

Return type:

ndarray

Jtor_unrefined(R, Z, psi, psi_bndry=None)[source]

Replaces the FreeGS4E call, while maintaining the same input structure.

Parameters:
  • R (np.ndarray) – R coordinates of the domain grid points

  • Z (np.ndarray) – Z coordinates of the domain grid points

  • psi (np.ndarray) – Poloidal field flux / 2*pi at each grid points (for example as returned by Equilibrium.psi())

  • psi_bndry (float, optional) – Value of the poloidal field flux at the boundary of the plasma (last closed flux surface), by default None

Returns:

2d map of toroidal current values

Return type:

ndarray

__init__(refine_jtor=False)[source]

Sets default unrefined Jtor.

select_refinement(eq, refine_jtor, nnx, nny)[source]

Initializes the object that handles the subgrid refinement of jtor

Parameters:
  • eq (freegs4e Equilibrium object) – Specifies the domain properties

  • refine_jtor (bool) – Flag to select whether to apply sug-grid refinement of plasma current distribution jtor

  • nnx (even integer) – refinement factor in the R direction

  • nny (even integer) – refinement factor in the Z direction

set_masks(eq)[source]

Universal function to set all masks related to the limiter.

Parameters:

eq (FreeGSNKE Equilibrium object) – Specifies the domain properties

set_refinement_thresholds(thresholds=(1.0, 1.0))[source]

Sets the default criteria for refinement – used when not directly set.

Parameters:

thresholds (tuple (threshold for jtor criterion, threshold for gradient criterion)) – tuple of values used to identify where to apply refinement

class freegsnke.jtor_update.Lao85(eq, *args, refine_jtor=False, nnx=None, nny=None, **kwargs)[source]

Bases: Lao85, Jtor_universal

FreeGSNKE profile class adapting the FreeGS4E object with the same name, with a few modifications, to: - retain memory of critical point calculation; - deal with limiter plasma configurations

Topeol_parameters(nn=100, max_it=100, tol=1e-05)[source]

Fids best combination of (alpha_m, alpha_n, beta_0) to instantiate a Topeol profile object as similar as possible to self

Parameters:
  • nn (int, optional) – number of points to sample 0,1 interval in the normalised psi, by default 100

  • max_it (int,) – maximum number of iterations in the optimization

  • tol (float) – iterations stop when change in the optimised parameters in smaller than tol

__init__(eq, *args, refine_jtor=False, nnx=None, nny=None, **kwargs)[source]

Instantiates the object.

Parameters:
  • eq (freegs4e Equilibrium object) – Specifies the domain properties

  • refine_jtor (bool) – Flag to select whether to apply sug-grid refinement of plasma current distribution jtor

  • nnx (even integer) – refinement factor in the R direction

  • nny (even integer) – refinement factor in the Z direction

class freegsnke.jtor_update.TensionSpline(eq, *args, **kwargs)[source]

Bases: TensionSpline, Jtor_universal

FreeGSNKE profile class adapting the FreeGS4E object with the same name, with a few modifications, to: - retain memory of critical point calculation; - deal with limiter plasma configurations

__init__(eq, *args, **kwargs)[source]

Instantiates the object.

Parameters:

eq (FreeGSNKE Equilibrium object) – Specifies the domain properties

assign_profile_parameter(pp_knots, pp_values, pp_values_2, pp_sigma, ffp_knots, ffp_values, ffp_values_2, ffp_sigma)[source]

Assigns to the profile object new values for the profile parameters