|
The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
|
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'`) | |
| __init__ | ( | self, | |
| OFT_env | |||
| ) |
Initialize TokaMaker object.
| OFT_env | OFT runtime environment object (See OFT_env) |
| __del__ | ( | self | ) |
Free Fortran-side objects by calling reset() before object is deleted or GC'd.
| abspsi_to_normalized | ( | self, | |
| psi_in | |||
| ) |
Convert unnormalized \( \psi \) values to normalized \( \hat{\psi} \) values.
| psi_in | Input \( \psi \) values |
| alam | ( | self | ) |
| alam | ( | self, | |
| value | |||
| ) |
| area_integral | ( | self, | |
| field, | |||
reg_mask = -1 |
|||
| ) |
Compute area integral of field over a specified region.
| field | Field to integrate [np,] |
| reg_mask | ID of region for integration (negative for whole mesh) |
| calc_loopvoltage | ( | self | ) |
Get plasma loop voltage.
| eta | Dictionary object containing resistivity profile ['y'] and sampled locations in normalized Psi ['x'] |
| ffp_NI | Dictionary object containing non-inductive FF' profile ['y'] and sampled locations in normalized Psi ['x'] |
| 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.
| coffs | Dictionary of coefficients \( \alpha \) (zero for unspecified coils) |
| target | Regularization target (default: 0.0) |
| weight | Weight for regularization term (default: 1.0) |
| 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.
| omega | Growth rate localization point (eigenvalues closest to this value will be found) |
| neigs | Number of eigenvalues to compute |
| include_bounds | Include bounding flux terms for constant normalized profiles? |
| pm | Print solver statistics and raw eigenvalues? |
| damping_scale | Scale factor for damping term to artificially limit growth rate (negative to disable)? |
| eig_wall | ( | self, | |
neigs = 4, |
|||
pm = False |
|||
| ) |
Compute eigenvalues ( \( 1 / \tau_{L/R} \)) for conducting structures.
| neigs | Number of eigenvalues to compute |
| pm | Print solver statistics and raw eigenvalues? |
| flux_integral | ( | self, | |
| psi_vals, | |||
| field_vals | |||
| ) |
Compute area integral of field over a specified region.
| field | Field to integrate [np,] |
| reg_mask | ID of region for integration (negative for whole mesh) |
| 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 \).
| psi | Psi corresponding to field with conductor currents (eg. from time-dependent simulation) |
| cell_centered | Get currents at cell centers |
| get_conductor_source | ( | self, | |
| dpsi_dt | |||
| ) |
Get toroidal current density in conducting regions for a \( d \psi / dt \) source.
| dpsi_dt | dPsi/dt source eddy currents (eg. from linear stability) |
| get_delstar_curr | ( | self, | |
| psi | |||
| ) |
Get toroidal current density from \( \psi \) through \( \Delta^{*} \) operator.
| psi | \( \psi \) corresponding to desired current density |
| get_field_eval | ( | self, | |
| field_type | |||
| ) |
Create field interpolator for vector potential.
| field_type | Field to interpolate, must be one of ("B", "psi", "F", "P", "dPSI", "dBr", "dBt", or "dBz") |
| 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.
| psi | Explicit sampling locations in \(\hat{\psi}\) |
| psi_pad | End padding (axis and edge) for uniform sampling (ignored if psi is not None) |
| npsi | Number of points for uniform sampling (ignored if psi is not None) |
| get_psi | ( | self, | |
normalized = True |
|||
| ) |
Get poloidal flux values on node points.
| normalized | Normalize (and offset) poloidal flux |
| get_q | ( | self, | |
psi = None, |
|||
psi_pad = 0.02, |
|||
npsi = 50, |
|||
compute_geo = False |
|||
| ) |
Get q-profile at specified or uniformly spaced points.
| psi | Explicit sampling locations in \(\hat{\psi}\) |
| psi_pad | End padding (axis and edge) for uniform sampling (ignored if psi is not None) |
| npsi | Number of points for uniform sampling (ignored if psi is not None) |
| compute_geo | Compute geometric values for LCFS |
| 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.
| lcfs_pad | Padding at LCFS for boundary calculations (default: 1.0 for limited; 0.99 for diverted) |
| li_normalization | Form of normalized \( l_i \) ('std', 'ITER') |
| geom_type | Method for computing geometric major/minor radius ('max': Use LCFS extrema, 'mid': Use axis plane extrema) |
| beta_Ip | Override \( I_p \) used for beta calculations |
| 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.
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.
| r0 | Major radial position for flux surface-based approach |
| z0 | Vertical position for flux surface-based approach |
| a | Minor radius for flux surface-based approach |
| kappa | Elongation for flux surface-based approach |
| delta | Triangularity for flux surface-based approach |
| curr_source | Current source for arbitrary current distribution |
| 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.
| f_file | File containing \(F*F'\) (or \(F'\) if mode=0) definition |
| foffset | Value of \(F0=R0*B0\) |
| p_file | File containing \(P'\) definition |
| eta_file | File containing $\eta$ definition |
| f_NI_file | File containing non-inductive \(F*F'\) definition |
| plot_constraints | ( | self, | |
| fig, | |||
| ax, | |||
isoflux_color = 'tab:red', |
|||
isoflux_marker = '+', |
|||
saddle_color = 'tab:green', |
|||
saddle_marker = 'x' |
|||
| ) |
Plot geometry constraints.
| fig | Figure to add to |
| ax | Axis to add to |
| isoflux_color | Color of isoflux points (None to disable) |
| saddle_color | Color of saddle points (None to disable) |
| 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}\).
| fig | Figure to add to |
| ax | Axis to add to |
| psi | Psi corresponding to eddy currents (eg. from time-dependent simulation) |
| dpsi_dt | dPsi/dt source eddy currents (eg. from linear stability) |
| nlevels | Number contour lines used for shading (with "psi" only) |
| colormap | Colormap to use for shadings |
| clabel | Label for colorbar (None to disable colorbar) |
| 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.
| fig | Figure to add to |
| ax | Axis to add to |
| vacuum_color | Color to shade vacuum region (None to disable) |
| cond_color | Color for conducting regions (None to disable) |
| limiter_color | Color for limiter contour (None to disable) |
| coil_color | Color for coil regions (None to disable) |
| coil_colormap | Colormap for coil current values |
| coil_symmap | Make coil current colorscale symmetric |
| coil_scale | Scale for coil currents when plotting |
| coil_clabel | Label for coil current colorbar (None to disable colorbar) |
| colorbar | Colorbar instance to overwrite (None to add) |
| 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}\).
| fig | Figure to add to |
| ax | Axis to add to |
| psi | Flux values to plot (otherwise self.get_psi() is called) |
| normalized | Retreive normalized flux, or assume normalized psi if passed as argument |
| plasma_color | Color for plasma contours |
| plasma_nlevels | Number of plasma contours |
| plasma_levels | Explicit levels for plasma contours |
| plasma_colormap | Colormap for plasma contours (cannot be specified with plasma_color) |
| plasma_linestyles | Linestyle for plasma contours |
| vacuum_color | Color for vacuum contours |
| vacuum_nlevels | Number of vacuum contours |
| vacuum_levels | Explicit levels for vacuum contours (cannot be specified with vacuum_color) |
| vacuum_colormap | Colormap for vacuum contours |
| vacuum_linestyles | Linestyle for vacuum contours |
| xpoint_color | Color for X-point markers (None to disable) |
| xpoint_marker | Marker style for X-points |
| xpoint_inactive_alpha | Alpha value for inactive X-points |
| opoint_color | Color for O-point markers (None to disable) |
| opoint_marker | Marker style for O-points |
| 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.
| lcfs_pad | Padding at LCFS for boundary calculations |
| li_normalization | Form of normalized \( l_i \) ('std', 'ITER') |
| geom_type | Method for computing geometric major/minor radius ('max': Use LCFS extrema, 'mid': Use axis plane extrema) |
| beta_Ip | Override \( I_p \) used for beta calculations |
| psinorm_to_absolute | ( | self, | |
| psi_in | |||
| ) |
Convert normalized \( \hat{\psi} \) values to unnormalized values \( \psi \).
| psi_in | Input \( \hat{\psi} \) values |
| 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.
| psi | Explicit sampling locations in \(\hat{\psi}\) |
| psi_pad | End padding (axis and edge) for uniform sampling (ignored if psi is not None) |
| npsi | Number of points for uniform sampling (ignored if psi is not None) |
| 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.
| filename | Filename to save equilibrium to |
| nr | Number of radial sampling points |
| nz | Number of vertical sampling points |
| rbounds | Extents of grid in R |
| zbounds | Extents of grid in Z |
| run_info | Run information for gEQDSK file (maximum of 40 characters) |
| lcfs_pad | Padding in normalized flux at LCFS |
| rcentr | RCENTR value for gEQDSK file (if None, geometric axis is used) |
| truncate_eq | Truncate equilibrium at lcfs_pad, if False \( q(\hat{\psi} > 1-pad) = q(1-pad) \) |
| limiter_file | File containing limiter contour to use instead of TokaMaker limiter |
| lcfs_pressure | Plasma pressure on the LCFS (zero by default) |
| 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.
| filename | Filename to save equilibrium to |
| npsi | Number of radial sampling points |
| ntheta | Number of vertical sampling points |
| lcfs_pad | Padding in normalized flux at LCFS |
| lcfs_pressure | Plasma pressure on the LCFS (zero by default) |
| pack_lcfs | Pack toward LCFS with quadraturic sampling? |
| single_precision | Save single precision file? (default: double precision) |
| save_mug | ( | self, | |
| filename | |||
| ) |
Save current equilibrium to MUG transfer format.
| filename | Filename to save equilibrium to |
| 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).
| coil_bounds | Minimum and maximum allowable coil currents [ncoils+1,2] |
| set_coil_current_dist | ( | self, | |
| coil_name, | |||
| curr_dist | |||
| ) |
Overwrite coil with non-uniform current distribution.
| coil_name | Name of coil to modify |
| curr_dist | Relative current density [self.np] |
| set_coil_currents | ( | self, | |
currents = None |
|||
| ) |
Set coil currents.
| currents | Current in each coil [A] |
| 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.
| reg_mat | Regularization matrix [nregularize,ncoils+1] |
| reg_targets | Regularization targets [nregularize] (default: 0) |
| reg_weights | Weights for regularization terms [nregularize] (default: 1) |
| reg_terms | List of regularization terms created with coil_reg_term |
| set_coil_vsc | ( | self, | |
| coil_gains | |||
| ) |
Define a vertical stability coil set from one or more coils.
| coil_gains | Gains for each coil (absolute scale is arbitrary) |
| set_dipole_a | ( | self, | |
| a_exp | |||
| ) |
Update anisotropy exponent a when dipole mode is used.
| a_exp | New value for a exponent |
| set_flux | ( | self, | |
| locations, | |||
| targets, | |||
weights = None |
|||
| ) |
Set explicit flux constraint points \( \psi(x_i) \).
| locations | List of points defining constraints [:,2] |
| targets | Target \( \psi \) value at each point [:] |
| weights | Weight to be applied to each constraint point [:] (default: 1) |
| 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)$
| isoflux | List of points defining constraints [:,2] |
| weights | Weight to be applied to each constraint point [:] (default: 1) |
| grad_wt_lim | Limit on gradient-based weighting (negative to disable) |
| ref_points | Reference points for each isoflux point [:,2] (default: isoflux[0,:] is used for all points) |
| 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.
| ffp_prof | Dictionary object containing FF' profile ['y'] and sampled locations in normalized Psi ['x'] |
| foffset | Value of \(F0=R0*B0\) |
| pp_prof | Dictionary object containing P' profile ['y'] and sampled locations in normalized Psi ['x'] |
| ffp_NI_prof | Dictionary object containing non-inductive FF' profile ['y'] and sampled locations in normalized Psi ['x'] |
| keep_files | Retain temporary profile files |
| set_psi | ( | self, | |
| psi | |||
| ) |
Set poloidal flux values on node points.
| psi | Poloidal flux values (should not be normalized!) |
| set_psi_dt | ( | self, | |
| psi0, | |||
| dt | |||
| ) |
Set reference poloidal flux and time step for eddy currents in .solve()
| psi0 | Reference poloidal flux at t-dt (unnormalized) |
| dt | Time since reference poloidal flux |
| 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]\).
| eta_prof | Values defining $\eta$ [:,2] |
| set_saddles | ( | self, | |
| saddles, | |||
weights = None |
|||
| ) |
Set saddle constraint points (poloidal field should vanish at each point)
| saddles | List of points defining constraints [:,2] |
| weights | Weight to be applied to each constraint point [:] (default: 1) |
| set_targets | ( | self, | |
Ip = None, |
|||
Ip_ratio = None, |
|||
pax = None, |
|||
estore = None, |
|||
R0 = None, |
|||
V0 = None, |
|||
retain_previous = False |
|||
| ) |
Set global target values.
retain_previous=True.| alam | Scale factor for \(F*F'\) term (disabled if Ip is set) |
| pnorm | Scale factor for \(P'\) term (disabled if pax, estore, or R0 are set) |
| Ip | Target plasma current [A] (disabled if OFT_env.float_disable_flag) |
| Ip_ratio | Amplitude of net plasma current contribution from FF' compared to P' (disabled if OFT_env.float_disable_flag) |
| pax | Target axis pressure [Pa] (disabled if OFT_env.float_disable_flag or if estore is set) |
| estore | Target sotred energy [J] (disabled if OFT_env.float_disable_flag) |
| R0 | Target major radius for magnetic axis (disabled if OFT_env.float_disable_flag or if pax or estore are set) |
| V0 | Target vertical position for magnetic axis (disabled if OFT_env.float_disable_flag) |
| retain_previous | Keep previously set targets unless explicitly updated? (default: False) |
| setup | ( | self, | |
order = 2, |
|||
F0 = 0.0, |
|||
full_domain = False |
|||
| ) |
Setup G-S solver.
| order | Order of FE representation to use |
| F0 | Vacuum \(F(\psi)\) value (B0*R0) |
| 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:
| r | Mesh point list [np,2] |
| lc | Mesh cell list [nc,3] (base one) |
| reg | Mesh region list [nc] (base one) |
| mesh_file | Filename containing mesh to load (native format only) |
| setup_regions | ( | self, | |
cond_dict = {}, |
|||
coil_dict = {} |
|||
| ) |
Define mesh regions (coils and conductors)
| cond_dict | Dictionary specifying conducting regions |
| setup_td | ( | self, | |
| dt, | |||
| lin_tol, | |||
| nl_tol, | |||
pre_plasma = False |
|||
| ) |
Setup the time-dependent G-S solver.
| dt | Starting time step |
| lin_tol | Tolerance for linear solver |
| nl_tol | Tolerance for non-linear solver |
| pre_plasma | Use plasma contributions in preconditioner (default: False) |
| 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.
| time | Growth rate enhancement point (should be approximately expected value) |
| dt | Number of eigenvalues to compute |
| trace_surf | ( | self, | |
| psi | |||
| ) |
Trace surface for a given poloidal flux.
| psi | Flux surface to trace \(\hat{\psi}\) |
| 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.
| psi | Boundary values for vacuum solve |
| rhs_source | Current source (optional) |
|
protected |
F*F' normalization value [1] (use alam property)
|
protected |
Coil definition dictionary.
|
protected |
Conductor definition dictionary.
|
protected |
Diverted (limited) flag [1] (use diverted property)
|
protected |
Stored energy target value (use set_targets)
|
protected |
Vacuum F value.
|
protected |
Flux constraint points (use set_flux)
|
protected |
Plasma current target ratio I_p(FF') / I_p(P') (use set_targets)
|
protected |
Plasma current target value (use set_targets)
|
protected |
Isoflux constraint points (use set_isoflux)
|
protected |
Internal mesh object.
|
protected |
|
protected |
Axis pressure target value (use set_targets)
|
protected |
Pressure normalization value [1] (use pnorm property)
|
protected |
R0 target value (use set_targets)
|
protected |
Saddle constraint points (use set_saddles)
|
protected |
Internal Grad-Shafranov object (gs_eq)
|
protected |
V0 target value (use set_targets)
|
protected |
Vacuum definition dictionary.
|
protected |
Virtual coils, if present (currently only ‘’#VSC'`)
| 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].