bfield

Magnetic field input.

Magnetic field is present in every simulation and plays for an important role in accuracy of the results. High-fidelity magnetic field is essential. Most of the time magnetic field is the most memory-intensive input (and also most of the CPU time is spent interpolating it).

The poloidal flux psi is also used as a radial coordinate in simulations e.g. to interpolate 1D data. More precisely, rho = sqrt( (psi - psi0) / ( psi1 - psi0 ) ), where psi0 is psi on axis and psi1 is psi on separatrix, is used. Note that psi0 and psi1 must be provided explicitly in all magnetic field inputs. It is a good idea to include some padding to psi0, as otherwise rho might have a complex value near separatrix (which leads to markers being aborted there) if the given psi0 differs from the actual extrema of the interpolation scheme.

a5py.ascot5io.bfield.B_2DS

Axisymmetric tokamak field interpolated with splines.

a5py.ascot5io.bfield.B_3DS

Non-axisymmetric tokamak field.

a5py.ascot5io.bfield.B_STS

Stellarator magnetic field.

a5py.ascot5io.bfield.B_GS

Analytical tokamak field.

a5py.ascot5io.bfield.B_TC

Magnetic field in Cartesian basis for testing purposes.

class a5py.ascot5io.bfield.B_2DS(root, path, **kwargs)

Bases: DataGroup

Axisymmetric tokamak field interpolated with splines.

This field is suitable for simulations where toroidal ripple or other 3D perturbations are not relevant (for MHD one can use the dedicated input). Using axisymmetric field is much faster than 3D fields, so using this input when applicable is strongly recommended.

The input consists of BR, Bphi, Bz, and psi tabulated on an uniform (R,z) grid which are interpolated with cubic splines in simulation. During the simulation, the total BR (and Bz) are computed as a sum of input BR (Bz) and the component calculated from the gradient of psi. In other words, the equilibrium BR and Bz are already contained in psi so in most cases the input BR and Bz should be set to zero. If for some reason psi has a poor quality, it can be scaled to an insignificant value and the field can be provided explicitly in BR and Bz. However, in such case the field is no longer divergence-free.

static convert_B_3DS(**kwargs)

Convert B_3DS input to B_2DS input.

This function takes a toroidal average of Bphi and sets BR and Bz to zero.

Parameters:
**kwargs

Arguments passed to B_3DS.write_hdf5() excluding fn and desc.

Returns:
outdict

B_3DS converted as an input for write_hdf5().

static convert_B_GS(rmin, rmax, nr, zmin, zmax, nz, **kwargs)

Convert B_GS input to B_2DS input.

Parameters:
rminfloat

R grid min edge [m].

rmaxfloat

R grid max edge [m].

nrint

Number of R grid points.

zminfloat

z grid min edge [m].

zmaxfloat

z grid max edge [m].

nzint

Number of z grid points.

**kwargs

Arguments passed to B_GS.write_hdf5() excluding fn and desc.

Returns:
outdict

B_GS converted as an input for write_hdf5().

static create_dummy()

Create dummy data that has correct format and is valid, but can be non-sensical.

This method is intended for testing purposes or to provide data whose presence is needed but which is not actually used in simulation.

The dummy output is an ITER-like but circular equilibrium.

Returns:
datadict

Input data that can be passed to write_hdf5 method of a corresponding type.

read()

Read data from HDF5 file.

Returns:
datadict

Data read from HDF5 stored in the same format as is passed to write_hdf5().

static write_hdf5(fn, rmin, rmax, nr, zmin, zmax, nz, axisr, axisz, psi, psi0, psi1, br, bphi, bz, desc=None)

Write input data to the HDF5 file.

Note that br and bz should not include the equilibrium component of the magnetic field as that is calculated from psi by ASCOT5 during the simulation.

Parameters:
fnstr

Full path to the HDF5 file.

rminfloat

R grid min edge [m].

rmaxfloat

R grid max edge [m].

nrint

Number of R grid points.

zminfloat

z grid min edge [m].

zmaxfloat

z grid max edge [m].

nzint

Number of z grid points.

axisrfloat

Magnetic axis R coordinate [m].

axiszfloat

Magnetic axis z coordinate [m].

psi0float

On-axis poloidal flux value [Vs/m].

psi1float

Separatrix poloidal flux value [Vs/m].

psiarray_like (nr, nz)

Poloidal flux values on the Rz grid [Vs/m].

brarray_like (nr,nz)

Magnetic field R component (excl. equilibrium comp.) on Rz grid [T].

bphiarray_like (nr,nz)

Magnetic field phi component on Rz grid [T].

bzarray_like (nr,nz)

Magnetic field z component (excl. equilibrium comp.) onRz grid [T].

descstr, optional

Input description.

Returns:
namestr

Name, i.e. “<type>_<qid>”, of the new input that was written.

Raises:
ValueError

If inputs were not consistent.

class a5py.ascot5io.bfield.B_3DS(root, path, **kwargs)

Bases: DataGroup

Non-axisymmetric tokamak field.

The input consists of BR, Bphi, and Bz tabulated on an uniform cylindrical grid and psi tabulated on (R,z) grid. Data is interpolated with cubic splines during the simulation, and the total BR (and Bz) are computed as a sum of input BR (Bz) and the component calculated from the gradient of psi. In other words, the equilibrium BR and Bz are already contained in psi so BR and Bz should only contain the components from 3D perturbation.

This input is suitable for artificial 3D perturbations. For MHD, consider the dedicated input.

Since this input interpolates BR and Bz with splines, the resulting B is not divergence free. This can be amended by using a dense grid for magnetic field evaluation. Note however that it can become memory-intensive.

static convert_B_GS(rmin, rmax, nr, zmin, zmax, nz, phimin, phimax, nphi, **kwargs)

Convert B_GS input to B_3DS input.

The resulting field is far from being divergence-free as the perturbation BR and Bz components are omitted.

Parameters:
rminfloat

Magnetic field data R grid min edge [m].

rmaxfloat

Magnetic field data R grid max edge [m].

nrint

Number of R grid points in magnetic field data.

zminfloat

Magnetic field data z grid min edge [m].

zmaxfloat

Magnetic field data z grid max edge [m].

nzint

Number of z grid points in magnetic field data.

phiminfloat

Beginning of the toroidal period [deg].

phimaxfloat

End of the toroidal period [deg].

nphiint

Number of phi grid points in magnetic field data.

**kwargs

Arguments passed to B_GS.write_hdf5() excluding fn and desc.

Returns:
outdict

B_GS converted as an input for write_hdf5().

static create_dummy()

Create dummy data that has correct format and is valid, but can be non-sensical.

This method is intended for testing purposes or to provide data whose presence is needed but which is not actually used in simulation.

Returns:
datadict

Input data that can be passed to write_hdf5 method of a corresponding type.

read()

Read data from HDF5 file.

Returns:
datadict

Data read from HDF5 stored in the same format as is passed to write_hdf5().

static write_hdf5(fn, b_rmin, b_rmax, b_nr, b_zmin, b_zmax, b_nz, b_phimin, b_phimax, b_nphi, axisr, axisz, psi, psi0, psi1, br, bphi, bz, psi_rmin=None, psi_rmax=None, psi_nr=None, psi_zmin=None, psi_zmax=None, psi_nz=None, desc=None)

Write input data to the HDF5 file.

It is possible to use different (R,z) grids for psi and magnetic field components by giving the (R,z)-grid for psi separately. 3D data can be memory intensive which necessitates sparser grid for B components, but psi can still be evaluated on a dense grid.

The toroidal angle phi is treated as a periodic coordinate, meaning A(phi=phimin) == A(phi=phimax). However, the phi grid, where input arrays are tabulated, is linspace(phimin, phimax, nphi+1)[:-1] to avoid storing duplicate data.

Parameters:
fnstr

Full path to the HDF5 file.

b_rminfloat

Magnetic field data R grid min edge [m].

b_rmaxfloat

Magnetic field data R grid max edge [m].

b_nrint

Number of R grid points in magnetic field data.

b_zminfloat

Magnetic field data z grid min edge [m].

b_zmaxfloat

Magnetic field data z grid max edge [m].

b_nzint

Number of z grid points in magnetic field data.

b_phiminfloat

Beginning of the toroidal period [deg].

b_phimaxfloat

End of the toroidal period [deg].

b_nphiint

Number of phi grid points in magnetic field data.

axisrfloat

Magnetic axis R coordinate [m].

axiszfloat

Magnetic axis z coordinate [m].

psi0float

On-axis poloidal flux value [Vs/m].

psi1float

Separatrix poloidal flux value [Vs/m].

psiarray_like (nr, nz)

Poloidal flux values on the Rz grid [Vs/m].

brarray_like (nr,nphi,nz)

Magnetic field R component (excl. equilibrium comp.) on Rz grid [T].

bphiarray_like (nr,nphi,nz)

Magnetic field phi component on Rz grid [T].

bzarray_like (nr,nphi,nz)

Magnetic field z component (excl. equilibrium comp.) onRz grid [T].

psi_rminfloat, optional

Psi data R grid min edge [m].

psi_rmaxfloat, optional

Psi data R grid max edge [m].

psi_nrint, optional

Number of R grid points in psi data.

psi_zminfloat, optional

Psi data z grid min edge [m].

psi_zmaxfloat, optional

Psi data z grid max edge [m].

psi_nzint, optional

Number of z grid points in psi data.

descstr, optional

Input description.

Returns:
namestr

Name, i.e. “<type>_<qid>”, of the new input that was written.

Raises:
ValueError

If inputs were not consistent.

class a5py.ascot5io.bfield.B_STS(root, path, **kwargs)

Bases: DataGroup

Stellarator magnetic field.

static convert_B_GS(R0, z0, B_phi0, psi_mult, psi_coeff, Rmin, Rmax, nR, zmin, zmax, nz, phimin, phimax, nphi, psi0=None, raxis=None, zaxis=None, desc=None, **kwargs)

Convert B_GS input to B_STS input.

The resulting field is far from being divergence-free as the perturbation BR and Bz components are omitted.

Parameters:
R0, z0real

Major axis Rz-coordinates.

B_phi0real

Toroidal field at axis.

psi_multreal

Scaling factor for psi.

psi_coeffreal 13 x 1 numpy array

Coefficients defining psi.

Rmin, Rmax, zmin, zmax, phimin, phimaxreal

Edges of the uniform Rphiz-grid.

nR, nz, nphiint

Number of Rphiz-grid points.

psi0float, optional <br>

Poloidal flux at magnetic axis.

raxisfloat, optional <br>

Magnetic axis R coordinate [m].

zaxisfloat, optional <br>

Magnetic axis z coordinate [m].

**kwargs

Arguments passed to B_GS.write_hdf5() excluding fn and desc.

Returns:
outdict

B_GS converted as an input for write_hdf5().

static create_dummy()

Create dummy data that has correct format and is valid, but can be non-sensical.

This method is intended for testing purposes or to provide data whose presence is needed but which is not actually used in simulation.

Returns:
datadict

Input data that can be passed to write_hdf5 method of a corresponding type.

read()

Read data from HDF5 file.

Returns:
datadict

Data read from HDF5 stored in the same format as is passed to write_hdf5().

static write_hdf5(fn, b_rmin, b_rmax, b_nr, b_zmin, b_zmax, b_nz, b_phimin, b_phimax, b_nphi, psi0, psi1, br, bphi, bz, psi, axis_phimin, axis_phimax, axis_nphi, axisr, axisz, psi_rmin=None, psi_rmax=None, psi_nr=None, psi_zmin=None, psi_zmax=None, psi_nz=None, psi_phimin=None, psi_phimax=None, psi_nphi=None, desc=None)

Write input data to the HDF5 file.

It is possible to use different (R,z) grids for psi and magnetic field components by giving the (R,z)-grid for psi separately. 3D data can be memory intensive which necessitates sparser grid for B components, but psi can still be evaluated on a dense grid.

The toroidal angle phi is treated as a periodic coordinate, meaning A(phi=phimin) == A(phi=phimax). However, the phi grid, where input arrays are tabulated, is linspace(phimin, phimax, nphi+1)[:-1] to avoid storing duplicate data.

Parameters:
fnstr

Full path to the HDF5 file.

b_rminfloat

Magnetic field data R grid min edge [m].

b_rmaxfloat

Magnetic field data R grid max edge [m].

b_nrint

Number of R grid points in magnetic field data.

b_zminfloat

Magnetic field data z grid min edge [m].

b_zmaxfloat

Magnetic field data z grid max edge [m].

b_nzint

Number of z grid points in magnetic field data.

b_phiminfloat

Beginning of the toroidal period [deg].

b_phimaxfloat

End of the toroidal period [deg].

b_nphiint

Number of phi grid points in magnetic field data.

axis_phiminfloat

Magnetic axis phi grid min value [deg].

axis_phimaxfloat

Magnetic axis phi grid max value [deg].

axis_nphifloat

Number of points in magnetic axis phi grid.

axisrfloat

Magnetic axis R coordinates on axis phi grid [m].

axiszfloat

Magnetic axis z coordinates on axis phi grid [m].

psi0float

On-axis poloidal flux value [Vs/m].

psi1float

Separatrix poloidal flux value [Vs/m].

psiarray_like (psi_nr,psi_nphi,psi_nz)

Poloidal flux values on the Rz grid [Vs/m].

brarray_like (b_nr,b_nphi,b_nz)

Magnetic field R component (excl. equilibrium comp.) on Rz grid [T].

bphiarray_like (b_nr,b_nphi,b_nz)

Magnetic field phi component on Rz grid [T].

bzarray_like (b_nr,b_nphi,b_nz)

Magnetic field z component (excl. equilibrium comp.) onRz grid [T].

psi_rminfloat, optional

Psi data R grid min edge [m].

psi_rmaxfloat, optional

Psi data R grid max edge [m].

psi_nrint, optional

Number of R grid points in psi data.

psi_zminfloat, optional

Psi data z grid min edge [m].

psi_zmaxfloat, optional

Psi data z grid max edge [m].

psi_nzint, optional

Number of z grid points in psi data.

psi_phiminfloat, optional

Psi data beginning of the toroidal period.

psi_phimaxfloat, optional

Psi data end of the toroidal period.

psi_nphiint, optional

Number of phi grid points in psi data.

descstr, optional

Input description.

Returns:
namestr

Name, i.e. “<type>_<qid>”, of the new input that was written.

Raises:
ValueError

If inputs were not consistent.

class a5py.ascot5io.bfield.B_GS(root, path, **kwargs)

Bases: DataGroup

Analytical tokamak field.

This field can either be axisymmetric or it can include a ripple-like perturbation. The field is very fast to interpolate and it is not memory intensive. However, it is still intended mostly for testing purposes.

The axisymmetric field is divergence-free whereas the 3D field is not as it included only the toroidal component of the perturbation.

See: A.J. Cerfon, J.P. Freidberg. “One size fits all” analytic solutions to the Grad-Shafranov equation. Physics of Plasmas 17 (3) (2010) 032502. http://scitation.aip.org/content/aip/journal/pop/17/3/10.1063/1.3328818

static create_dummy()

Create dummy data that has correct format and is valid, but can be non-sensical.

This method is intended for testing purposes or to provide data whose presence is needed but which is not actually used in simulation.

Returns:
datadict

Input data that can be passed to write_hdf5 method of a corresponding type.

read()

Read data from HDF5 file.

Returns:
datadict

Data read from HDF5 stored in the same format as is passed to write_hdf5().

static write_hdf5(fn, r0, z0, bphi0, psimult, coefficients, psi0=None, psi1=0, raxis=None, zaxis=None, nripple=0, a0=2, alpha0=2, delta0=0.05, desc=None)

Write input data to the HDF5 file.

If psi0 is None, magnetic axis location is found numerically and psi0 value is evaluated at that point. The analytical field is defined by assuming midplane is at z=0, but it is moved to z=z0 here.

Setting number of toroidal field coils to non-zero makes the field 3D.

Parameters:
fnstr

Full path to the HDF5 file.

r0float

Major radius R coordinate [m].

z0float

Distance by which midplane is moved from z=0 plane [m].

bphi0float

Toroidal field at axis [T].

psimultfloat

Scaling factor for psi.

coefficientsarray_like (13,1)

Coefficients defining psi [c0, c1, …, c11, A].

psi0float, optional

Poloidal flux at magnetic axis.

psi1float, optional

Poloidal flux at the separatrix.

raxisfloat, optional

Magnetic axis R coordinate [m].

zaxisfloat, optional

Magnetic axis z coordinate [m].

nripplefloat, optional

Number of TF coils.

a0float, optional

Minor radius. [m]

alpha0float, optional

Ripple penetration.

delta0float, optional

Ripple (maximum) strength.

descstr, optional

Input description.

Returns:
namestr

Name, i.e. “<type>_<qid>”, of the new input that was written.

Raises:
ValueError

If inputs were not consistent.

class a5py.ascot5io.bfield.B_TC(root, path, **kwargs)

Bases: DataGroup

Magnetic field in Cartesian basis for testing purposes.

This input defines the magnetic field vector on the Cartesian basis, and uses constant Jacobian to make the field non-uniform if needed. This field is only used to test the validity of the orbit-integrators.

static create_dummy()

Create dummy data that has correct format and is valid, but can be non-sensical.

This method is intended for testing purposes or to provide data whose presence is needed but which is not actually used in simulation.

Returns:
datadict

Input data that can be passed to write_hdf5 method of a corresponding type.

read()

Read data from HDF5 file.

Returns:
datadict

Data read from HDF5 stored in the same format as is passed to write_hdf5().

static write_hdf5(fn, bxyz, jacobian, rhoval, psival=None, axisr=1, axisz=0, desc=None)

Write input data to the HDF5 file.

Parameters:
fnstr

Full path to the HDF5 file.

bxyzarray_like (3,1)

Magnetic field in cartesian coordinates at origo.

jacobianarray_like (3,3)

Magnetic field Jacobian, jacobian[i,j] = dB_i/dx_j.

rhoval: float

Constant rho value.

psival: float, optional

Constant psi value. If None, same as rhoval.

axisr: float, optional

Magnetic axis R coordinate.

axisz: real, optional

Magnetic axis z coordinate.

descstr, optional

Input description.

Returns:
namestr

Name, i.e. “<type>_<qid>”, of the new input that was written.

Raises:
ValueError

If inputs were not consistent.