Python API
a5py
Package for processing ASCOT5 data and generating inputs.
Interface for accessing data in ASCOT5 HDF5 files. |
|
Python interface to libascot.so. |
|
Package for generating inputs from templates and imported data. |
|
Tools for assessing physical quantities and units. |
|
Unit tests, physics tests and regressions runs. |
|
Advanced pre- and postprocessing routines. |
|
Graphical user interface for ASCOT5. |
Important classes
Primary tool for processing ASCOT5 data. |
|
Entry node for accessing data in the HDF5 file. |
|
Node containing results and methods to process them. |
|
Node containing AFSI results and methods to process them. |
|
Node containing BBNBI results and methods to process them. |
|
ASCOT Fusion Source Integrator AFSI. |
|
Biot-Savart solver for computing magnetic field from a coil geometry. |
|
Tool for processing distributions into markers. |
Ascot
- class a5py.Ascot(inputfile=None, create=False, mute='err', mpirank=0, mpisize=1, condor=False)
Primary tool for processing ASCOT5 data.
This object is supposed to be all you need to work with ASCOT5 data unless one is doing something exotic.
The data is accessed via
dataattribute and its attributes that form a treeview. For example, active magnetic field is accessed asAscot.data.bfield.active. SeeAscot5IOfor details.The active run is accessed as Ascot5.data.active and its methods contain most of the plotting and data access routines.
The methods in this class provide interface to ASCOT5 C-routines and allow simulating markers in Python or evaluating data exactly as it is evaluated run-time.
- Attributes:
- data
Ascot5IO Container for the HDF5 data.
- biosaw
BioSaw Tool for calculating magnetic field from coils.
- afsi
Afsi Tool for calculating fusion source from thermal plasma and fast ion distributions.
- markergen
MarkerGenerator Tool for generating markers from distributions.
- data
- DUMMY_QID = b''
Flag to use in _sim.qid_* to indicate that the input is not initialized.
- file_getpath()
Return name of the HDF5 file from which this instance reads the data.
- Returns:
- path
str Absolute path or None if this is an empty instance.
- path
- file_load(inputfile)
Open a new HDF5 file and free resources from the old one.
- Parameters:
- inputfile
str Path to the new HDF5 file.
- inputfile
- get_plasmaquantities()
Return species present in plasma input.
- input_check_bfield_seam(qnt, r=None, z=None, axes=None, centremostpercentage=10, diff=True, forceplot=False, n_centremost_elements_for_fit=None, plot_spline_interpolated=False)
Tries to check if an (R,z) plane is missing near the “seam” of the bfield input.
This function is intended to be used to verify that the correct number of (R,z) planes is provided when using a stellarator bfield input. To recap, psi and Bfield components should be given at the points np.linspace(phimin, phimax, nphi+1)[:-1]. (It is possible to independently specify phimin, phimax and nphi for both psi and the bfield components.)
If an (R,z) plane might be missing, plots qnt as a function of phi, evaluated only at the grid points where the input is given.
- Parameters:
- qnt
str One of the following: {“psi”, “br”, “bphi”, “bz”}
- r
float,optional R coordinate where data is evaluated. (The functions uses the R value of the input that is closest to this requested one) [m]
- z
float,optional z coordinate where data is evaluated. (The functions uses the z value of the input that is closest to this requested one) [m]
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- centremostpercentage
float,optional Percentage of the grid closest to the requested r and z. (Default is 10, i.e. 5 % of the grid on either side)
- diffbool,
optional If True, plots np.diff(qnt)/np.diff(phi) instead of qnt. A possible kink may be more clearly visible in the diff(qnt) plot.
- forceplotbool,
optional If True, plots the requested quantity even if no (R,z) plane is deemed missing.
- n_centremost_elements_for_fit
int,optional Number of elements from the centre of the quantity array to use for least square sum fitting a polynomial. Default is None, in which case the whole array, consisting of ‘centremostpercentage’, is used.
- plot_spline_interpolatedbool,
optional If True, plots the spline interpolated version of qnt.
- qnt
- input_eval(r, phi, z, t, *qnt, grid=False)
Evaluate input quantities at given coordinates.
This method uses ASCOT5 C-routines for evaluating inputs exactly as they are evaluated during a simulation. See
input_eval_list()for a complete list of available input and derived quantities.- Parameters:
- rarray_like
R coordinates where data is evaluated.
- phiarray_like
phi coordinates where data is evaluated.
- zarray_like
z coordinates where data is evaluated.
- tarray_like
Time coordinates where data is evaluated.
- *qnt
str Names of evaluated quantities.
- gridbool,
optional Treat input coordinates as abscissae and return the evaluated quantities on a grid instead.
- Returns:
- valarray_like (n,)
or(nr,nphi,nz,nt) Quantity evaluated in each queried point.
NaN is returned at those points where the evaluation failed.
4D matrix is returned when grid=``True``.
- valarray_like (n,)
- Raises:
AssertionErrorIf required data has not been initialized.
RuntimeErrorIf evaluation in libascot.so failed.
- input_eval_atomiccoefs(ma, anum, znum, r, phi, z, t, va, reaction)
Evaluate atomic reaction rates for a given test particle.
- Parameters:
- ma
float Test particle mass
- anum
int Test particle atomic mass number.
- znum
int Test particle charge number.
- rarray_like, (n,)
R coordinates where data is evaluated [m].
- phiarray_like (n,)
phi coordinates where data is evaluated [rad].
- zarray_like (n,)
z coordinates where data is evaluated [m].
- tarray_like (n,)
Time coordinates where data is evaluated [s].
- vaarray_like (n,)
Test particle velocities where data is evaluated in each grid point.
- reaction{“ionization”, “recombination”, “charge-exchange”,
- “beamstopping”}
Reaction whose cross-section is computed.
- ma
- Returns:
- sigmavarray_like, (n,nv)
Reaction cross-section.
- Raises:
AssertionErrorIf required data has not been initialized.
- input_eval_atomicsigma(ma, anum, znum, r, phi, z, t, va, ion, reaction)
Deprecated.
- input_eval_collcoefs(ma, qa, r, phi, z, t, va, *coefs, grid=True)
Evaluate Coulomb collision coefficients for a given test particle.
Collision coefficients are evaluated by interpolating the plasma parameters on given coordinates. The coefficients are returned as a function of the test particle velocity.
- Parameters:
- ma
float Test particle mass.
- qa
float Test particle charge.
- rarray_like, (n,)
R coordinates where data is evaluated.
- phiarray_like, (n,)
phi coordinates where data is evaluated.
- zarray_like, (n,)
z coordinates where data is evaluated.
- tarray_like, (n,)
Time coordinates where data is evaluated.
- vaarray_like, (nv,)
or(n,) Test particle velocity.
- *coefs
str Names of the coefficients to be evaluated.
“clog” - Coulomb logarithm. “f” - Drag in particle Fokker-Planck equation. “k” - Drag in guiding-center Fokker-Planck equation. “nu” - Pitch collision frequency. “dpara” - Parallel diffusion. “dperp” - Perpendicular diffusion. “ddpara” - d(Dpara) / dv. “q” - From K = Q + dDpara + 2*Dpara / va. “dq” - d(Q) / dv. “mu0” - One of the special functions needed in evaluation. “mu1” - One of the special functions needed in evaluation. “dmu0” - One of the special functions needed in evaluation.
- gridbool,
optional If True, all velocity components are evaluated at each (R,phi,z) position, i.e., the returned values are of shape (nion+1,n,nv).
If False, the coefficients are evaluated at positions (R,phi,z,va), e.g. along an orbit, and the returned values have shape (nion+1,n).
- ma
- Returns:
- *outarray_like, (nion+1,
n,nv)or(nion+1,n) Evaluated collision coefficients in same order as declared in
*coefs.The first dimension is the background plasma species where the first index is for the electrons followed by ions.
- *outarray_like, (nion+1,
- Raises:
AssertionErrorIf required data has not been initialized.
RuntimeErrorIf evaluation in libascot.so failed.
- input_eval_list(show=True)
Return quantities that can be evaluated from inputs.
The returned quantities are those that can be evaluated from the inputs with input_eval method. Units are not listed as that function returns quantities with units included.
The number ion species in plasma input is not fixed so “ni” refers to i’th plasma species as listed by input_getplasmaspecies (index starts from 1).
- input_eval_orbitresonance(rhogrid, xigrid, egrid, species, plot=False, n=1, p=0, omega=0.0, axes=None, cax=None)
Evaluate the resonance condition between particles and perturbations.
Orbit resonance refers to the phenomenon where the frequency of a particle’s orbit in a magnetic field matches the frequency of a perturbation, such as a wave or an external perturbation. This resonance can lead to increased particle transport.
The resonance condition for orbit resonance is given by:
f_{tor} / f_{pol} = m / n
where: - f_{tor} is the toroidal orbit frequency of the particle, - f_{pol} is the poloidal orbit frequency of the particle, - n is the toroidal mode number of the perturbation. - m is an arbitrary integer.
This function can evaluate and plot the resonance in 1D or 2D space. Note that one of the (‘rhogrid’, ‘xigrid’, ‘egrid’) must be a scalar.
- Parameters:
- rhogridarray_like, (
nrho) Value of the radial coordinate of radial grid for which the resonance is evaluated.
- xigrid
floator array_like, (nxi) Value of the pitch of pitch grid for which the resonance is evaluated.
- egrid
floator array_like, (ne) Value of the energy of energy grid for which the resonance is evaluated.
- species
str Name of the particle species.
- plotbool,
optional If True, the results are not returned but visualized instead.
- resonance
floator array_like,optional Resonance(s) to be indicated on the plot e.g. m/n=2/1.
- omega
float,optional Perturbation mode frequency.
The default value assumes stationary perturbation.
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- cax
Axes,optional The color bar axes or otherwise taken from the main axes.
- rhogridarray_like, (
- Returns:
- torfreqarray_like (
nrho,nxi)or(nrho,ne)or(nxi,ne) Evaluated toroidal orbit frequencies if ‘plot’ = False.
- polfreqarray_like (
nrho,nxi)or(nrho,ne)or(nxi,ne) Evaluated poloidal orbit frequencies if ‘plot’ = False.
- qeffarray_like (
nrho,nxi)or(nrho,ne)or(nxi,ne) Evaluated effective qfactor = torfreq / polfreq if ‘plot’ = False.
- torfreqarray_like (
- input_eval_rfof(m, q, vpar, r, phi, z, t)
Evaluate Evaluate ICRH electric field and the resonance condition.
The evaluated electric field consists of left-hand (-) and right-hand (+) circularly polarized components. The resonance condition is given by
omega_wave - n * omega_gyro - k_parallel * v_parallel - k_perp dot v_drift = 0.
- Parameters:
- m
float Test particle mass (for calculating the resonance).
- q
float Test particle charge (for calculating the resonance).
- vpar
float Test particle parallel velocity (for calculating the resonance).
- rarray_like, (n,)
R coordinates where data is evaluated.
- phiarray_like (n,)
phi coordinates where data is evaluated.
- zarray_like (n,)
z coordinates where data is evaluated.
- tarray_like (n,)
Time coordinates where data is evaluated.
- m
- Returns:
- eplus_realarray_like, (n,)
- eminus_realarray_like, (n,)
- eplus_imagarray_like, (n,)
- eminus_imagarray_like, (n,)
- rescondarray_like, (n,)
- Raises:
AssertionErrorIf required data has not been initialized.
- input_eval_ripple(rgrid, zgrid, rlarmor, rho0=unyt_quantity(1., 'dimensionless'), theta0=unyt_quantity(0, 'rad'), t=unyt_quantity(0, 's'), nphi=362, plot=True, axes1=None, axes2=None)
Evaluate various ripple quantities.
The variation of toroidal field strength as a function of toroidal angle is called the toroidal field (TF) ripple. The TF ripple can lead to ripple-induced losses provided that i) the particle is trapped in a ripple well in which case it is promptly lost or ii) the ripple causes particle motion to become stochastic which leads to radial diffusion.
This function evaluates the quantities relevant for ripple transport. Note that
deltacritandripplewellare approximate quantities based on analytical estimates.- Parameters:
- rgridarray_like
R points where
delta,deltacrit, andripplewellare evaluated.- zgridarray_like
z points where
delta,deltacrit, andripplewellare evaluated.- rlarmor
floatortuple(str,int,float) Larmor radius of the test particle or tuple with (species, charge, energy).
- rho0
float,optional Radial coordinate where
amplitudeis evaluated.- theta0
float,optional Poloidal coordinate where
amplitudeis evaluated.- t
float,optional Time instance for time-dependent fields.
- nphi
int,optional Number of toroidal grid points in ripple evaluation.
There is probably never a need to change this number.
- plotbool,
optional If True, the evaluated quantities are also plotted.
- axes1,
optional Axes for plotting
amplitudeifplot= True.- axes2,
optional Axes for plotting
delta,deltacrit, andripplewellifplot= True.
- Returns:
- phigridarray_like, (
nphi) Toroidal grid [-2, 362] where
amplitudeis given.This grid is not limited to [0,360] as usual so that one can also use it to spot any discontinuities in field if the input was poorly constructed.
- amplitudearray_like, (nphi,)
Variation in toroidal field strength, i.e. Bphi(phi) - Bphi_axisymmetric.
- deltaarray_like, (nr,nz)
Ripple magnitude defined as (Bmax-Bmin) / (Bmax+Bmin), where Bmax and Bmin are toroidal extrema of Bphi.
- deltacritarray_like, (nr,nz)
Critical ripple magnitude divided by
delta.Particles are subject to stochastic ripple diffusion when
deltacrit< 1.- ripplewellarray_like, (nr,nz)
Ripple well parameter.
Particles are bound to be ripple-trapped when
ripplewell< 1. However, this does not necessarily mean that the particle is lost as it may escape the well when ripple-trapping causes it to drift outwards. As particles that are ripple-trapped drift vertically, a good indication that a trapped particle will be lost if the ripple well region covers everything vertically below/above that point.
- phigridarray_like, (
- input_eval_safetyfactor(rho, nth=10000, return_ftrap=False)
Evaluate safety factor and associated quantities.
This function most likely works for tokamaks only.
- Parameters:
- rhoarray_like, (n,)
Radial positions where the safety factor is evaluated.
Note that evaluating it near the axis or exactly at the separatrix may yield mess.
- return_bminmaxbool,
optional Return also minimum and maximum values of B along the flux surface.
- Returns:
- qprofarray_like, (n,)
Safety factor.
- Iprofarray_like, (n,)
Toroidal current term which, when multiplied with 2pi/mu0, gives the enclosed toroidal current.
- gprofarray_like, (n,)
This is just R * Bphi which is constant since Bphi ~ 1/R.
- ftraparray_like, (n,)
Trapped particle fraction if return_ftrap is true.
- input_findpsi0(psi1, phi=unyt_quantity(0., 'rad'), nphi=None, phimin=None, phimax=None, manual_pad=1e-09)
Find poloidal flux on axis value numerically.
Before this function is called, the magnetic field data should contain initial guess for the position of the magnetic axis. The algorithm then uses the gradient descent method to find psi0 and its location (either 2D or 3D). The interpolation is done using Ascot’s magnetic field interpolation and a little bit of padding is added to psi0 so the value can be used as an input parameter for the magnetic field without any errors.
- Parameters:
- psi1
float Poloidal flux at the separatrix.
This value is used to deduce whether the algorithm searches minimum or maximum value when finding psi0.
- phi
float,optional Can be used to manually set the phi co-ordinate when using the 2D search. Zero by default.
- nphi
int,optional Number of B field grid points in the phi direction between phimin and phimax including the end points of this interval.
Needed for 3D fields. For clarity: if the first and the last grid point in the interval are the same point, this point is counted twice. Example: you have three grid points at 0 deg, 120 deg and 240 deg and you use phimin=0, phimax=2pi. The input argument nphi is then four instead of three because at 2pi you count again the first grid point.
- phimin
float,optional Minimum of the phi interval.
Needed for 3D fields.
- phimax
float,optional Maximum of the phi interval.
Needed for 3D fields.
- manual_pad
float,optional The magnitude of manual paddin that is applied to psi0 afterwards. None by default but can be used when the default padding has not been enough. This could manifest itself, e.g., as errors when converting (rho,phi,theta) to (R,z) near the magnetic axis.
- psi1
- Returns:
- r
floator array_like (nphi-1,) Axis R-coordinate. If nphi given, returns an array.
- phiarray_like (nphi-1,),
optional Axis phi-coordinate. Returned only if nphi is given.
- z
floator array_like (nphi-1,) Axis z-coordinate. If nphi given, returns an array.
- psi0
floator array_like (nphi-1,) Poloidal flux on axis. If nphi given, returns an array.
- r
- Raises:
AssertionErrorIf required data has not been initialized.
RuntimeErrorIf evaluation in libascot.so failed.
Notes
In 2D geometry (no nphi specified), returns (r, z). In 3D geometry (when nhpi specified), returns (r, phi, z). These points are generally speaking not uniformly spaced.
- input_free(bfield=False, efield=False, plasma=False, wall=False, neutral=False, boozer=False, mhd=False, asigma=False)
Free input used by the Python interface.
Arguments toggle which input fields are free’d. If called without providing any arguments, then all input fields are free’d.
- Parameters:
- bfieldbool,
optional Toggle to free the magnetic field.
- efieldbool,
optional Toggle to free the electric field input.
- plasmabool,
optional Toggle to free the plasma input.
- wallbool,
optional Toggle to free the wall input.
- neutralbool,
optional Toggle to free the neutral input.
- boozerbool,
optional Toggle to free the boozer input.
- mhdbool,
optional Toggle to free the MHD input.
- asigmabool,
optional Toggle to free the atomic data input.
- bfieldbool,
- input_free_rfof()
- input_getplasmaspecies()
Get species present in plasma input (electrons first).
- Returns:
- nspecies
int Number of species.
- massarray_like
Species’ mass.
- chargearray_like
Species’ charge state.
- anumarray_like
Species’ atomic mass number.
- znumarray_like
Species’ charge number.
- nspecies
- Raises:
AssertionErrorIf required data has not been initialized.
RuntimeErrorIf evaluation in libascot.so failed.
- input_init(run=False, bfield=False, efield=False, plasma=False, wall=False, neutral=False, boozer=False, mhd=False, asigma=False, switch=False)
Initialize input so it can be accessed via Python interface.
This method can be used in three ways:
Without any arguments in which case all inputs are initialized. The initialized input is the one that has been set as active.
Provide one of the input arguments (bfield, efield, plasma, wall, neutral, boozer, mhd) in which case only that input is initialized. The argument can either be QID of the input to be initialized, or True in which case the active input of that type is initialized.
Set run = True or QID of a run group, in which case the inputs of that run group are initialized. If you only want a specific input to be initialized, set that argument as True, e.g. bfield = True would only initialize the magnetic field used by the given run.
If the input is already initialized, nothing is done. In case there is already different input of same type initialized, an error is raised.
If the input argument is a dict, it is used instead of reading the data from hdf5.
- Parameters:
- run
stror bool,optional The run group or True to use the one that is active.
- bfield
stror bool ordict,optional The magnetic field input or True to use the one that is active.
- efield
stror bool ordict,optional The electric field input or True to use the one that is active.
- plasma
stror bool ordict,optional The plasma input or True to use the one that is active.
- wall
stror bool ordict,optional The wall input or True to use the one that is active.
- neutral
stror bool ordict,optional The neutral input or True to use the one that is active.
- boozer
stror bool ordict,optional The boozer input or True to use the one that is active.
- mhd
stror bool ordict,optional The MHD input or True to use the one that is active.
- asigma
stror bool ordict,optional The atomic data input or True to use the one that is active.
- switchbool,
optional If
True, no error is raised if there is already a different input initialized and instead the old input is free’d and the new one initialized.
- run
- input_init_rfof()
- input_initialized()
Get inputs that are currently initialized.
- input_plotalongaxis(qnt, nphi=None, phimin=0, phimax=360, t=unyt_quantity(0., 's'), axes=None)
Plot qnt on the magnetic axis as a function of phi.
By default plots qnt in the phigrid of the axis inputs. However, if nphi is given, qnt is plotted against np.linspace(phimin, phimax, nphi)[:-1]
- input_plotcontoursurface3d(qnt, isovalue, xmin, xmax, ymin, ymax, zmin, zmax, pointspermeterx=10, pointspermetery=10, pointspermeterz=10, t=unyt_quantity(0, 's'), max_rho=1.0, ax=None, c=0.135, delta=None, tri_lc=None, tri_lw=0.4, opacity=0.8, reltol=0.01, axscatter=None, surfcolor='blue')
Plot a 3D contour surface of the specified quantity at a given isovalue.
If there are holes in you surface, try lowering input parameter c or increasing reltol.
- Parameters:
- qnt
str Name of the quantity whose contour surface is drawn (e.g., “bnorm”).
- isovalue
float Value of the quantity at which to generate the isosurface.
- xmin, xmax
float Minimum and maximum bounds along the x-axis.
- ymin, ymax
float Minimum and maximum bounds along the y-axis.
- zmin, zmax
float Minimum and maximum bounds along the z-axis.
- pointspermeterx, pointspermetery, pointspermeterz
int,optional Number of grid points along the x, y, and z directions per meter. Ideally, use uniform grid spacing in all directions.
- t
floatorunyt_quantity,optional Time coordinate at which the data is evaluated. Default is 0 seconds.
- max_rho
float,optional Maximum rho until (surface outside this is not drawn)
- ax
matplotlib.axes.Axes,optional Matplotlib 3D axes to plot on. If None, a new figure and axes are created.
- c
float,optional Scaling parameter for the grid spacing or contour thresholding. Default is 0.135.
- delta
floatorNone,optional Characteristic grid spacing. Ideally, use uniform grid of points.
- tri_lc
strorNone,optional Color of the mesh triangle edges. If None, edges are invisible.
- tri_lw
float,optional Width of the mesh triangle edges. Default is 0.4.
- opacity
float,optional Opacity of the drawn surface. (0: transparent, 1: opaque)
- reltol
float,optional Relative tolerance for isosurface extraction. Default is 0.01.
- axscatter3D
axes,optional If given, plots the grid points of the surface on these axes. This can be useful especially when verifying that the the plotted surface is sensible.
- qnt
- input_plotfluxsurface3d(rho, phi, theta, qnt=None, t=unyt_quantity(0, 's'), clim=None, cmap=None, axes=None, cax=None, meshalpha=0.6, tri_lc=None, tri_lw=0.4, opacity=0.8)
Plot a flux surface or, given qnt, a quantity on a flux surface in 3D.
To plot quantity on a logarithmic scale (base 10), add “log” in the name of the quantity e.g. “log ne”.
- Parameters:
- rho
float Rho co-ordinate of the flux surface.
- phiarray_like (nphi,)
phi abscissa where data is evaluated and plotted.
- thetaarray_like (ntheta,)
theta abscissa where data is evaluated and plotted.
- qnt
str,optional Name of the plotted quantity. E.g. “bnorm”
- t
float,optional Time coordinate where data is evaluated.
- clim
list[float,float],optional Minimum and maximum color range values.
- cmap
str,optional Colormap.
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- cax
Axes,optional The color bar axes or otherwise taken from the main axes.
- meshalpha: float, optional
A parameter used when generating the triangles for the surface. For more info, see “alphashape” and “alpha parameter”.
- tri_lc
str,optional Color of the mesh triangle edges. By default, edges are invisible.
- tri_lw
float,optional Width of the mesh triangle edges. To enable edges, set also tri_lc.
- opacity
float,optional Opacity of the drawn surface (0: transparent, 1: opaque).
- rho
- input_plotphitheta(phi, theta, rho, qnt, t=unyt_quantity(0, 's'), clim=None, cmap=None, axes=None, cax=None, drawcontours=False, contourvalues=None, contourcolors=None, contourlinestyles='-', contourlinewidths=None, labelfontsize=10, drawbfielddir=False, bfieldlinecolor='black', bfieldlinedensity=0.5, blinewidth=1.0, blinelength=50.0)
Plot input quantity on a (phi, theta) surface at given rho.
With drawcontours=True, draws also the contourlines. To plot quantity on a logarithmic scale (base 10), add “log” in the name of the quantity e.g. “log ne”.
- Parameters:
- phiarray_like (nphi,)
phi abscissa where data is evaluated and plotted.
- thetaarray_like (ntheta,)
theta abscissa where data is evaluated and plotted.
- rho
float Rho co-ordinate of the flux surface.
- qnt
str Name of the plotted quantity.
- t
float,optional Time coordinate where data is evaluated.
- clim
list[float,float],optional Minimum and maximum color range values.
- cmap
str,optional Colormap.
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- cax
Axes,optional The color bar axes or otherwise taken from the main axes.
- drawcontoursbool,
optional If true, draws contourlines on top of the heatmap.
- contourvalues
list[float],optional Which contour lines to draw. Leave empty for automatic selection
- contourcolors
strorlist[str],optional Colors of the contour lines. By giving a list, one can match the colors with the contourvalues (see above).
- contourlinestyles
strorlist[str],optional Linestyle of the contour lines – as the reader might have guessed. By giving a list, one can match the styles with the contourvalues ( see above).
- contourlinewidths
floatorlist[float],optional You can probably guess. By giving a list, one can match the widths with the contourvalues (see above).
- labelfontsize
float,optional Font size of the labels that indicate contour values.
- drawbfielddirbool,
optional Draw the direction of the bfield on top of the figure as streamlines. This makes sense only in a tokamak geometry since the poloidal plane is a contant phi plane. If drawfielddir true, one must have theta and phi go from 0 to 360.
- bfieldlinecolor
str,optional Color of the bfield lines if drawbfielddir=True.
- bfieldlinedensity
float,optional Density of the bfield lines.
- blinewidth
float,optional Width of the bfield lines.
- blinelength
float,optional Length of the bfield lines.
- input_plotrhocontour(rho=1.0, phi=unyt_quantity(0, 'degree'), t=unyt_quantity(0, 's'), ntheta=180, axes=None, **kwargs)
Plot contours of constant rho.
Calling this without arguments plots the separatrix.
- Parameters:
- rho
floator array_like,optional Value(s) of square root of normalized poloidal flux whose contours are plotted.
- phi
float,optional Toroidal angle at which contours are plotted.
- t
float,optional Time instance when the contours are plotted.
- ntheta
int,optional Poloidal resolution of the contours.
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- **kwargs
Arguments passed to
matplotlib.pyplot.plot().
- rho
- input_plotrz(r, z, qnt, phi=unyt_quantity(0, 'degree'), t=unyt_quantity(0, 's'), rhomax=None, clim=None, cmap=None, axes=None, cax=None, drawcontours=False, contourvalues=None, contourcolors=None, contourlinestyles='-', contourlinewidths=None, labelfontsize=10)
Plot input quantity on a (R, z) plane at given coordinates.
With drawcontours=True, draws also the contourlines. To plot quantity on a logarithmic scale (base 10), add “log” in the name of the quantity e.g. “log ne”.
- Parameters:
- rarray_like (nr,)
R abscissa where data is evaluated and plotted.
- zarray_like (nz,)
z abscissa where data is evaluated and plotted.
- qnt
str Name of the plotted quantity.
- phi
float,optional Toroidal coordinate where data is evaluated.
- t
float,optional Time coordinate where data is evaluated.
- rhomax
float,optional Do not plot (R, z) points for which rho > rhomax
- clim
list[float,float],optional Minimum and maximum color range values.
- cmap
str,optional Colormap.
- axes
Axes,optional The axes where figure is plotted or otherwise new figure is created.
- cax
Axes,optional The color bar axes or otherwise taken from the main axes.
- drawcontoursbool,
optional If true, draws contourlines on top of the heatmap.
- contourvalues
list[float],optional Which contour lines to draw. Leave empty for automatic selection
- contourcolors
strorlist[str],optional Colors of the contour lines. By giving a list, one can match the colors with the contourvalues (see above).
- contourlinestyles
strorlist[str],optional Linestyle of the contour lines – as the reader might have guessed. By giving a list, one can match the styles with the contourvalues ( see above).
- contourlinewidths
floatorlist[float],optional You can probably guess. By giving a list, one can match the widths with the contourvalues (see above).
- labelfontsize
float,optional Font size of the labels that indicate contour values.
- input_plotwallcontour(phi=unyt_quantity(0, 'degree'), axes=None, color='black', label=None)
Plot intersection of the wall and the poloidal plane at the given toroidal angle.
- input_recalculate_axis(psi1, nphi, phimin, phimax)
Calculates the (R, phi, z) co-ordinates of the magnetic axis.
Before this function is called, the magnetic field data should contain initial guess for the position of the magnetic axis. (see input_findpsi0)
- Parameters:
- psi1
float Poloidal flux at the separatrix.
This value is used to deduce whether the algorithm searches minimum or maximum value of psi when finding the magnetic axis.
- nphi
float Number of grid points in the phi direction.
- phimin
float Minimum of the phi co-ordinate
- phimax
float Maximum of the phi co-ordinate
- psi1
- Returns:
- rarray_like (nphi,)
Axis R-coordinate.
- phiarray_like (nphi,)
Axis phi-coordinate.
- zarray_like (nphi,)
Axis z-coordinate.
- input_rhotheta2rz(rho, theta, phi, time, maxiter=100, tol=1e-05)
Convert (rho, theta, phi) coordinates to (R,z) positions.
- Parameters:
- rhovalsarray_like (n,)
Normalized poloidal flux coordinates to be converted.
- thetaarray_like (n,)
Poloidal angle coordinates to be converted.
- phiarray_like (n,)
Toroidal angle coordinates to be converted.
- time
float Time slice (same for all).
- Returns:
- rarray_like (n,)
R coordinates at (rhovals, theta, phi)
- zarray_like (n,)
z coordinates at (rhovals, theta, phi)
- Raises:
AssertionErrorIf required data has not been initialized.
RuntimeErrorIf evaluation in libascot.so failed.
- input_rhovolume(nrho=10, ntheta=360, nphi=10, method='mc', tol=0.1, t=unyt_quantity(0, 's'), return_area=False, return_coords=False, minrho=0, maxrho=1, minphi=0, maxphi=360, mintheta=0, maxtheta=360)
Evaluate flux surface volumes.
- Parameters:
- nrho
int,optional Number of radial grid edges between [0, 1].
- ntheta
int,optional Number of poloidal grid edges.
- nphi
int,optional Number of toroidal grid edges.
- method{“prism”, “mc”},
optional Method used to evaluate the volumes.
“prism” divides the volume into structured 3D grid consisting of triangular prims where the triangular plane is poloidal and one of the vertices is at magnetic axis while the other two are along a given flux surface. The volume is the estimated as a sum of these prisms.
“mc” divides the volume into toroidal sectors of same size, and on each sector fins the separatrix bounding box. Bounding box is assumed to be fixed along the sector and markers are randomly sampled within it. The volumes are then estimated via the Monte Carlo method.
- tol
float,optional If relative difference in volume in subsequent iterations falls below this value, the algorithm finishes.
- t
float,optional Time slice when the volumes are computed.
- return_areabool,
optional Return also the area on (R, z) plane.
- return_coordsbool,
optional Return also the (R, phi, z) coordinates for grid center points.
- minrho
float,optional minimum value of rho
- maxrho
float,optional maximum value of rho
- minphi
float,optional minimum value of phi (in degrees)
- maxphi
float,optional maximum value of phi (in degrees)
- mintheta
float,optional minimum value of theta (in degrees)
- maxtheta
float,optional maximum value of theta (in degrees)
- nrho
- Returns:
- volumearray_like (nrho-1, ntheta-1, nphi-1)
Volume of each bin.
- areaarray_like (nrho-1, ntheta-1, nphi-1)
Area of the poloidal cross section of each bin if
return_area=True.- rarray_like (nrho-1, ntheta-1, nphi-1)
R coordinates of bin centers if
return_coords=True.- parray_like (nrho-1, ntheta-1, nphi-1)
phi coordinates of bin centers if
return_coords=True.- zarray_like (nrho-1, ntheta-1, nphi-1)
z coordinates of bin centers if
return_coords=True.
- preflight_bfieldpsi0()
Checks whether psi0 given in input is actually extreme value.
Because psi is interpolated with splines, there might be numerical error that causes psi (near the axis) to have more extreme value than psi0 which is given in input. This leads to imaginary rho and termination of the simulation if marker ends up there.
This check uses Monte Carlo method to 1. Draw phi 2. Evaluate axis (R,z) 3. Draw random (R,z) coordinates within 10 cm of the axis. 4. Evaluate psi at that point and compare to psi0. Process repeats N times. Check passes if all evaluations are valid.
- preflight_inputspresent()
Check required inputs are present for this run.
- preflight_optionsconsistent()
Check that options are consistent with inputs and each other.
- preflight_plottopview(hidewall=False, hidemarkers=False, axes=None)
Plot top view of the machine showing Ip, Bphi, and possibly markers and wall if present.
Assumes bfield is initialized in ascotpy.
- simulation_bbnbi(nprt, t1=0, t2=0, printsummary=True)
Run BBNBI simulation.
- Parameters:
- nprt
Numberofmarkerstobelaunchedintotal. This number is the combined total for all injectors which are then distributed among the injectors in proportion to the injector power. For example, suppose there are 3 injectors with P_1 = 10 MW and P_2 = 5 MW and P_3 = 5 MW. If nprt = 10000, then 5000 markers are generated for injector 1 and 2500 markers for both injectors 2 and 3.
The number of ionized markers is either equal or smaller than nprt as some of the markers are lost as shinethrough.
- t1
float,optional The time instant at which the injector is turned on.
In the usual case where the background is time independent, the parameters t1 and t2 only affect the marker initial time which will be uniformly distributed in [t1, t2].
- printsummarybool,
optional If True, summary of marker endstates is printed after the simulation completes.
- nprt
- Returns:
- run
VirtualBBNBI An object that acts almost exactly as
BBNBIGroupexcept that the data is read from C arrays in the memory instead of HDF5 file.The run can be used until the output is freed with
simulation_free(). Previous data must be freed before rerunning the simulation.
- run
- Raises:
AscotInitExceptionIf inputs are not packed, markers are not initialized or previous results have not been freed.
- simulation_free(inputs=False, markers=False, diagnostics=False)
Free resources used by the interactive simulation.
If called without arguments, everything will be freed.
- simulation_initbbnbi(bfield=True, plasma=True, neutral=True, wall=True, asigma=True, nbi=True, switch=True)
Initialize inputs for BBNBI simulation.
This method must be called before running BBNBI simulation.
- simulation_initinputs(bfield=True, efield=True, plasma=True, neutral=True, wall=True, boozer=True, mhd=True, asigma=True, switch=True)
Prepare input fields for the interactive simulation.
Initializes simulation inputs. The inputs used in the simulation are those that are active.
This method differs from
input_init()in that here the inputs are “packed” internally in a single (offload) array as this is what ASCOT5 does. Inputs cannot be changed while the they are packed. The unpacking is done withsimulation_free().This method must be called before running the simulation.
- simulation_initmarkers(**mrk)
Create markers for the interactive simulations.
Any existing markers are deallocated when new ones are created.
- Parameters:
- **mrk
Marker input (all marker types are supported) as accepted by
Prt.write_hdf5(),GC.write_hdf5(), orFL.write_hdf5().If not given, markers are read from the HDF5 file.
- simulation_initoptions(**opt)
Set simulation options for the interactive simulation.
There is no need to free options before setting new ones.
- Parameters:
- **opt
Options in same format as accepted by
Opt.write_hdf5().If not given, the options are read from the HDF5 file.
- simulation_run(printsummary=True)
Run the interactive simulation using inputs, options and markers that were set.
- Parameters:
- printsummarybool,
optional If True, summary of marker endstates is printed after the simulation completes.
- printsummarybool,
- Returns:
- run
VirtualRun An object that acts almost exactly as
RunGroupexcept that the data is read from C arrays in the memory instead of HDF5 file.The run can be used until the output is freed with
simulation_free(). Previous data must be freed before rerunning the simulation.
- run
- Raises:
AscotInitExceptionIf inputs are not packed, markers are not initialized or previous results have not been freed.
- class a5py.routines.afsi5.Afsi(ascot)
ASCOT Fusion Source Integrator AFSI.
ASCOT Fusion Source Integrator (AFSI) is a tool for calculating fusion rates and fusion products for arbitrary particle populations. It can be used either as a preprocessing tool e.g. to generate source distribution of alphas or as a postprocessing tool e.g. to evaluate neutron source from beam-thermal fusion.
AFSI can be used in three ways:
thermal: Calculates fusion rates between two Maxwellian populations.
beam-thermal: Calculates fusion rates between a Maxwellian and arbitrary population e.g. beam ions as obtained from ASCOT5 simulation.
beam-beam: Calculates fusion rates between two arbitrary populations.
Reference: https://iopscience.iop.org/article/10.1088/1741-4326/aa92e9
- Attributes:
- _ascot
Ascot Ascot object used to run AFSI.
- _ascot
- beambeam(reaction, beam1, beam2=None, nmc=1000, ppar1=None, pperp1=None, ppar2=None, pperp2=None)
Calculate beam-beam fusion.
- Parameters:
- reaction
int Fusion reaction index.
- beam1
dict Beam1 distribution.
- beam2
dict,optional Beam2 distribution or None to calculate fusion generation with beam1 itself.
- nmc
int,optional Number of MC samples used in each (R, phi, z) bin.
- ppar1array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the first reaction product.
- pperp1array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the first reaction product.
- ppar2array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the second reaction product.
- pperp2array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the second reaction product.
- reaction
- Returns:
- prod1array_like
Fusion product 1 distribution.
- prod2array_like
Fusion product 2 distribution.
- beamthermal(reaction, beam, swap=False, nmc=1000, ppar1=None, pperp1=None, ppar2=None, pperp2=None)
Calculate beam-thermal fusion.
- Parameters:
- reaction
int Fusion reaction index
- beam
dict Beam distribution that acts as the first reactant.
- swapbool,
optional If True, beam distribution acts as the second reactant and the first reactant is a background species.
- nmc
int,optional Number of MC samples used in each (R, phi, z) bin.
- ppar1array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the first reaction product.
- pperp1array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the first reaction product.
- ppar2array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the second reaction product.
- pperp2array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the second reaction product.
- reaction
- Returns:
- prod1array_like
Fusion product 1 distribution.
- prod2array_like
Fusion product 2 distribution.
- reactions(reaction=None)
Return reaction data for a given reaction.
- Parameters:
- reaction
str,optional Reaction or None to return list of available reactions.
- reaction
- Returns:
- reactions[
str] List of available reactions if
reaction=None.- m1
float Mass of reactant 1.
- q1
float Charge of reactant 1.
- m2
float Mass of reactant 2.
- q2
float Charge of reactant 2.
- mprod1
float Mass of product 1.
- qprod1
float Charge of product 1.
- mprod2
float Mass of product 2.
- qprod2
float Charge of product 2.
- q
float Energy released in the reaction.
- reactions[
- thermal(reaction, r=None, phi=None, z=None, rho=None, theta=None, ppar1=None, pperp1=None, ekin1=None, pitch1=None, ppar2=None, pperp2=None, ekin2=None, pitch2=None, nmc=1000)
Calculate thermonuclear fusion between two thermal (Maxwellian) species.
- Parameters:
- reaction
intorstr Fusion reaction index or name.
- rarray_like
Abscissa for the radial coordinate in (R,phi,z) basis.
- phiarray_like
Abscissa for the toroidal coordinate in (R,phi,z) and (rho,theta,phi) basis.
- zarray_like
Abscissa for the axial coordinate in (R,phi,z) basis.
- rhoarray_like
Abscissa for the radial coordinate in (rho,theta,phi) basis.
- thetaarray_like
Abscissa for the poloidal coordinate in (rho,theta,phi) basis.
- ppar1array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the first reaction product.
- pperp1array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the first reaction product.
- ekin1array_like
Abscissa for the kinetic energy in (E,pitch) basis for the first reaction product.
- pitch1array_like
Abscissa for the pitch in (E,pitch) basis for the first reaction product.
- ppar2array_like
Abscissa for the parallel momentum in (ppar,pperp) basis for the second reaction product.
- pperp2array_like
Abscissa for the perpendicular momentum in (ppar,pperp) basis for the second reaction product.
- ekin2array_like
Abscissa for the kinetic energy in (E,pitch) basis for the second reaction product.
- pitch2array_like
Abscissa for the pitch in (E,pitch) basis for the second reaction product.
- nmc
int,optional Number of MC samples used in each (R, phi, z) bin.
- reaction
- Returns:
- prod1array_like
Source distribution of the first fusion product.
- prod2array_like
Source distribution of the second fusion product.
- class a5py.routines.biosaw5.BioSaw(ascot)
Biot-Savart solver for computing magnetic field from a coil geometry.
- Attributes:
- _ascot
Ascot Ascot object used to run BioSaw.
- _ascot
- addto2d(coilxyz, phimin=0, phimax=360, nphi=360, rmin=None, rmax=None, nr=None, zmin=None, zmax=None, nz=None, current=1.0, revolve=None, exclude_bphi_from_eqdsk=False, b0=None, b2d=None)
- Generate 3D tokamak magnetic field input from coil geometry and
2D input.
This function uses the active
B_2DSinput to generate theB_3DSinput.Note that the toroidal grid is defined as:
np.linspace(phimin, phimax, nphi+1)[:-1]- Parameters:
- coilxyz
listor array_like, (ncoil,3) xyz coordinates defining the coil segments.
This argument can also be a list of coordinate arrays, in which case each list element is assumed to be an individual coil and the resulting field is summed together.
- phimin
float Toroidal grid minimum value.
- phimax
float Toroidal grid maximum value.
- nphi
int Number of toroidal grid points.
- rmin
float,optional Minimum value in R grid unless taken from the input data.
- rmax
float,optional Maximum value in R grid unless taken from the input data.
- nr
int,optional Number of R grid points unless taken from the input data.
- zmin
float,optional Minimum value in z grid unless taken from the input data.
- zmax
float,optional Maximum value in z grid unless taken from the input data.
- nz
int,optional Number of z grid points unless taken from the input data.
- current
listorfloat Coil current in Ampères either same for all coils or given individually.
Note that the absolute value does not matter here (only relative values between the coils) if
b0is provided.- revolve
int,optional If set, the field is rotated repeatedly in toroidal direction by a given amount of phi indices and summed, as if the coils were revolved.
In other words, the total field is given by:
B_tot = B_coil(phi0 at 0) + B_coil(phi0 at phi[revolve]) + B_coil(phi0 at phi[2*revolve]) + …,
where the field is rotated through the whole
phiarray.- exclude_bphi_from_eqdskbool,
optional Exclude diamagnetic contribution to Bphi that comes from the equilibrium.
If this flag is set to True, the toroidal magnetic field is perfectly 1/R, but for realistic simulations this flag should be False.
- b0{“axis”}
orfloat,optional If given, the current in the coils is scaled (after revolving them if requested) so that the toroidal field on axis has requested value on average.
If this parameter is float, the toroidal field will have this value on axis on average. If this parameter is “axis”, the average value is the same as what EQDSK has.
- b2d
dict,optional B_2DSdata to be used or otherwise read from active bfield input.
- coilxyz
- Returns:
- addto3d(coilxyz, current=1.0, revolve=None, b3d=None)
- Generate new 3D tokamak magnetic field input where contribution from
given coils are included.
- Parameters:
- coilxyz
listor array_like, (ncoil,3) xyz coordinates defining the coil segments.
This argument can also be a list of coordinate arrays, in which case each list element is assumed to be an individual coil and the resulting field is summed together.
- current
listorfloat Coil current in Ampères either same for all coils or given individually.
- revolve
int,optional If set, the field is rotated repeatedly in toroidal direction by a given amount of phi indices and summed, as if the coils were revolved.
In other words, the total field is given by:
B_tot = B_coil(phi0 at 0) + B_coil(phi0 at phi[revolve]) + B_coil(phi0 at phi[2*revolve]) + …,
where the field is rotated through the whole
phiarray.- b3d
dict,optional B_3DSdata to be used or otherwise read from active bfield input.
- coilxyz
- Returns:
- calculate(r, phi, z, coilxyz, current=1.0, revolve=None)
Calculate magnetic field at given points due to given coil(s).
- Parameters:
- rarray_like (nr,)
R grid for the coordinates where the field is evaluated.
- phiarray_like (nphi,)
phi grid for the coordinates where the field is evaluated.
- zarray_like (nz,)
z grid for the coordinates where the field is evaluated.
- coilxyz
listor array_like, (ncoil,3) xyz coordinates defining the coil segments.
This argument can also be a list of coordinate arrays, in which case each list element is assumed to be an individual coil and the resulting field is summed together.
- current
listorfloat Coil current in Ampères either same for all coils or given individually.
- revolve
int,optional If set, the field is rotated repeatedly in toroidal direction by a given amount of phi indices and summed, as if the coils were revolved.
In other words, the total field is given by:
B_tot = B_coil(phi0 at 0) + B_coil(phi0 at phi[revolve]) + B_coil(phi0 at phi[2*revolve]) + …,
where the field is rotated through the whole
phiarray.
- Returns:
- brarray_like, (nr,nphi,nz)
Magnetic field R component.
- bphiarray_like, (nr,nphi,nz)
Magnetic field phi component.
- bzarray_like, (nr,nphi,nz)
Magnetic field z component.
- class a5py.routines.markergen.MarkerGenerator(ascot)
Tool for processing distributions into markers.
- Attributes:
- _ascot
Ascot The Ascot object this instance belongs to.
- _ascot
- generate(nmrk, mass, charge, anum, znum, particledist, markerdist=None, mode='gc', minweight=0, return_dists=False)
Generate weighted markers from marker and particle distributions.
This function takes two 5D distributions that must have identical grids. The marker distribution is a (normalized) probability distribution from which markers are sampled. The particle distribution presents physical particles, e.g. it could be an alpha particle birth rate. Once markers are sampled from the marker distribution, then each marker in a same cell is given same weight, so that those markers represent all the physical particles in that particular cell.
As markers are sampled from the marker distribution, their initial coordinates are randomly chosen from the phase space region contained within the cell.
- Parameters:
- nmrk
int Number of markers to be generated.
- mass
float Particle species` mass.
- charge
float Particle species` charge.
- anum
int Particle species` atomic mass number.
- znum
int Particle species` charge number.
- particledist
Dist 5D distribution from which marker weight is sampled.
- markerdist
Dist,optional 5D distribution from which marker initial coordinates are sampled.
If not given, the coordinates are sampled from
particledist.- mode{‘prt’, ‘gc’},
optional Decides whether the distributions (and returned markeres) are in particle or guiding center phase-space.
- minweight
float,optional Minimum weight a marker must have or else it is rejected.
- return_distsbool,
optional Return
markerdistandparticledistcalculated from the generated markers.These distributions hould converge to inputs as the marker number is increased. The exception is
particledistwhich will always get zero values in regions wheremarkerdistis zero.
- nmrk
- Returns:
- mrk
dict Marker input data on a dictionary.
- mrkdist
Dist Distribution of initialized markers which ideally should be close to the one given as an input unless markers were rejected by
minweight.- prtdist
Dist Weighted distribution of initialized markers which ideally should be close to the one given as an input, unless some phase-space regions in
markerdistwere left intentionally empty even though they would contain physical particles.
- mrk
- rhoto5d(rho, prob, r_edges, phi_edges, z_edges, mom1_edges, mom2_edges)
Maps a 1D profile to a normalized 5D dist.
The resulting probability density has only spatial variance.
- Parameters:
- rhoarray_like, (nrho,)
The radial rho grid where the profile is given.
- probarray_like, (nrho,)
Distribution value at the grid centers.
- r_edgesarray_like
R abscissa edges for the output distribution.
- phi_edgesarray_like
Phi abscissa edges for the output distribution.
- z_edgesarray_like
Z abscissa edges for the output distribution.
- mom1_edgesarray_like
Either ppar or ekin abscissa edges for the output distribution.
- mom2_edgesarray_like
Either pperp or pitch abscissa edges for the output distribution.
- Returns:
- markerdist
Dist Normalized 5D distribution from which markers can be sampled.
- markerdist
exceptions
Contains definitions of this package’s exceptions.
- exception a5py.exceptions.AscotIOException
Raised when there is an internal error in Ascot5IO.
This exception should be raised in cases where the issue is most certainly a bug in code or inconsistent HDF5 file.
- exception a5py.exceptions.AscotInitException
Raised when there is an issue with data initialized by Ascot.
This exception should be raised whenever there is an issue with libascot or any related data.
- exception a5py.exceptions.AscotNoDataException
Raised when required input or output is not present.
This exception should be raised when user tries to query input or results that is not present. Other case is when using methods that require the presence of certain input or output data (e.g. when plotting orbits but orbit data is not present).
- exception a5py.exceptions.AscotUnitWarning
Warning raised when quantities are provided without units.