Python API

a5py

Package for processing ASCOT5 data and generating inputs.

a5py.ascot5io

Interface for accessing data in ASCOT5 HDF5 files.

a5py.ascotpy

Python interface to libascot.so.

a5py.templates

Package for generating inputs from templates and imported data.

a5py.physlib

Tools for assessing physical quantities and units.

a5py.testascot

Unit tests, physics tests and regressions runs.

a5py.routines

Advanced pre- and postprocessing routines.

a5py.gui

Graphical user interface for ASCOT5.

Important classes

Ascot

Primary tool for processing ASCOT5 data.

Ascot5IO

Entry node for accessing data in the HDF5 file.

RunGroup

Node containing results and methods to process them.

AfsiGroup

Node containing AFSI results and methods to process them.

BBNBIGroup

Node containing BBNBI results and methods to process them.

Afsi

ASCOT Fusion Source Integrator AFSI.

BioSaw

Biot-Savart solver for computing magnetic field from a coil geometry.

MarkerGenerator

Tool for processing distributions into markers.

Ascot

class a5py.Ascot(inputfile=None, create=False, mute='err')

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 data attribute and its attributes that form a treeview. For example, active magnetic field is accessed as Ascot.data.bfield.active. See Ascot5IO for 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:
dataAscot5IO

Container for the HDF5 data.

biosawBioSaw

Tool for calculating magnetic field from coils.

afsiAfsi

Tool for calculating fusion source from thermal plasma and fast ion distributions.

markergenMarkerGenerator

Tool for generating markers from distributions.

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:
pathstr

Absolute path or None if this is an empty instance.

file_load(inputfile)

Open a new HDF5 file and free resources from the old one.

Parameters:
inputfilestr

Path to the new HDF5 file.

get_plasmaquantities()

Return species present in plasma input.

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.

*qntstr

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``.

Raises:
AssertionError

If required data has not been initialized.

RuntimeError

If evaluation in libascot.so failed.

input_eval_atomicsigma(ma, anum, znum, r, phi, z, t, va, ion, reaction)

Evaluate atomic reaction rate cross-sections for a given test particle.

Parameters:
mafloat

Test particle mass

anumint

Test particle atomic mass number.

znumint

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.

ionint

Index number of the background ion species in plasma input.

reaction{“ionization”, “recombination”, “charge-exchange”,
“beamstopping”}

Reaction whose cross-section is computed.

Returns:
sigmavarray_like, (n,nv)

Reaction cross-section.

Raises:
AssertionError

If required data has not been initialized.

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:
mafloat

Test particle mass.

qafloat

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.

*coefsstr

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).

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.

Raises:
AssertionError

If required data has not been initialized.

RuntimeError

If 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).

Parameters:
showbool, optional

Print output.

Returns:
quantitiesdict [str, str]

Name of each quantity and a short description.

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 deltacrit and ripplewell are approximate quantities based on analytical estimates.

Parameters:
rgridarray_like

R points where delta, deltacrit, and ripplewell are evaluated.

zgridarray_like

z points where delta, deltacrit, and ripplewell are evaluated.

rlarmorfloat or tuple (str, int, float)

Larmor radius of the test particle or tuple with (species, charge, energy).

rho0float, optional

Radial coordinate where amplitude is evaluated.

theta0float, optional

Poloidal coordinate where amplitude is evaluated.

tfloat, optional

Time instance for time-dependent fields.

nphiint, 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 amplitude if plot = True.

axes2, optional

Axes for plotting delta, deltacrit, and ripplewell if plot = True.

Returns:
phigridarray_like, (nphi)

Toroidal grid [-2, 362] where amplitude is 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.

input_eval_safetyfactor(rho, nth=10000)

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.

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.

input_findpsi0(psi1)

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. 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:
psi1float

Poloidal flux at the separatrix.

This value is used to deduce whether the algorithm searches minimum or maximum value when finding psi0.

Returns:
rfloat

Axis R-coordinate.

zfloat

Axis z-coordinate.

psi0float

Poloidal flux on axis.

Raises:
AssertionError

If required data has not been initialized.

RuntimeError

If evaluation in libascot.so failed.

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.

input_getplasmaspecies()

Get species present in plasma input (electrons first).

Returns:
nspeciesint

Number of species.

massarray_like

Species’ mass.

chargearray_like

Species’ charge state.

anumarray_like

Species’ atomic mass number.

znumarray_like

Species’ charge number.

Raises:
AssertionError

If required data has not been initialized.

RuntimeError

If 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:

  1. Without any arguments in which case all inputs are initialized. The initialized input is the one that has been set as active.

  2. 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.

  3. 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:
runstr or bool, optional

The run group or True to use the one that is active.

bfieldstr or bool or dict, optional

The magnetic field input or True to use the one that is active.

efieldstr or bool or dict, optional

The electric field input or True to use the one that is active.

plasmastr or bool or dict, optional

The plasma input or True to use the one that is active.

wallstr or bool or dict, optional

The wall input or True to use the one that is active.

neutralstr or bool or dict, optional

The neutral input or True to use the one that is active.

boozerstr or bool or dict, optional

The boozer input or True to use the one that is active.

mhdstr or bool or dict, optional

The MHD input or True to use the one that is active.

asigmastr or bool or dict, 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.

input_initialized()

Get inputs that are currently initialized.

Returns:
outputsdict [str, str]

Name of the input parent groups that are initialized and QID of the initialized input.

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:
rhofloat or array_like, optional

Value(s) of square root of normalized poloidal flux whose contours are plotted.

phifloat, optional

Toroidal angle at which contours are plotted.

tfloat, optional

Time instance when the contours are plotted.

nthetaint, optional

Poloidal resolution of the contours.

axesAxes, optional

The axes where figure is plotted or otherwise new figure is created.

**kwargs

Arguments passed to matplotlib.pyplot.plot().

input_plotrz(r, z, qnt, phi=unyt_quantity(0, 'degree'), t=unyt_quantity(0, 's'), clim=None, cmap=None, axes=None, cax=None)

Plot input quantity on a (R, z) plane at given coordinates.

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.

qntstr

Name of the plotted quantity.

phifloat, optional

Toroidal coordinate where data is evaluated.

tfloat, optional

Time coordinate where data is evaluated.

climlist [float, float], optional

Minimum and maximum color range values.

axesAxes, optional

The axes where figure is plotted or otherwise new figure is created.

caxAxes, optional

The color bar axes or otherwise taken from the main axes.

input_plotwallcontour(phi=unyt_quantity(0, 'degree'), axes=None)

Plot intersection of the wall and the poloidal plane at the given toroidal angle.

Parameters:
phifloat

Toroidal angle of the poloidal plane.

axesAxes, optional

The axes where figure is plotted or otherwise new figure is created.

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.

timefloat

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:
AssertionError

If required data has not been initialized.

RuntimeError

If 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)

Evaluate flux surface volumes.

Parameters:
nrhoint, optional

Number of radial grid edges between [0, 1].

nthetaint, optional

Number of poloidal grid edges.

nphiint, 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.

tolfloat, optional

If relative difference in volume in subsequent iterations falls below this value, the algorithm finishes.

tfloat, 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.

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.

Parameters:
hidewallbool, optional

Don’t show wall even if present.

hidemarkersbool, optional

Don’t show markers even if present.

axesAxes, optional

The axes where figure is plotted or otherwise new figure is created.

simulation_free(inputs=False, markers=False, diagnostics=False)

Free resources used by the interactive simulation.

If called without arguments, everything will be freed.

Parameters:
inputsbool, optional

If True, inputs are unpacked and freed.

markersbool, optional

If True, markers are freed.

diagnosticsbool, optional

If True, diagnostics data is freed.

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 with simulation_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(), or FL.write_hdf5().

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().

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.

Returns:
runVirtualRun

An object that acts almost exactly as RunGroup except 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.

Raises:
AscotInitException

If 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:
_ascotAscot

Ascot object used to run AFSI.

beambeam(reaction, beam1, beam2=None, nmc=1000, minppara=-1.3e-19, maxppara=1.3e-19, nppara=80, minpperp=0, maxpperp=1.3e-19, npperp=40)

Calculate beam-beam fusion.

Parameters:
reactionint

Fusion reaction index.

beam1dict

Beam1 distribution.

beam2dict, optional

Beam2 distribution or None to calculate fusion generation with beam1 itself.

nmcint, optional

Number of MC samples used in each (R, phi, z) bin.

minpparafloat, optional

Minimum ppara value in the output distribution.

maxpparafloat, optional

Maximum ppara value in the output distribution.

npparaint, optional

Number of pperp bins in the output distribution.

minpperpfloat, optional

Minimum pperp value in the output distribution.

maxpperpfloat, optional

Maximum ppara value in the output distribution.

npperpint, optional

Number of pperp bins in the output distribution.

Returns:
prod1array_like

Fusion product 1 distribution.

prod2array_like

Fusion product 2 distribution.

beamthermal(reaction, beam, swap=False, nmc=1000, minppara=-1.3e-19, maxppara=1.3e-19, nppara=80, minpperp=0, maxpperp=1.3e-19, npperp=40)

Calculate beam-thermal fusion.

Parameters:
reactionint

Fusion reaction index

beamdict

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.

nmcint, optional

Number of MC samples used in each (R, phi, z) bin.

minpparafloat, optional

Minimum ppara value in the output distribution.

maxpparafloat, optional

Maximum ppara value in the output distribution.

npparaint, optional

Number of pperp bins in the output distribution.

minpperpfloat, optional

Minimum pperp value in the output distribution.

maxpperpfloat, optional

Maximum ppara value in the output distribution.

npperpint, optional

Number of pperp bins in the output distribution.

Returns:
prod1array_like

Fusion product 1 distribution.

prod2array_like

Fusion product 2 distribution.

reactions(reaction=None)

Return reaction data for a given reaction.

Parameters:
reactionstr, optional

Reaction or None to return list of available reactions.

Returns:
reactions[str]

List of available reactions if reaction=None.

m1float

Mass of reactant 1.

q1float

Charge of reactant 1.

m2float

Mass of reactant 2.

q2float

Charge of reactant 2.

mprod1float

Mass of product 1.

qprod1float

Charge of product 1.

mprod2float

Mass of product 2.

qprod2float

Charge of product 2.

qfloat

Energy released in the reaction.

thermal(reaction, minr, maxr, nr, minz, maxz, nz, minphi=0, maxphi=6.283185307179586, nphi=1, nmc=1000, minppara=-1.3e-19, maxppara=1.3e-19, nppara=80, minpperp=0, maxpperp=1.3e-19, npperp=40)

Calculate thermonuclear fusion between two thermal (Maxwellian) species.

Parameters:
reactionint or str

Fusion reaction index or name.

minrfloat

Minimum R value in the output distribution.

maxrfloat

Maximum R value in the output distribution.

nrint

Number of R bins in the output distribution.

minzfloat

Minimum z value in the output distribution.

maxzfloat

Maximum z value in the output distribution.

nzint

Number of z bins in the output distribution.

minphifloat, optional

Minimum phi value in the output distribution.

maxphifloat, optional

Maximum phi value in the output distribution.

nphiint, optional

Number of phi bins in the output distribution.

nmcint, optional

Number of MC samples used in each (R, phi, z) bin.

minpparafloat, optional

Minimum ppara value in the output distribution.

maxpparafloat, optional

Maximum ppara value in the output distribution.

npparaint, optional

Number of pperp bins in the output distribution.

minpperpfloat, optional

Minimum pperp value in the output distribution.

maxpperpfloat, optional

Maximum ppara value in the output distribution.

npperpint, optional

Number of pperp bins in the output distribution.

Returns:
prod1array_like

Fusion product 1 distribution.

prod2array_like

Fusion product 2 distribution.

class a5py.routines.biosaw5.BioSaw(ascot)

Biot-Savart solver for computing magnetic field from a coil geometry.

Attributes:
_ascotAscot

Ascot object used to run BioSaw.

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_2DS input to generate the B_3DS input.

Note that the toroidal grid is defined as:

np.linspace(phimin, phimax, nphi+1)[:-1]

Parameters:
coilxyzlist or 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.

phiminfloat

Toroidal grid minimum value.

phimaxfloat

Toroidal grid maximum value.

nphiint

Number of toroidal grid points.

rminfloat, optional

Minimum value in R grid unless taken from the input data.

rmaxfloat, optional

Maximum value in R grid unless taken from the input data.

nrint, optional

Number of R grid points unless taken from the input data.

zminfloat, optional

Minimum value in z grid unless taken from the input data.

zmaxfloat, optional

Maximum value in z grid unless taken from the input data.

nzint, optional

Number of z grid points unless taken from the input data.

currentlist or float

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 b0 is provided.

revolveint, 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 phi array.

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”} or float, 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.

b2ddict, optional

B_2DS data to be used or otherwise read from active bfield input.

Returns:
outdict

B_3DS data.

addto3d(coilxyz, current=1.0, revolve=None, b3d=None)
Generate new 3D tokamak magnetic field input where contribution from

given coils are included.

Parameters:
coilxyzlist or 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.

currentlist or float

Coil current in Ampères either same for all coils or given individually.

revolveint, 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 phi array.

b3ddict, optional

B_3DS data to be used or otherwise read from active bfield input.

Returns:
outdict

B_3DS data.

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.

coilxyzlist or 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.

currentlist or float

Coil current in Ampères either same for all coils or given individually.

revolveint, 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 phi array.

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:
_ascotAscot

The Ascot object this instance belongs to.

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:
nmrkint

Number of markers to be generated.

massfloat

Particle species` mass.

chargefloat

Particle species` charge.

anumint

Particle species` atomic mass number.

znumint

Particle species` charge number.

particledistDist

5D distribution from which marker weight is sampled.

markerdistDist, 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.

minweightfloat, optional

Minimum weight a marker must have or else it is rejected.

return_distsbool, optional

Return markerdist and particledist calculated from the generated markers.

These distributions hould converge to inputs as the marker number is increased. The exception is particledist which will always get zero values in regions where markerdist is zero.

Returns:
mrkdict

Marker input data on a dictionary.

mrkdistDist

Distribution of initialized markers which ideally should be close to the one given as an input unless markers were rejected by minweight.

prtdistDist

Weighted distribution of initialized markers which ideally should be close to the one given as an input, unless some phase-space regions in markerdist were left intentionally empty even though they would contain physical particles.

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:
markerdistDist

Normalized 5D distribution from which markers can be sampled.

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.