The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Attributes | List of all members
TokaMaker Class Reference

Detailed Description

TokaMaker G-S solver class.

Public Member Functions

 __del__ (self)
 Free Fortran-side objects by calling reset() before object is deleted or GC'd.
 
 __init__ (self, OFT_env)
 Initialize TokaMaker object.
 
 abspsi_to_normalized (self, psi_in)
 Convert unnormalized \( \psi \) values to normalized \( \hat{\psi} \) values.
 
 alam (self)
 
 alam (self, value)
 
 area_integral (self, field, reg_mask=-1)
 Compute area integral of field over a specified region.
 
 calc_loopvoltage (self)
 Get plasma loop voltage.
 
 coil_reg_term (self, coffs, target=0.0, weight=1.0)
 Define coil current regularization term for the form \( target = \Sigma_i \alpha_i I_i \) to be used in set_coil_reg.
 
 diverted (self)
 
 eig_td (self, omega=-1.E4, neigs=4, include_bounds=True, pm=False, damping_scale=-1.0)
 Compute eigenvalues for the linearized time-dependent system.
 
 eig_wall (self, neigs=4, pm=False)
 Compute eigenvalues ( \( 1 / \tau_{L/R} \)) for conducting structures.
 
 flux_integral (self, psi_vals, field_vals)
 Compute area integral of field over a specified region.
 
 get_coil_currents (self)
 Get currents in each coil [A] and coil region [A-turns].
 
 get_coil_Lmat (self)
 Get mutual inductance matrix between coils.
 
 get_conductor_currents (self, psi, cell_centered=False)
 Get toroidal current density in conducting regions for a given \( \psi \).
 
 get_conductor_source (self, dpsi_dt)
 Get toroidal current density in conducting regions for a \( d \psi / dt \) source.
 
 get_delstar_curr (self, psi)
 Get toroidal current density from \( \psi \) through \( \Delta^{*} \) operator.
 
 get_field_eval (self, field_type)
 Create field interpolator for vector potential.
 
 get_globals (self)
 Get global plasma parameters.
 
 get_jtor_plasma (self)
 Get plasma toroidal current density for current equilibrium.
 
 get_profiles (self, psi=None, psi_pad=1.E-8, npsi=50)
 Get G-S source profiles.
 
 get_psi (self, normalized=True)
 Get poloidal flux values on node points.
 
 get_q (self, psi=None, psi_pad=0.02, npsi=50, compute_geo=False)
 Get q-profile at specified or uniformly spaced points.
 
 get_stats (self, lcfs_pad=None, li_normalization='std', geom_type='max', beta_Ip=None)
 Get information (Ip, q, kappa, etc.) about current G-S equilbirium.
 
 get_targets (self)
 Get global target values.
 
 get_vfixed (self)
 Get required vacuum flux values to balance fixed boundary equilibrium.
 
 get_xpoints (self)
 Get X-points.
 
 init_psi (self, r0=-1.0, z0=0.0, a=0.0, kappa=0.0, delta=0.0, curr_source=None)
 Initialize \(\psi\) using uniform current distributions.
 
 load_profiles (self, f_file='none', foffset=None, p_file='none', eta_file='none', f_NI_file='none')
 Load flux function profiles ( \(F*F'\) and \(P'\)) from files.
 
 plot_constraints (self, fig, ax, isoflux_color='tab:red', isoflux_marker='+', saddle_color='tab:green', saddle_marker='x')
 Plot geometry constraints.
 
 plot_eddy (self, fig, ax, psi=None, dpsi_dt=None, nlevels=40, colormap='jet', clabel=r ' $J_w$[$A/m^2$]', symmap=False)
 Plot contours of \(\hat{\psi}\).
 
 plot_machine (self, fig, ax, vacuum_color='whitesmoke', cond_color='gray', limiter_color='k', coil_color='gray', coil_colormap=None, coil_symmap=False, coil_scale=1.0, coil_clabel=r ' $I_C$[A]', colorbar=None)
 Plot machine geometry.
 
 plot_psi (self, fig, ax, psi=None, normalized=True, plasma_color=None, plasma_nlevels=8, plasma_levels=None, plasma_colormap=None, plasma_linestyles=None, vacuum_color='darkgray', vacuum_nlevels=8, vacuum_levels=None, vacuum_colormap=None, vacuum_linestyles=None, xpoint_color='k', xpoint_marker='x', xpoint_inactive_alpha=0.5, opoint_color='k', opoint_marker=' *')
 Plot contours of \(\hat{\psi}\).
 
 pnorm (self)
 
 pnorm (self, value)
 
 print_info (self, lcfs_pad=0.01, li_normalization='std', geom_type='max', beta_Ip=None)
 Print information (Ip, q, etc.) about current G-S equilbirium.
 
 psinorm_to_absolute (self, psi_in)
 Convert normalized \( \hat{\psi} \) values to unnormalized values \( \psi \).
 
 reset (self)
 Reset G-S object to enable loading a new mesh and coil configuration.
 
 sauter_fc (self, psi=None, psi_pad=0.02, npsi=50)
 Evaluate Sauter trapped particle fractions at specified or uniformly spaced points.
 
 save_eqdsk (self, filename, nr=65, nz=65, rbounds=None, zbounds=None, run_info='', lcfs_pad=0.01, rcentr=None, truncate_eq=True, limiter_file='', lcfs_pressure=0.0)
 Save current equilibrium to gEQDSK format.
 
 save_ifile (self, filename, npsi=65, ntheta=65, lcfs_pad=0.01, lcfs_pressure=0.0, pack_lcfs=True, single_precision=False)
 Save current equilibrium to iFile format.
 
 save_mug (self, filename)
 Save current equilibrium to MUG transfer format.
 
 set_coil_bounds (self, coil_bounds=None)
 Set hard constraints on coil currents.
 
 set_coil_current_dist (self, coil_name, curr_dist)
 Overwrite coil with non-uniform current distribution.
 
 set_coil_currents (self, currents=None)
 Set coil currents.
 
 set_coil_reg (self, reg_mat=None, reg_targets=None, reg_weights=None, reg_terms=None)
 Set regularization matrix for coil currents when isoflux and/or saddle constraints are used.
 
 set_coil_vsc (self, coil_gains)
 Define a vertical stability coil set from one or more coils.
 
 set_dipole_a (self, a_exp)
 Update anisotropy exponent a when dipole mode is used.
 
 set_flux (self, locations, targets, weights=None)
 Set explicit flux constraint points \( \psi(x_i) \).
 
 set_isoflux (self, isoflux, weights=None, grad_wt_lim=-1.0, ref_points=None)
 Set isoflux constraint points (all points lie on a flux surface)
 
 set_profiles (self, ffp_prof=None, foffset=None, pp_prof=None, ffp_NI_prof=None, keep_files=False)
 Set flux function profiles ( \(F*F'\) and \(P'\)) using a piecewise linear definition.
 
 set_psi (self, psi)
 Set poloidal flux values on node points.
 
 set_psi_dt (self, psi0, dt)
 Set reference poloidal flux and time step for eddy currents in .solve()
 
 set_resistivity (self, eta_prof=None)
 Set flux function profile $\eta$ using a piecewise linear definition.
 
 set_saddles (self, saddles, weights=None)
 Set saddle constraint points (poloidal field should vanish at each point)
 
 set_targets (self, Ip=None, Ip_ratio=None, pax=None, estore=None, R0=None, V0=None, retain_previous=False)
 Set global target values.
 
 setup (self, order=2, F0=0.0, full_domain=False)
 Setup G-S solver.
 
 setup_mesh (self, r=None, lc=None, reg=None, mesh_file=None)
 Setup mesh for static and time-dependent G-S calculations.
 
 setup_regions (self, cond_dict={}, coil_dict={})
 Define mesh regions (coils and conductors)
 
 setup_td (self, dt, lin_tol, nl_tol, pre_plasma=False)
 Setup the time-dependent G-S solver.
 
 solve (self, vacuum=False)
 Solve G-S equation with specified constraints, profiles, etc.
 
 step_td (self, time, dt)
 Compute eigenvalues for the time-dependent system.
 
 trace_surf (self, psi)
 Trace surface for a given poloidal flux.
 
 update_settings (self)
 Update settings after changes to values in python.
 
 vac_solve (self, psi=None, rhs_source=None)
 Solve for vacuum solution (no plasma), with present coil currents and optional other currents.
 

Public Attributes

 alam
 
 coil_set_names
 Coil set names in order of id number.
 
 coil_sets
 Coil set definitions, including sub-coils.
 
 dist_coils
 Distribution coils, only (currently) saved for plotting utility.
 
 lc
 Mesh triangles [nc,3].
 
 lim_contour
 Limiting contour.
 
 lim_contours
 
 lim_point
 Limiting point (limter or active X-point) [2].
 
 lim_pts
 
 nc
 Number of cells in mesh.
 
 ncoils
 Number of coils in mesh.
 
 np
 Number of points in mesh.
 
 nregs
 Number of regions in mesh.
 
 nvac
 Number of vacuum regions in mesh.
 
 o_point
 Location of O-point (magnetic axis) [2].
 
 pnorm
 
 psi_bounds
 Bounding values for \(\psi\) ( \(\psi_a\), \(\psi_0\)) [2].
 
 psi_convention
 Normalized flux convention (0 -> tokamak, 1 -> spheromak)
 
 r
 Mesh vertices [np,3] (last column should be all zeros)
 
 reg
 Mesh regions [nc].
 
 settings
 General settings object.
 
 x_points
 Location of X-points [20,2].
 

Protected Attributes

 _alam
 F*F' normalization value [1] (use alam property)
 
 _coil_dict
 Coil definition dictionary.
 
 _cond_dict
 Conductor definition dictionary.
 
 _diverted
 Diverted (limited) flag [1] (use diverted property)
 
 _estore_target
 Stored energy target value (use set_targets)
 
 _F0
 Vacuum F value.
 
 _flux_targets
 Flux constraint points (use set_flux)
 
 _Ip_ratio_target
 Plasma current target ratio I_p(FF') / I_p(P') (use set_targets)
 
 _Ip_target
 Plasma current target value (use set_targets)
 
 _isoflux_targets
 Isoflux constraint points (use set_isoflux)
 
 _mesh_ptr
 Internal mesh object.
 
 _oft_env
 
 _pax_target
 Axis pressure target value (use set_targets)
 
 _pnorm
 Pressure normalization value [1] (use pnorm property)
 
 _R0_target
 R0 target value (use set_targets)
 
 _saddle_targets
 Saddle constraint points (use set_saddles)
 
 _tMaker_ptr
 Internal Grad-Shafranov object (gs_eq)
 
 _V0_target
 V0 target value (use set_targets)
 
 _vac_dict
 Vacuum definition dictionary.
 
 _virtual_coils
 Virtual coils, if present (currently only ‘’#VSC'`)
 

Constructor & Destructor Documentation

◆ __init__()

__init__ (   self,
  OFT_env 
)

Initialize TokaMaker object.

Parameters
OFT_envOFT runtime environment object (See OFT_env)

◆ __del__()

__del__ (   self)

Free Fortran-side objects by calling reset() before object is deleted or GC'd.

Member Function Documentation

◆ abspsi_to_normalized()

abspsi_to_normalized (   self,
  psi_in 
)

Convert unnormalized \( \psi \) values to normalized \( \hat{\psi} \) values.

Parameters
psi_inInput \( \psi \) values
Returns
Normalized \( \hat{\psi} \) values

◆ alam() [1/2]

alam (   self)

◆ alam() [2/2]

alam (   self,
  value 
)

◆ area_integral()

area_integral (   self,
  field,
  reg_mask = -1 
)

Compute area integral of field over a specified region.

Parameters
fieldField to integrate [np,]
reg_maskID of region for integration (negative for whole mesh)
Returns
\( \int f dA \)

◆ calc_loopvoltage()

calc_loopvoltage (   self)

Get plasma loop voltage.

Parameters
etaDictionary object containing resistivity profile ['y'] and sampled locations in normalized Psi ['x']
ffp_NIDictionary object containing non-inductive FF' profile ['y'] and sampled locations in normalized Psi ['x']
Returns
Vloop [Volts]

◆ coil_reg_term()

coil_reg_term (   self,
  coffs,
  target = 0.0,
  weight = 1.0 
)

Define coil current regularization term for the form \( target = \Sigma_i \alpha_i I_i \) to be used in set_coil_reg.

Parameters
coffsDictionary of coefficients \( \alpha \) (zero for unspecified coils)
targetRegularization target (default: 0.0)
weightWeight for regularization term (default: 1.0)

◆ diverted()

diverted (   self)

◆ eig_td()

eig_td (   self,
  omega = -1.E4,
  neigs = 4,
  include_bounds = True,
  pm = False,
  damping_scale = -1.0 
)

Compute eigenvalues for the linearized time-dependent system.

Parameters
omegaGrowth rate localization point (eigenvalues closest to this value will be found)
neigsNumber of eigenvalues to compute
include_boundsInclude bounding flux terms for constant normalized profiles?
pmPrint solver statistics and raw eigenvalues?
damping_scaleScale factor for damping term to artificially limit growth rate (negative to disable)?
Returns
eigenvalues[neigs], eigenvectors[neigs,:]

◆ eig_wall()

eig_wall (   self,
  neigs = 4,
  pm = False 
)

Compute eigenvalues ( \( 1 / \tau_{L/R} \)) for conducting structures.

Parameters
neigsNumber of eigenvalues to compute
pmPrint solver statistics and raw eigenvalues?
Returns
eigenvalues[neigs], eigenvectors[neigs,:]

◆ flux_integral()

flux_integral (   self,
  psi_vals,
  field_vals 
)

Compute area integral of field over a specified region.

Parameters
fieldField to integrate [np,]
reg_maskID of region for integration (negative for whole mesh)
Returns
\( \int f dA \)

◆ get_coil_currents()

get_coil_currents (   self)

Get currents in each coil [A] and coil region [A-turns].

Returns
Coil currents [ncoils], Coil currents by region [nregs]

◆ get_coil_Lmat()

get_coil_Lmat (   self)

Get mutual inductance matrix between coils.

Note
This is the inductance in terms of A-turns. To get in terms of current in a single of the \(n\) windings you must multiply by \(n_i*n_j\).
Returns
L[ncoils+1,ncoils+1]

◆ get_conductor_currents()

get_conductor_currents (   self,
  psi,
  cell_centered = False 
)

Get toroidal current density in conducting regions for a given \( \psi \).

Parameters
psiPsi corresponding to field with conductor currents (eg. from time-dependent simulation)
cell_centeredGet currents at cell centers

◆ get_conductor_source()

get_conductor_source (   self,
  dpsi_dt 
)

Get toroidal current density in conducting regions for a \( d \psi / dt \) source.

Parameters
dpsi_dtdPsi/dt source eddy currents (eg. from linear stability)

◆ get_delstar_curr()

get_delstar_curr (   self,
  psi 
)

Get toroidal current density from \( \psi \) through \( \Delta^{*} \) operator.

Parameters
psi\( \psi \) corresponding to desired current density
Returns
\( J_{\phi} = \textrm{M}^{-1} \Delta^{*} \psi \)

◆ get_field_eval()

get_field_eval (   self,
  field_type 
)

Create field interpolator for vector potential.

Parameters
field_typeField to interpolate, must be one of ("B", "psi", "F", "P", "dPSI", "dBr", "dBt", or "dBz")
Returns
Field interpolation object

◆ get_globals()

get_globals (   self)

Get global plasma parameters.

Returns
Ip, [R_Ip, Z_Ip], \(\int dV\), \(\int P dV\), diamagnetic flux, enclosed toroidal flux

◆ get_jtor_plasma()

get_jtor_plasma (   self)

Get plasma toroidal current density for current equilibrium.

Returns
\( J_{\phi} \) by evalutating RHS source terms

◆ get_profiles()

get_profiles (   self,
  psi = None,
  psi_pad = 1.E-8,
  npsi = 50 
)

Get G-S source profiles.

Parameters
psiExplicit sampling locations in \(\hat{\psi}\)
psi_padEnd padding (axis and edge) for uniform sampling (ignored if psi is not None)
npsiNumber of points for uniform sampling (ignored if psi is not None)
Returns
\(\hat{\psi}\), \(F(\hat{\psi})\), \(F'(\hat{\psi})\), \(P(\hat{\psi})\), \(P'(\hat{\psi})\)

◆ get_psi()

get_psi (   self,
  normalized = True 
)

Get poloidal flux values on node points.

Parameters
normalizedNormalize (and offset) poloidal flux
Returns
\(\hat{\psi} = \frac{\psi-\psi_0}{\psi_a-\psi_0} \) or \(\psi\)

◆ get_q()

get_q (   self,
  psi = None,
  psi_pad = 0.02,
  npsi = 50,
  compute_geo = False 
)

Get q-profile at specified or uniformly spaced points.

Parameters
psiExplicit sampling locations in \(\hat{\psi}\)
psi_padEnd padding (axis and edge) for uniform sampling (ignored if psi is not None)
npsiNumber of points for uniform sampling (ignored if psi is not None)
compute_geoCompute geometric values for LCFS
Returns
\(\hat{\psi}\), \(q(\hat{\psi})\), \([<R>,<1/R>,dV/dPsi]\), length of last surface, [r(R_min),r(R_max)], [r(z_min),r(z_max)]

◆ get_stats()

get_stats (   self,
  lcfs_pad = None,
  li_normalization = 'std',
  geom_type = 'max',
  beta_Ip = None 
)

Get information (Ip, q, kappa, etc.) about current G-S equilbirium.

See eq. 1 for ‘li_normalization='std’and eq 2. forli_normalization='iter'` in Jackson et al.

Parameters
lcfs_padPadding at LCFS for boundary calculations (default: 1.0 for limited; 0.99 for diverted)
li_normalizationForm of normalized \( l_i \) ('std', 'ITER')
geom_typeMethod for computing geometric major/minor radius ('max': Use LCFS extrema, 'mid': Use axis plane extrema)
beta_IpOverride \( I_p \) used for beta calculations
Returns
Dictionary of equilibrium parameters

◆ get_targets()

get_targets (   self)

Get global target values.

Returns
Dictionary of global target values

◆ get_vfixed()

get_vfixed (   self)

Get required vacuum flux values to balance fixed boundary equilibrium.

Returns
sampling points [:,2], flux values [:]

◆ get_xpoints()

get_xpoints (   self)

Get X-points.

Returns
X-points, is diverted?

◆ init_psi()

init_psi (   self,
  r0 = -1.0,
  z0 = 0.0,
  a = 0.0,
  kappa = 0.0,
  delta = 0.0,
  curr_source = None 
)

Initialize \(\psi\) using uniform current distributions.

If r0>0 then a uniform current density inside a surface bounded by a curve of the form defined in oftpy.create_isoflux is used. If r0<0 then a uniform current density over the entire plasma region is used.

Parameters
r0Major radial position for flux surface-based approach
z0Vertical position for flux surface-based approach
aMinor radius for flux surface-based approach
kappaElongation for flux surface-based approach
deltaTriangularity for flux surface-based approach
curr_sourceCurrent source for arbitrary current distribution

◆ load_profiles()

load_profiles (   self,
  f_file = 'none',
  foffset = None,
  p_file = 'none',
  eta_file = 'none',
  f_NI_file = 'none' 
)

Load flux function profiles ( \(F*F'\) and \(P'\)) from files.

Parameters
f_fileFile containing \(F*F'\) (or \(F'\) if mode=0) definition
foffsetValue of \(F0=R0*B0\)
p_fileFile containing \(P'\) definition
eta_fileFile containing $\eta$ definition
f_NI_fileFile containing non-inductive \(F*F'\) definition

◆ plot_constraints()

plot_constraints (   self,
  fig,
  ax,
  isoflux_color = 'tab:red',
  isoflux_marker = '+',
  saddle_color = 'tab:green',
  saddle_marker = 'x' 
)

Plot geometry constraints.

Parameters
figFigure to add to
axAxis to add to
isoflux_colorColor of isoflux points (None to disable)
saddle_colorColor of saddle points (None to disable)

◆ plot_eddy()

plot_eddy (   self,
  fig,
  ax,
  psi = None,
  dpsi_dt = None,
  nlevels = 40,
  colormap = 'jet',
  clabel = r'$J_w$ [$A/m^2$]',
  symmap = False 
)

Plot contours of \(\hat{\psi}\).

Parameters
figFigure to add to
axAxis to add to
psiPsi corresponding to eddy currents (eg. from time-dependent simulation)
dpsi_dtdPsi/dt source eddy currents (eg. from linear stability)
nlevelsNumber contour lines used for shading (with "psi" only)
colormapColormap to use for shadings
clabelLabel for colorbar (None to disable colorbar)
Returns
Colorbar object

◆ plot_machine()

plot_machine (   self,
  fig,
  ax,
  vacuum_color = 'whitesmoke',
  cond_color = 'gray',
  limiter_color = 'k',
  coil_color = 'gray',
  coil_colormap = None,
  coil_symmap = False,
  coil_scale = 1.0,
  coil_clabel = r'$I_C$ [A]',
  colorbar = None 
)

Plot machine geometry.

Parameters
figFigure to add to
axAxis to add to
vacuum_colorColor to shade vacuum region (None to disable)
cond_colorColor for conducting regions (None to disable)
limiter_colorColor for limiter contour (None to disable)
coil_colorColor for coil regions (None to disable)
coil_colormapColormap for coil current values
coil_symmapMake coil current colorscale symmetric
coil_scaleScale for coil currents when plotting
coil_clabelLabel for coil current colorbar (None to disable colorbar)
colorbarColorbar instance to overwrite (None to add)
Returns
Colorbar instance for coil colors or None

◆ plot_psi()

plot_psi (   self,
  fig,
  ax,
  psi = None,
  normalized = True,
  plasma_color = None,
  plasma_nlevels = 8,
  plasma_levels = None,
  plasma_colormap = None,
  plasma_linestyles = None,
  vacuum_color = 'darkgray',
  vacuum_nlevels = 8,
  vacuum_levels = None,
  vacuum_colormap = None,
  vacuum_linestyles = None,
  xpoint_color = 'k',
  xpoint_marker = 'x',
  xpoint_inactive_alpha = 0.5,
  opoint_color = 'k',
  opoint_marker = '*' 
)

Plot contours of \(\hat{\psi}\).

Parameters
figFigure to add to
axAxis to add to
psiFlux values to plot (otherwise self.get_psi() is called)
normalizedRetreive normalized flux, or assume normalized psi if passed as argument
plasma_colorColor for plasma contours
plasma_nlevelsNumber of plasma contours
plasma_levelsExplicit levels for plasma contours
plasma_colormapColormap for plasma contours (cannot be specified with plasma_color)
plasma_linestylesLinestyle for plasma contours
vacuum_colorColor for vacuum contours
vacuum_nlevelsNumber of vacuum contours
vacuum_levelsExplicit levels for vacuum contours (cannot be specified with vacuum_color)
vacuum_colormapColormap for vacuum contours
vacuum_linestylesLinestyle for vacuum contours
xpoint_colorColor for X-point markers (None to disable)
xpoint_markerMarker style for X-points
xpoint_inactive_alphaAlpha value for inactive X-points
opoint_colorColor for O-point markers (None to disable)
opoint_markerMarker style for O-points

◆ pnorm() [1/2]

pnorm (   self)

◆ pnorm() [2/2]

pnorm (   self,
  value 
)

◆ print_info()

print_info (   self,
  lcfs_pad = 0.01,
  li_normalization = 'std',
  geom_type = 'max',
  beta_Ip = None 
)

Print information (Ip, q, etc.) about current G-S equilbirium.

Parameters
lcfs_padPadding at LCFS for boundary calculations
li_normalizationForm of normalized \( l_i \) ('std', 'ITER')
geom_typeMethod for computing geometric major/minor radius ('max': Use LCFS extrema, 'mid': Use axis plane extrema)
beta_IpOverride \( I_p \) used for beta calculations

◆ psinorm_to_absolute()

psinorm_to_absolute (   self,
  psi_in 
)

Convert normalized \( \hat{\psi} \) values to unnormalized values \( \psi \).

Parameters
psi_inInput \( \hat{\psi} \) values
Returns
Unnormalized \( \psi \) values

◆ reset()

reset (   self)

Reset G-S object to enable loading a new mesh and coil configuration.

◆ sauter_fc()

sauter_fc (   self,
  psi = None,
  psi_pad = 0.02,
  npsi = 50 
)

Evaluate Sauter trapped particle fractions at specified or uniformly spaced points.

Parameters
psiExplicit sampling locations in \(\hat{\psi}\)
psi_padEnd padding (axis and edge) for uniform sampling (ignored if psi is not None)
npsiNumber of points for uniform sampling (ignored if psi is not None)
Returns
\( f_c \), [ \(<R>,<1/R>,<a>\)], [ \(<|B|>,<|B|^2>\)]

◆ save_eqdsk()

save_eqdsk (   self,
  filename,
  nr = 65,
  nz = 65,
  rbounds = None,
  zbounds = None,
  run_info = '',
  lcfs_pad = 0.01,
  rcentr = None,
  truncate_eq = True,
  limiter_file = '',
  lcfs_pressure = 0.0 
)

Save current equilibrium to gEQDSK format.

Parameters
filenameFilename to save equilibrium to
nrNumber of radial sampling points
nzNumber of vertical sampling points
rboundsExtents of grid in R
zboundsExtents of grid in Z
run_infoRun information for gEQDSK file (maximum of 40 characters)
lcfs_padPadding in normalized flux at LCFS
rcentrRCENTR value for gEQDSK file (if None, geometric axis is used)
truncate_eqTruncate equilibrium at lcfs_pad, if False \( q(\hat{\psi} > 1-pad) = q(1-pad) \)
limiter_fileFile containing limiter contour to use instead of TokaMaker limiter
lcfs_pressurePlasma pressure on the LCFS (zero by default)

◆ save_ifile()

save_ifile (   self,
  filename,
  npsi = 65,
  ntheta = 65,
  lcfs_pad = 0.01,
  lcfs_pressure = 0.0,
  pack_lcfs = True,
  single_precision = False 
)

Save current equilibrium to iFile format.

Parameters
filenameFilename to save equilibrium to
npsiNumber of radial sampling points
nthetaNumber of vertical sampling points
lcfs_padPadding in normalized flux at LCFS
lcfs_pressurePlasma pressure on the LCFS (zero by default)
pack_lcfsPack toward LCFS with quadraturic sampling?
single_precisionSave single precision file? (default: double precision)

◆ save_mug()

save_mug (   self,
  filename 
)

Save current equilibrium to MUG transfer format.

Parameters
filenameFilename to save equilibrium to

◆ set_coil_bounds()

set_coil_bounds (   self,
  coil_bounds = None 
)

Set hard constraints on coil currents.

Can be used with or without regularization terms (see set_coil_reg).

Parameters
coil_boundsMinimum and maximum allowable coil currents [ncoils+1,2]

◆ set_coil_current_dist()

set_coil_current_dist (   self,
  coil_name,
  curr_dist 
)

Overwrite coil with non-uniform current distribution.

Parameters
coil_nameName of coil to modify
curr_distRelative current density [self.np]

◆ set_coil_currents()

set_coil_currents (   self,
  currents = None 
)

Set coil currents.

Parameters
currentsCurrent in each coil [A]

◆ set_coil_reg()

set_coil_reg (   self,
  reg_mat = None,
  reg_targets = None,
  reg_weights = None,
  reg_terms = None 
)

Set regularization matrix for coil currents when isoflux and/or saddle constraints are used.

Can be used to enforce "soft" constraints on coil currents. For hard constraints see set_coil_bounds.

Parameters
reg_matRegularization matrix [nregularize,ncoils+1]
reg_targetsRegularization targets [nregularize] (default: 0)
reg_weightsWeights for regularization terms [nregularize] (default: 1)
reg_termsList of regularization terms created with coil_reg_term

◆ set_coil_vsc()

set_coil_vsc (   self,
  coil_gains 
)

Define a vertical stability coil set from one or more coils.

Parameters
coil_gainsGains for each coil (absolute scale is arbitrary)

◆ set_dipole_a()

set_dipole_a (   self,
  a_exp 
)

Update anisotropy exponent a when dipole mode is used.

Parameters
a_expNew value for a exponent

◆ set_flux()

set_flux (   self,
  locations,
  targets,
  weights = None 
)

Set explicit flux constraint points \( \psi(x_i) \).

Parameters
locationsList of points defining constraints [:,2]
targetsTarget \( \psi \) value at each point [:]
weightsWeight to be applied to each constraint point [:] (default: 1)

◆ set_isoflux()

set_isoflux (   self,
  isoflux,
  weights = None,
  grad_wt_lim = -1.0,
  ref_points = None 
)

Set isoflux constraint points (all points lie on a flux surface)

To constraint points more uniformly in space additional weighting based on the gradient of $\psi$ at each point can also be included by setting grad_wt_lim>0. When set the actual weight will be $w_i * min(grad_wt_lim,|\nabla \psi|_{max} / |\nabla \psi|_i)$

Parameters
isofluxList of points defining constraints [:,2]
weightsWeight to be applied to each constraint point [:] (default: 1)
grad_wt_limLimit on gradient-based weighting (negative to disable)
ref_pointsReference points for each isoflux point [:,2] (default: isoflux[0,:] is used for all points)

◆ set_profiles()

set_profiles (   self,
  ffp_prof = None,
  foffset = None,
  pp_prof = None,
  ffp_NI_prof = None,
  keep_files = False 
)

Set flux function profiles ( \(F*F'\) and \(P'\)) using a piecewise linear definition.

Parameters
ffp_profDictionary object containing FF' profile ['y'] and sampled locations in normalized Psi ['x']
foffsetValue of \(F0=R0*B0\)
pp_profDictionary object containing P' profile ['y'] and sampled locations in normalized Psi ['x']
ffp_NI_profDictionary object containing non-inductive FF' profile ['y'] and sampled locations in normalized Psi ['x']
keep_filesRetain temporary profile files

◆ set_psi()

set_psi (   self,
  psi 
)

Set poloidal flux values on node points.

Parameters
psiPoloidal flux values (should not be normalized!)

◆ set_psi_dt()

set_psi_dt (   self,
  psi0,
  dt 
)

Set reference poloidal flux and time step for eddy currents in .solve()

Parameters
psi0Reference poloidal flux at t-dt (unnormalized)
dtTime since reference poloidal flux

◆ set_resistivity()

set_resistivity (   self,
  eta_prof = None 
)

Set flux function profile $\eta$ using a piecewise linear definition.

Arrays should have the form array[i,:] = ( \(\hat{\psi}_i\), \(f(\hat{\psi}_i)\)) and span \(\hat{\psi}_i = [0,1]\).

Parameters
eta_profValues defining $\eta$ [:,2]

◆ set_saddles()

set_saddles (   self,
  saddles,
  weights = None 
)

Set saddle constraint points (poloidal field should vanish at each point)

Parameters
saddlesList of points defining constraints [:,2]
weightsWeight to be applied to each constraint point [:] (default: 1)

◆ set_targets()

set_targets (   self,
  Ip = None,
  Ip_ratio = None,
  pax = None,
  estore = None,
  R0 = None,
  V0 = None,
  retain_previous = False 
)

Set global target values.

Note
Values that are not specified are reset to their defaults on each call unless retain_previous=True.
Parameters
alamScale factor for \(F*F'\) term (disabled if Ip is set)
pnormScale factor for \(P'\) term (disabled if pax, estore, or R0 are set)
IpTarget plasma current [A] (disabled if OFT_env.float_disable_flag)
Ip_ratioAmplitude of net plasma current contribution from FF' compared to P' (disabled if OFT_env.float_disable_flag)
paxTarget axis pressure [Pa] (disabled if OFT_env.float_disable_flag or if estore is set)
estoreTarget sotred energy [J] (disabled if OFT_env.float_disable_flag)
R0Target major radius for magnetic axis (disabled if OFT_env.float_disable_flag or if pax or estore are set)
V0Target vertical position for magnetic axis (disabled if OFT_env.float_disable_flag)
retain_previousKeep previously set targets unless explicitly updated? (default: False)

◆ setup()

setup (   self,
  order = 2,
  F0 = 0.0,
  full_domain = False 
)

Setup G-S solver.

Parameters
orderOrder of FE representation to use
F0Vacuum \(F(\psi)\) value (B0*R0)

◆ setup_mesh()

setup_mesh (   self,
  r = None,
  lc = None,
  reg = None,
  mesh_file = None 
)

Setup mesh for static and time-dependent G-S calculations.

A mesh should be specified by passing "r", "lc", and optionally "reg" or using a "mesh_file". When a region is specified the following ordering should apply:

  • 1: Plasma region
  • 2: Vacuum/air regions
  • 3+: Conducting regions and coils
Parameters
rMesh point list [np,2]
lcMesh cell list [nc,3] (base one)
regMesh region list [nc] (base one)
mesh_fileFilename containing mesh to load (native format only)

◆ setup_regions()

setup_regions (   self,
  cond_dict = {},
  coil_dict = {} 
)

Define mesh regions (coils and conductors)

Parameters
cond_dictDictionary specifying conducting regions

◆ setup_td()

setup_td (   self,
  dt,
  lin_tol,
  nl_tol,
  pre_plasma = False 
)

Setup the time-dependent G-S solver.

Parameters
dtStarting time step
lin_tolTolerance for linear solver
nl_tolTolerance for non-linear solver
pre_plasmaUse plasma contributions in preconditioner (default: False)

◆ solve()

solve (   self,
  vacuum = False 
)

Solve G-S equation with specified constraints, profiles, etc.

◆ step_td()

step_td (   self,
  time,
  dt 
)

Compute eigenvalues for the time-dependent system.

Parameters
timeGrowth rate enhancement point (should be approximately expected value)
dtNumber of eigenvalues to compute
Returns
new time, new dt, # of NL iterations, # of linear iterations, # of retries

◆ trace_surf()

trace_surf (   self,
  psi 
)

Trace surface for a given poloidal flux.

Parameters
psiFlux surface to trace \(\hat{\psi}\)
Returns
\(r(\hat{\psi})\)

◆ update_settings()

update_settings (   self)

Update settings after changes to values in python.

◆ vac_solve()

vac_solve (   self,
  psi = None,
  rhs_source = None 
)

Solve for vacuum solution (no plasma), with present coil currents and optional other currents.

Parameters
psiBoundary values for vacuum solve
rhs_sourceCurrent source (optional)

Member Data Documentation

◆ _alam

_alam
protected

F*F' normalization value [1] (use alam property)

◆ _coil_dict

_coil_dict
protected

Coil definition dictionary.

◆ _cond_dict

_cond_dict
protected

Conductor definition dictionary.

◆ _diverted

_diverted
protected

Diverted (limited) flag [1] (use diverted property)

◆ _estore_target

_estore_target
protected

Stored energy target value (use set_targets)

◆ _F0

_F0
protected

Vacuum F value.

◆ _flux_targets

_flux_targets
protected

Flux constraint points (use set_flux)

◆ _Ip_ratio_target

_Ip_ratio_target
protected

Plasma current target ratio I_p(FF') / I_p(P') (use set_targets)

◆ _Ip_target

_Ip_target
protected

Plasma current target value (use set_targets)

◆ _isoflux_targets

_isoflux_targets
protected

Isoflux constraint points (use set_isoflux)

◆ _mesh_ptr

_mesh_ptr
protected

Internal mesh object.

◆ _oft_env

_oft_env
protected

◆ _pax_target

_pax_target
protected

Axis pressure target value (use set_targets)

◆ _pnorm

_pnorm
protected

Pressure normalization value [1] (use pnorm property)

◆ _R0_target

_R0_target
protected

R0 target value (use set_targets)

◆ _saddle_targets

_saddle_targets
protected

Saddle constraint points (use set_saddles)

◆ _tMaker_ptr

_tMaker_ptr
protected

Internal Grad-Shafranov object (gs_eq)

◆ _V0_target

_V0_target
protected

V0 target value (use set_targets)

◆ _vac_dict

_vac_dict
protected

Vacuum definition dictionary.

◆ _virtual_coils

_virtual_coils
protected

Virtual coils, if present (currently only ‘’#VSC'`)

◆ alam

alam

◆ coil_set_names

coil_set_names

Coil set names in order of id number.

◆ coil_sets

coil_sets

Coil set definitions, including sub-coils.

◆ dist_coils

dist_coils

Distribution coils, only (currently) saved for plotting utility.

◆ lc

lc

Mesh triangles [nc,3].

◆ lim_contour

lim_contour

Limiting contour.

◆ lim_contours

lim_contours

◆ lim_point

lim_point

Limiting point (limter or active X-point) [2].

◆ lim_pts

lim_pts

◆ nc

nc

Number of cells in mesh.

◆ ncoils

ncoils

Number of coils in mesh.

◆ np

np

Number of points in mesh.

◆ nregs

nregs

Number of regions in mesh.

◆ nvac

nvac

Number of vacuum regions in mesh.

◆ o_point

o_point

Location of O-point (magnetic axis) [2].

◆ pnorm

pnorm

◆ psi_bounds

psi_bounds

Bounding values for \(\psi\) ( \(\psi_a\), \(\psi_0\)) [2].

◆ psi_convention

psi_convention

Normalized flux convention (0 -> tokamak, 1 -> spheromak)

◆ r

r

Mesh vertices [np,3] (last column should be all zeros)

◆ reg

reg

Mesh regions [nc].

◆ settings

settings

General settings object.

◆ x_points

x_points

Location of X-points [20,2].


The documentation for this class was generated from the following file: