API Reference
Benchmark
IMAS.Benchmark — Type
Benchmark{T}Struct for benchmarking time-dependent dd data structures
IMAS.benchmark — Function
benchmark(dd_ben::IMAS.dd{T}, times::AbstractVector{Float64}) where {T<:Real}Evaluate self-consistency of a given IDS (flux-matching errors, etc.)
benchmark(dd_ben::IMAS.dd{T}, dd_ref::IMAS.dd{T}, times::AbstractVector{Float64}; self_benchmark::Bool=true) where {T<:Real}Benchmark two dd structures at specified times
benchmark(x::AbstractVector{T}, y::AbstractVector{T}, x_ref::AbstractVector{T}, y_ref::AbstractVector{T}) where {T<:Real}evaluate error between two 1D arrays
Control
IMAS.controller — Function
controller(ct::IMAS.controllers, name::String)Pick a controller based on its name, returns nothing if not found
NOTE: for now only looks under dd.controllers.linear_controller
IMASdd.controllers__linear_controller — Type
(controller::controllers__linear_controller{T})(control_name::Any, control_algorithm::Val{:PID}, setpoint::T, value::T, time0::Float64) where {T<:Real}Operates PID controllers
Extract
IMAS.extract — Function
extract(dd::IMAS.dd, library::Symbol=:extract)library can be one of:
:extract=>ExtractFunctionsLibrary:moopt=>ConstraintFunctionsLibrary+ObjectiveFunctionsLibrary:all=>ExtractFunctionsLibrary+ConstraintFunctionsLibrary+ObjectiveFunctionsLibrary
extract(dd::IMAS.dd, xtract::AbstractDict{Symbol,<:ExtractFunction})Extract data from dd. Each of the ExtractFunction should accept dd as input, like this:
xtract = IMAS.ExtractFunction[
:κ => ExtractFunction(:equilibrium, :κ, "-", dd -> dd.equilibrium.time_slice[].boundary.elongation)
:Te0 => ExtractFunction(:profiles, :Te0, "keV", dd -> dd.core_profiles.profiles_1d[].electrons.temperature[1] / 1E3)
]extract(dd::IMAS.dd, filename::AbstractString, args...; kw...)Like extract(dd) but saves the data to HDF5 file
Functions library
IMAS.ConstraintFunction — Type
name::Symbol
units::String
func::Function
operation::Function
limit::Float64
tolerance::Float64IMAS.ConstraintFunctionsLibrary — Constant
ConstraintFunctionsLibrary::OrderedCollections.OrderedDict{Symbol,ConstraintFunction}Collection of ConstraintFunction
power_electric_net → == 0.0 ± 0.2 [%]min_power_electric_net → > 0.0 [%]flattop_duration → == 0.0 ± 0.01 [%]min_flattop_duration → > 0.0 [%]min_required_B0 → > 0.0 [%]zero_ohmic → == 0.0 ± 0.1 [MA]max_ne_peaking → < 0.0 [%]min_lh_power_threshold_fraction → > 1.0 [%]max_ωpe_ωce → < 1.0 [%]max_qpol_omp → < 0.0 [%]max_tf_coil_j → < 1.0 [%]max_oh_coil_j → < 1.0 [%]max_pl_stress → < 1.0 [%]max_tf_coil_stress → < 1.0 [%]max_oh_coil_stress → < 1.0 [%]max_hds03 → < 0.0 [%]min_q95 → > 0.0 [%]min_qmin → > 0.0 [%]max_beta_normal → < 0.0 []max_Psol_R → < 0.0 [%]max_transport_error → < 0.1 []
IMAS.ObjectiveFunction — Type
name::Symbol
units::String
func::Function
target::Float64IMAS.ObjectiveFunctionsLibrary — Constant
ObjectiveFunctionsLibrary::OrderedCollections.OrderedDict{Symbol,ObjectiveFunction}Collection of ObjectiveFunction
min_levelized_CoE → -Inf [$/kWh]min_log10_levelized_CoE → -Inf [log₁₀($/kW)]min_capital_cost → -Inf [$B]max_fusion → Inf [MW]max_power_electric_net → Inf [MW]req_power_electric_net → 0.0 [ΔMW]max_flattop → Inf [hours]max_q95 → Inf []req_flattop → 0.0 [Δhours]max_log10_flattop → Inf [log₁₀(hours)]min_βn → -Inf []min_R0 → -Inf [m]max_zeff → Inf []
IMAS.ExtractFunction — Type
group::Symbol
name::Symbol
units::String
func::Function
error::Union{Nothing,Exception}
value::AnyIMAS.ExtractFunctionsLibrary — Constant
ExtractFunctionsLibrary::OrderedCollections.OrderedDict{Symbol,ExtractFunction}Collection of ExtractFunction
geometry.R0 → [m]geometry.a → [m]geometry.1/ϵ → [-]geometry.κ → [-]geometry.δ → [-]geometry.ζ → [-]geometry.Volume → [m³]geometry.Surface → [m²]equilibrium.B0 → [T]equilibrium.ip → [MA]equilibrium.q95 → [-]equilibrium.<Bpol> → [T]equilibrium.βpol_MHD → [-]equilibrium.βtor_MHD → [-]equilibrium.βn_MHD → [-]temperatures.Te0 → [keV]temperatures.Ti0 → [keV]temperatures.<Te> → [keV]temperatures.<Ti> → [keV]temperatures.Te0/<Te> → [-]temperatures.Ti0/<Ti> → [-]densities.ne0 → [m⁻³]densities.ne_ped → [m⁻³]densities.ne_line → [m⁻³]densities.<ne> → [m⁻³]densities.ne0/<ne> → [-]densities.fGW → [-]densities.zeff_ped → [-]densities.<zeff> → [-]densities.impurities → [-]pressures.P0 → [MPa]pressures.<P> → [MPa]pressures.P0/<P> → [-]pressures.βn → [-]pressures.βn_th → [-]transport.τe → [s]transport.τe_exp → [s]transport.H98y2 → [-]transport.H98y2_exp → [-]transport.Hds03 → [-]transport.Hds03_exp → [-]transport.τα_thermalization → [s]transport.τα_slowing_down → [s]sources.Pec → [MW]sources.rho0_ec → [MW]sources.Pnbi → [MW]sources.Enbi1 → [MeV]sources.Pic → [MW]sources.Plh → [MW]sources.Paux_tot → [MW]sources.Pα → [MW]sources.Pohm → [MW]sources.Pheat → [MW]sources.Prad_tot → [MW]exhaust.Psol → [MW]exhaust.PLH → [MW]exhaust.Bpol_omp → [T]exhaust.λq → [mm]exhaust.qpol → [MW/m²]exhaust.qpar → [MW/m²]exhaust.P/R0 → [MW/m]exhaust.PB/R0 → [MW T/m]exhaust.PBp/R0 → [MW T/m]exhaust.PBϵ/R0q95 → [MW T/m]exhaust.neutrons_peak → [MW/m²]currents.ip_bs_aux_ohm → [MA]currents.ip_ni → [MA]currents.ip_bs → [MA]currents.ip_aux → [MA]currents.ip_ohm → [MA]currents.ejima → [-]currents.flattop → [Hours]bop.Pfusion → [MW]bop.Qfusion → [-]bop.thermal_cycle_type → [-]bop.thermal_efficiency_plant → [%]bop.thermal_efficiency_cycle → [%]bop.power_electric_generated → [MW]bop.Pelectric_net → [MW]bop.Qplant → [-]bop.TBR → [-]build.PF_material → [-]build.TF_material → [-]build.OH_material → [-]build.TF_max_b → [T]build.OH_max_b → [T]build.TF_j_margin → [-]build.OH_j_margin → [-]build.TF_stress_margin → [-]build.OH_stress_margin → [-]costing.capital_cost → [$B]costing.levelized_CoE → [$/kWh]costing.TF_of_total → [%]costing.BOP_of_total → [%]costing.blanket_of_total → [%]costing.cryostat_of_total → [%]
Geometry
IMAS.centroid — Function
centroid(x::AbstractVector{<:T}, y::AbstractVector{<:T}) where {T<:Real}Calculate centroid of polygon
IMAS.perimeter — Function
perimeter(r::AbstractVector{T}, z::AbstractVector{T})::T where {T<:Real}Calculate the perimeter of a polygon
IMAS.revolution_volume — Function
revolution_volume(x::AbstractVector{<:T}, y::AbstractVector{<:T}) where {T<:Real}Calculate volume of polygon revolved around x=0
IMAS.intersection_angles — Function
intersection_angles(
path1_r::AbstractVector{T},
path1_z::AbstractVector{T},
path2_r::AbstractVector{T},
path2_z::AbstractVector{T},
intersection_indexes::Vector{StaticArrays.SVector{2,Int}};
mod_pi::Bool=true
) where {T<:Real}returns angles of intersections between two paths and intersection_indexes given by intersection() function
IMAS.intersection — Function
intersection(
l1_x::AbstractVector{T},
l1_y::AbstractVector{T},
l2_x::AbstractVector{T},
l2_y::AbstractVector{T}
) where {T<:Real}Intersections between two 2D paths, returns list of (x,y) intersection indexes and crossing points
intersection(
l1_x::AbstractVector{T},
l1_y::AbstractVector{T},
l2_x::AbstractVector{T},
l2_y::AbstractVector{T},
tolerance::Float64) where {T<:Real}Intersections between two 2D paths, returns list of (x,y) intersection indexes and crossing points
Endpoints crossings are checked with some tolerance
IMAS.intersection_split — Function
intersection_split(
l1_x::AbstractVector{T},
l1_y::AbstractVector{T},
l2_x::AbstractVector{T},
l2_y::AbstractVector{T}) where {T<:Real}Returns vector of segments of l1x,l1y split at the intersections with l2x,l2y
IMAS.point_to_line_distance — Function
point_to_line_distance(x0::Real, y0::Real, x1::Real, y1::Real, x2::Real, y2::Real)Distance of point (x0,y0) from line defined by points (x1,y1) and (x2,y2)
IMAS.closest_point_to_segment — Function
closest_point_to_segment(x0::Real, y0::Real, x1::Real, y1::Real, x2::Real, y2::Real)Closest point on segment defined by points (x1,y1) and (x2,y2) to point (x0,y0)
IMAS.point_to_segment_distance — Function
point_to_segment_distance(x0::Real, y0::Real, x1::Real, y1::Real, x2::Real, y2::Real)Distance of point (x0,y0) from segment defined by points (x1,y1) and (x2,y2)
IMAS.point_to_path_distance — Function
point_to_path_distance(x0::Real, y0::Real, x::AbstractVector{<:Real}, y::AbstractVector{<:Real})Distance of point (x0,y0) from path defined by vectors x and y
IMAS.rdp_simplify_2d_path — Function
rdp_simplify_2d_path(x::AbstractArray{T}, y::AbstractArray{T}, epsilon::T) where {T<:Real}Simplifies a 2D line represented by arrays of x and y coordinates using the Ramer-Douglas-Peucker algorithm. The epsilon parameter controls the maximum distance allowed between a point on the original line and its simplified representation.
IMAS.rwa_simplify_2d_path — Function
rwa_simplify_2d_path(x::AbstractArray{T}, y::AbstractArray{T}, epsilon::T) where {T<:Real}Simplifies a 2D line represented by arrays of x and y coordinates using the Reumann-Witkam Algorithm algorithm. This algorithm uses a threshold value to determine which points to keep in the path. Points are kept if the angle between the previous and next line segments is greater than the threshold, and removed if it is less than or equal to the threshold.
IMAS.calculate_angle — Function
calculate_angle(p1::T, p2::T, p3::T) where {T}Calculate the angle between three points
IMAS.simplify_2d_path — Function
simplify_2d_path(x::AbstractArray{T}, y::AbstractArray{T}, simplification_factor::T; model::Symbol=:distance)Simplify 2D path by :curvature (Reumann-Witkam Algorithm) or :distance (Ramer-Douglas-Peucker) algorithms
IMAS.resample_2d_path — Function
resample_2d_path(
x::AbstractVector{T},
y::AbstractVector{T};
step::Float64=0.0,
n_points::Integer=0,
curvature_weight::Float64=0.0,
retain_extrema::Bool=false,
retain_original_xy::Bool=false,
method::Symbol=:cubic) where {T<:Real}Resample 2D line with uniform stepping (or number of points) with option to add more points where curvature is highest and option to retain extrema in x and y (in these cases stepping is not constant anymore!)
IMAS.resample_plasma_boundary — Function
resample_plasma_boundary(
x::AbstractVector{T},
y::AbstractVector{T};
step::Float64=0.0,
n_points::Integer=0,
curvature_weight::Float64=0.0,
retain_extrema::Bool=true,
retain_original_xy::Bool=false,
method::Symbol=:linear) where {T<:Real}Like resample2dpath but with retain_extrema=true and method=:linear as defaults
IMAS.is_z_offset — Function
is_z_offset(pr::Vector{T}, pz::Vector{T}; order::Int=4, precision::Float64=1E-3) where {T<:Real}Returns true the shape is offset from z=0 (ie. does not have mxh.z0=0)
is_z_offset(mxh::MXH; precision::Float64=1E-3)IMAS.is_updown_symmetric — Function
is_updown_symmetric(pr::Vector{T}, pz::Vector{T}; order::Int=4, precision::Float64=1E-3) where {T<:Real}Returns true if boundary is updown symmetric (independent of mxh.z0)
is_updown_symmetric(mxh::MXH; precision::Float64=1E-3)IMAS.minimum_distance_polygons_vertices — Function
minimum_distance_polygons_vertices(
R_obj1::AbstractVector{<:T},
Z_obj1::AbstractVector{<:T},
R_obj2::AbstractVector{<:T},
Z_obj2::AbstractVector{<:T};
return_index::Bool=false) where {T<:Real}Returns minimum distance between two polygons vertices and index of points on the two polygons
IMAS.minimum_distance_polygons — Function
minimum_distance_polygons(
R_obj1::AbstractVector{<:T},
Z_obj1::AbstractVector{<:T},
R_obj2::AbstractVector{<:T},
Z_obj2::AbstractVector{<:T}) where {T<:Real}Returns minimum distance between two polygons
NOTE: this is the actual distance, not the distance between the vertices
IMAS.min_mean_distance_polygons — Function
min_mean_distance_polygons(
R_obj1::AbstractVector{<:T},
Z_obj1::AbstractVector{<:T},
R_obj2::AbstractVector{<:T},
Z_obj2::AbstractVector{<:T}) where {T<:Real}Calculate the minimum and mean distances between two polygons in 2D space.
NOTE: this is the actual distance, not the distance between the vertices
IMAS.curvature — Function
curvature(pr::AbstractVector{T}, pz::AbstractVector{T}) where {T<:Real}Calculate the curvature of a 2D path defined by pr and pz using a finite difference approximation.
The path is assumed to be closed if the first and last points are the same, and open otherwise.
Arguments
pr: Real abstract vector representing the r-coordinates of the path.pz: Real abstract vector representing the z-coordinates of the path.
Returns
- A vector of the same length as
prandpzwith the calculated curvature values.
IMAS.angle_between_two_vectors — Function
angle_between_two_vectors(
v1_p1::Tuple{T,T},
v1_p2::Tuple{T,T},
v2_p1::Tuple{T,T},
v2_p2::Tuple{T,T}) where {T<:Real}Returns angle in radiants between two vectors defined by their start and end points
IMAS.bisector — Function
bisector(v1, v2, v3)Returns signed unitary bisector given three vertices
IMAS.polygon_rays — Function
polygon_rays(vertices::AbstractVector, extent_a::Float64, extent_b::Float64)Returns bisecting "rays" (lines) that radiate from vertices of a polygon
IMAS.split_long_segments — Function
split_long_segments(R::AbstractVector{T}, Z::AbstractVector{T}, max_length::Float64) where {T<:Real}Split long segments of a polygon so that each resulting segment is always <= max_length
split_long_segments(R::AbstractVector{T}, Z::AbstractVector{T}, n_points::Int) where {T<:Real}Split long segments of a polygon so that there are at least n_points in it
IMAS.thick_line_polygon — Function
thick_line_polygon(r1, z1, r2, z2, thickness1, thickness2)Generates a closed polygon from a thick line. Returns points of the quadrilateral polygon
thick_line_polygon(pr::AbstractVector{Float64}, pz::AbstractVector{Float64}, thickness::AbstractVector{Float64})IMAS.convex_hull — Function
convex_hull(xy_points::AbstractVector; closed_polygon::Bool)Compute the convex hull of a set of 2D points, sorted in counter-clockwise order. The resulting convex hull forms a closed polygon by appending the first point at the end.
convex_hull(x::AbstractVector{T}, y::AbstractVector{T}; closed_polygon::Bool) where {T}IMAS.convex_hull_indices — Function
convex_hull_indices(xy_points::AbstractVector)Compute the indices of points that form the convex hull of a set of 2D points. Returns the indices in counter-clockwise order.
This is the internal function that works with indices only.
Math
IMAS.norm01 — Function
norm01(x::T)::T where {T<:AbstractVector{<:Real}}Normalize a vector so that the first item in the array is 0 and the last one is 1 This is handy where psinorm should be used (and IMAS does not define a psinorm array)
IMAS.to_range — Function
to_range(vector::AbstractVector)Turn a vector into a range (if possible)
IMAS.meshgrid — Function
meshgrid(x::AbstractVector{T}, y::AbstractVector{T}) where {T}Return coordinate matrices from coordinate vectors
IMAS.calc_z — Function
calc_z(x::AbstractVector{<:Real}, f::AbstractVector{<:Real}, method::Symbol)Returns the gradient scale lengths of vector f on x
The finite difference method of the gradient can be one of [:backward, :central, :forward, :secondorder, :thirdorder]
NOTE: the inverse scale length is NEGATIVE for typical density/temperature profiles
IMAS.integ_z — Function
integ_z(rho::AbstractVector{<:Real}, z_profile::AbstractVector{<:Real}, bc::Real)Backward integration of inverse scale length vector with given edge boundary condition
IMAS.pack_grid_gradients — Function
pack_grid_gradients(x::AbstractVector{T}, y::AbstractVector{T}; n_points::Int=length(x), l::Float64=1E-2) where {T<:Float64}Returns grid between minimum(x) and maximum(x) with n_points points positioned to sample y(x) in such a way to pack more points where gradients are greates.
l controls how much the adaptive gradiant sampling should approach linear sampling.
IMAS.unique_indices — Function
unique_indices(vec::AbstractVector)::Vector{Int}Return the indices of the first occurrence of each unique element in the input vector vec
IMAS.index_circular — Function
index_circular(N::Int, idx::Int)If idx is beyond N or less than 1, it wraps around in a circular manner.
IMAS.getindex_circular — Function
getindex_circular(vec::AbstractVector{T}, idx::Int)::T where {T}Return the element of the vector vec at the position idx.
If idx is beyond the length of vec or less than 1, it wraps around in a circular manner.
IMAS.chunk_indices — Function
chunk_indices(dims::Tuple{Vararg{Int}}, N::Int)Split the indices of an array with dimensions dims into N chunks of similar size. Each chunk is a generator of CartesianIndex objects.
Arguments
dims::Tuple{Vararg{Int}}: A tuple specifying the dimensions of the array. For a 2D array, this would be(rows, cols).N::Int: The number of chunks to split the indices into.
Returns
- Vector with chunks of cartesian indices. Each chunk can be iterated over to get the individual
CartesianIndexobjects.
Outlines
IMAS.OutlineClosedVector — Type
OutlineClosedVector{T,A<:AbstractVector{T}} <: AbstractVector{T}Outlines in IMAS v3 ontology are always open (NOTE: in v4 they will always be closed!)
Use OutlineClosedVector to have a vector that behaves like a it's a closed outline, without creating a new array (allocation free)
When OutlineClosedVector is put into the dd, it will automatically save it there as an open outline
IMAS.OutlineOpenVector — Type
OutlineOpenVector{T,A<:AbstractVector{T}} <: AbstractVector{T}Outlines in IMAS v3 ontology are always open (NOTE: in v4 they will always be closed!)
Use OutlineOpenVector to have a vector that behaves like a it's a open outline, without creating a new array (allocation free)
When OutlineOpenVector is put into the dd, it will automatically save it there as an open outline
IMAS.CircularVector — Type
CircularVector{T,A<:AbstractVector{T}} <: AbstractVector{T}Vector with circular indexes, without creating a new array (allocation free)
NOTE: It can be combined with OutlineClosedVector and/or OutlineOpenVector, like this: CircularVector(OutlineClosedVector(x, x_is_closed))
IMAS.is_open_polygon — Function
is_open_polygon(R::AbstractVector{T}, Z::AbstractVector{T})::Bool where {T<:Real}Determine if a polygon, defined by separate vectors for R and Z coordinates, is open.
is_open_polygon(vertices::AbstractVector)::BoolIMAS.is_closed_polygon — Function
is_closed_polygon(R::AbstractVector{T}, Z::AbstractVector{T})::Bool where {T<:Real}Determine if a polygon, defined by separate vectors for R and Z coordinates, is closed.
is_closed_polygon(vertices::AbstractVector)::BoolIMAS.open_polygon — Function
open_polygon(R::AbstractVector{T}, Z::AbstractVector{T}, args...) where {T<:Real}Returns a view of the vectors R and Z such that they are a open polygon
Returns a named tuple containing the status of the polygon (wasclosed, wasopen) and the views of the R and Z vectors.
IMAS.closed_polygon — Function
closed_polygon(R::AbstractVector{T}, Z::AbstractVector{T}, args...) where {T<:Real}Returns a view of the vectors R and Z such that they are a closed polygon
Returns a named tuple containing the status of the polygon (wasclosed, wasopen) and the views of the R and Z vectors.
closed_polygon(R::AbstractVector{T}, Z::AbstractVector{T}, closed::Bool, args...) where {T<:Real}Returns a closed polygon depending on closed otherwise returns an open polygon
IMAS.is_clockwise — Function
is_clockwise(r::AbstractVector{T}, z::AbstractVector{T})::Bool where {T<:Real}Returns true/false if polygon is defined clockwise
IMAS.is_counterclockwise — Function
is_counterclockwise(r::AbstractVector{T}, z::AbstractVector{T})::Bool where {T<:Real}Returns true/false if polygon is defined counterclockwise
Physics boundary
IMAS.boundary_shape — Function
boundary_shape(;
a::T,
eps::T,
kapu::T,
kapl::T,
delu::T,
dell::T,
zetaou::T,
zetaiu::T,
zetail::T,
zetaol::T,
zoffset::T,
upnull::Bool=false,
lonull::Bool=false,
npts::Int=90
) where {T<:Real}Function used to generate boundary shapes based on T. C. Luce, PPCF, 55 9 (2013)
a: minor radiuseps: aspect ratiokapu: upper elongationlkap: lower elongationdelu: upper triangularitydell: lower triangularityzetaou: upper outer squarenesszetaiu: upper inner squarenesszetail: lower inner squarenesszetaol: lower outer squarenesszoffset: z-offsetupnull: toggle upper x-pointlonull: toggle lower x-pointnpts: number of points (per quadrant)
returns tuple with arrays of (r, z, zref)
>> boundary_shape(;a=0.608,eps=0.374,kapu=1.920,kapl=1.719,delu=0.769,dell=0.463,zetaou=-0.155,zetaiu=-0.255,zetail=-0.174,zetaol=-0.227,zoffset=0.000,upnull=true,lonull=false)IMAS.boundary — Function
boundary(pc::IMAS.pulse_schedule__position_control{T}, time0::Float64)return boundary from pulseschedule.positioncontrol at a given time0
boundary(pc::IMAS.pulse_schedule__position_control{T}, time_index::Int)returns boundary from pulseschedule.positioncontrol at a given time_index
boundary(pc::IMAS.pulse_schedule__position_control; time0::Float64=global_time(pc))Beturns r,z vectors from pulseschedule.positioncontrol.equilibriumtime_slice___boundaryoutline
IMAS.x_points — Function
x_points(x_points::IMAS.IDSvector{<:IMAS.pulse_schedule__position_control__x_point{T}}; time0::Float64=global_time(x_points)) where {T<:Real}Beturns vector with tuples of R,Z coordinates of x-points in pulse_schedule at time0
IMAS.strike_points — Function
strike_points(strike_points::IMAS.IDSvector{<:IMAS.pulse_schedule__position_control__strike_point{T}}; time0::Float64=global_time(strike_points)) where {T<:Real}Beturns vector with tuples of R,Z coordinates of x-points in pulse_schedule at time0
Physics build
IMAS.BuildLayerType — Type
Enum BuildLayerTypeUsed for dd.build.layer[:].type
_plasma_-> -1_gap_-> 0_oh_-> 1_tf_-> 2_shield_-> 3_blanket_-> 4_wall_-> 5_vessel_-> 6_cryostat_-> 7_divertor_-> 8_port_-> 9
IMAS.BuildLayerSide — Type
Enum BuildLayerSideUsed for dd.build.layer[:].side
_lfs_-> -1_lhfs_-> 0_hfs_-> 1_in_-> 2_out_-> 3
IMAS.BuildLayerShape — Type
Enum BuildLayerShapeUsed for dd.build.layer[:].shape
_offset_-> 0_negative_offset_-> 1_convex_hull_-> 2_mirror_princeton_D_exact_-> 3_princeton_D_-> 4_mirror_princeton_D_-> 5_princeton_D_scaled_-> 6_mirror_princeton_D_scaled_-> 7_rectangle_-> 8_double_ellipse_-> 9_mirror_double_ellipse_-> 10_rectangle_ellipse_-> 11_mirror_rectangle_ellipse_-> 12_circle_ellipse_-> 13_mirror_circle_ellipse_-> 14_triple_arc_-> 15_mirror_triple_arc_-> 16_miller_-> 17_silo_-> 18_racetrack_-> 19_undefined_-> 20
IMAS.build_radii — Function
build_radii(layers::IMAS.IDSvector{<:IMAS.build__layer{T}}) where {T<:Real}Convert thicknesses to absolute radii in the build layers
IMAS.get_build_indexes — Function
get_build_indexes(
layers::IMAS.IDSvector{<:IMAS.build__layer};
type::Union{Nothing,BuildLayerType}=nothing,
name::Union{Nothing,String}=nothing,
identifier::Union{Nothing,Integer}=nothing,
fs::Union{Nothing,BuildLayerSide,AbstractVector{BuildLayerSide}}=nothing)Returns indexes of layer(s) in build based on a series of selection criteria
IMAS.get_build_index — Function
get_build_index(
layers::IMAS.IDSvector{<:IMAS.build__layer};
type::Union{Nothing,BuildLayerType}=nothing,
name::Union{Nothing,String}=nothing,
identifier::Union{Nothing,Integer}=nothing,
fs::Union{Nothing,BuildLayerSide,AbstractVector{BuildLayerSide}}=nothing)Returns index of layer in build based on a series of selection criteria
It raises an error if none or more than one layer matches.
IMAS.get_build_layers — Function
get_build_layers(
layers::IMAS.IDSvector{<:IMAS.build__layer};
type::Union{Nothing,BuildLayerType}=nothing,
name::Union{Nothing,String}=nothing,
identifier::Union{Nothing,Integer}=nothing,
fs::Union{Nothing,BuildLayerSide,AbstractVector{BuildLayerSide}}=nothing)Select layer(s) in build based on a series of selection criteria
IMAS.get_build_layer — Function
get_build_layer(
layers::IMAS.IDSvector{<:IMAS.build__layer};
type::Union{Nothing,BuildLayerType}=nothing,
name::Union{Nothing,String}=nothing,
identifier::Union{Nothing,Integer}=nothing,
fs::Union{Nothing,BuildLayerSide,AbstractVector{BuildLayerSide}}=nothing)Select layer in build based on a series of selection criteria
It raises an error if none or more than one layer matches.
IMAS.opposite_side_layer — Function
opposite_side_layer(layer::IMAS.build__layer)Returns corresponding layer on the high/low field side
IMAS.func_nested_layers — Function
func_nested_layers(layer::IMAS.build__layer{D}, func::Function)::D where {D<:Real}Apply function func to a layer, then subtract func applied to the layer inside of it.
This is used to caclulate area, volume, etc.. of each layer.
IMAS.volume — Function
volume(layer::IMAS.build__layer)Calculate volume of a build layer revolved around z axis
volume(layer::IMAS.build__structure)Calculate volume of a build structure outline revolved around z axis
volume(coil::IMAS.pf_active__coil)Returns volume of PF coils
IMAS.first_wall — Function
first_wall(wall::IMAS.wall{T}; simplify_to_inscribed_fractional_area::Float64=1 - 1E-6) where {T<:Real}Returns named tuple with outline of the official contiguous first wall limiter contour, or an empty outline if not present
NOTE: in IMAS wall.description_2d[].limiter.type.index == 0 indicates an official contiguous limiter contour
NOTE: By default the first wall simplifies points on straight lines
first_wall(pf_active::IMAS.pf_active{T}) where {T<:Real}Returns named tuple with outline of the first wall, defined by the pf_active.coils
first_wall(eqt::IMAS.equilibrium__time_slice; precision::Float=1E-3) where {T<:Real}Returns named tuple with outline of the first wall, defined by equilibrium computation domain
first_wall(eqt::IMAS.equilibrium__time_slice, pf_active::IMAS.pf_active{T}) where {T<:Real}Returns named tuple with outline of the first wall, defined as the equilibrium computational domain with cutouts for the pf_active.coils that fall in it
IMAS.first_wall! — Function
first_wall!(wall::IMAS.wall{T}, r::AbstractVector{T}, z::AbstractVector{T}) where {T<:Real}Set wall.description_2d[?].limiter.unit[1].outline from input r and z
Returns the limiter.unit with the new outline
IMAS.contiguous_limiter_from_open_limiters! — Function
contiguous_limiter_from_open_limiters!(wall::IMAS.wall{T}) where {T<:Real}Combines multiple open limiter units into a single contiguous closed limiter by connecting their outlines in sequence. The function chooses the connection order to minimize gaps between adjacent limiter segments.
Modifies the wall structure in-place by replacing existing limiters with one closed contiguous limiter named "contiguous first wall" and limiter.type.index=0
Returns the limiter.unit with the contiguous outline.
IMAS.build_max_R0_B0 — Function
build_max_R0_B0(bd::IMAS.build)Returns the plasma geometric center (r0) and the maximum vacuum toroidal magnetic field (b0) evaluated at (r0) that the TF build allows
IMAS.vertical_maintenance — Function
vertical_maintenance(bd::IMAS.build; tor_modularity::Int=2, pol_modularity::Int=1, gap_VV_BL::Float64=0.1)Returns the radial dimensions of the vertical vacuum port for blanket maintenance
gapVVBL is the margin between the blanket module and the wall (10cm seems reasonable)
IMAS.outline — Function
outline(layer::Union{IMAS.build__layer, IMAS.build__structure})
outline(out::Union{IMAS.build__layer___outline,IMAS.build__structure___outline})Returns a closed polygon as a named tuple with (r,z) of a dd.build.layer or dd.build.structure
outline(element::Union{IMAS.pf_active__coil___element{T},IMAS.pf_passive__loop___element{T}}) where {T<:Real}Returns named tuple with r and z arrays outline, independent of geometry_type used to describe the element
Physics collisions
IMAS.lnΛ_ee — Function
lnΛ_ee(ne::Real, Te::Real)Calculate Couloumb logarithm (lnΛ) for thermal electron-electron collisions [NRL Plasma Formulary]
ne: electron density [m^-3]Te: electron temperature [eV]
IMAS.lnΛ_ei — Function
nΛ_ei(ne::S, Te::P, ni::AbstractVector{Q}, Ti::AbstractVector{R}, mi::AbstractVector{T}, Zi::AbstractVector{Int}) where {S<:Real,P<:Real,Q<:Real,R<:Real,T<:Real}Calculate Couloumb logarithm (lnΛ) for thermal electron-ion collisions [NRL Plasma Formulary]
ne: electron density [m^-3]Te: electron temperature [eV]ni: list of ion densities [m^-3]Ti: list of ion temperaturs [eV]mi: list of ion masses [amu]Zi: list of ion charges
lnΛ_ei(ne::Real, Te::Real)Calculate Couloumb logarithm (lnΛ) for thermal electron-ion collisions where Ti*me/mi < 10Zi^2 eV < Te [NRL Plasma Formulary]
ne: electron density [m^-3]Te: electron temperature [eV]
IMAS.lnΛ_ii — Function
lnΛ_ii(ne::Real, Te::Real, ni1::Real, Ti1::Real, mi1::Real, Zi1::Real, ni2::Real, Ti2::Real, mi2::Real, Zi2::Real; beta_D::Real=0.0)Calculate Couloumb logarithm (lnΛ) for mixed thermal ion1-ion2 collisions [NRL Plasma Formulary]
ne: electron density [m^-3]Te: electron temperaturs [eV]ni1: ion1 density [m^-3]Ti1: ion1 temperature [eV]mi1: ion1 mass [amu]Zi1: ion1 chargeni1: ion2 density [m^-3]Ti1: ion2 temperature [eV]mi1: ion2 mass [amu]Zi1: ion2 chargebeta_D: relative drift velocities between ion speciesv_D = beta_D*c
IMAS.lnΛ_fi — Function
lnΛ_fi(ne::Real, Te::Real, n_i::Real, T_i::Real, m_i::Real, Z_i::Real, beta_f::Real, mf::Real, Zf::Real; verbose=true)Calculate Couloumb logarithm (lnΛ) for beam/fast ion in the presence of warm electrons/ions [NRL Plasma Formulary]
ne: electron density [m^-3]Te: electron temperaturs [eV]n_i: ion density [m^-3]Ti: ion temperature [eV]mi: ion mass [amu]Zi: ion chargebeta_f: relative fast ion velocityv_f = beta_f*cmf: mass of fast ion [amu]Zf: charge of fast ion
Physics constants
IMAS.mks — Constant
Named tuple with physics constants in mks (and eV):
μ_0 = 1.25663706212e-6 [N A^-2] # Vacuum permeability
c = 2.99792458e8 [m s^-1] # Speed of light in vacuum
ϵ_0 = 8.8541878128e-12 [F m^-1] # Vacuum permittivity
k_B = 1.380649e-23 [J K^-1] # Boltzmann constant
e = 1.602176634e-19 [C] # Elementary charge
m_e = 9.1093837015e-31 [kg] # Electron mass
m_p = 1.67262192369e-27 [kg] # Proton mass
m_n = 1.67492749804e-27 [kg] # Neutron mass
m_d = 3.3435837768e-27 [kg] # Deuteron mass
atm = 101325.0 [Pa] # Standard atmosphere
m_u = 1.6605390666e-27 [kg] # Atomic mass constant
avog = 6.02214076e23 [mol^-1] # Avogadro constant
E_α = 3.518e6 [eV] # Alpha particle energy
E_n = 14.072e6 [eV] # Neutron energyPhysics currents
IMAS.j_ohmic_steady_state — Function
j_ohmic_steady_state(eqt::IMAS.equilibrium__time_slice{T}, cp1d::IMAS.core_profiles__profiles_1d{T}, ip::T, j_ohmic_shape::AbstractVector{T}=cp1d.conductivity_parallel) where {T<:Real}Sets j_ohmic parallel current density to what it would be at steady-state, based on conductivity_parallel, j_non_inductive and a target total ip
Requires constant loop voltage: Vl = 2π * η * <J_oh⋅B> / (F * <R⁻²>) = constant
IMAS.JtoR_2_JparB — Function
JtoR_2_JparB(rho_tor_norm::Vector{<:Real}, JtoR::Vector{<:Real}, includes_bootstrap::Bool, eqt::IMAS.equilibrium__time_slice)Given <Jt/R> returns <J⋅B>
Transformation obeys <J⋅B> = (1/f)*(<B^2>/<1/R^2>)*(<Jt/R> + dp/dpsi*(1 - f^2*<1/R^2>/<B^2>))
Includes_bootstrap set to true if input current includes bootstrap
NOTE: Jtor ≂̸ JtoR
<Jt/R> = <Jt/R>/<1/R> * <1/R> = Jtor * <1/R> = Jtor * gm9NOTE: Jpar ≂̸ JparB
JparB = Jpar * B0IMAS.JparB_2_JtoR — Function
JparB_2_JtoR(rho_tor_norm::Vector{<:Real}, JparB::Vector{<:Real}, includes_bootstrap::Bool, eqt::IMAS.equilibrium__time_slice)Given <J⋅B> returns <Jt/R>
Transformation obeys <J⋅B> = (1/f)*(<B^2>/<1/R^2>)*(<Jt/R> + dp/dpsi*(1 - f^2*<1/R^2>/<B^2>))
Includes_bootstrap set to true if input current includes bootstrap
NOTE: Jtor ≂̸ JtoR
<Jt/R> = <Jt/R>/<1/R> * <1/R> = Jtor * <1/R> = Jtor * gm9NOTE: Jpar ≂̸ JparB
JparB = Jpar * B0IMAS.Jtor_2_Jpar — Function
Jtor_2_Jpar(rho_tor_norm::Vector{<:Real}, Jtor::Vector{<:Real}, includes_bootstrap::Bool, eqt::IMAS.equilibrium__time_slice)Given Jtor returns Jpar
Includes_bootstrap set to true if input current includes bootstrap
NOTE: Jtor ≂̸ JtoR
<Jt/R> = <Jt/R>/<1/R> * <1/R> = Jtor * <1/R> = Jtor * gm9NOTE: Jpar ≂̸ JparB
JparB = Jpar * B0IMAS.Jpar_2_Jtor — Function
Jpar_2_Jtor(rho_tor_norm::Vector{<:Real}, Jpar::Vector{<:Real}, includes_bootstrap::Bool, eqt::IMAS.equilibrium__time_slice)Given Jpar returns Jtor
Includes_bootstrap set to true if input current includes bootstrap
NOTE: Jtor ≂̸ JtoR
<Jt/R> = <Jt/R>/<1/R> * <1/R> = Jtor * <1/R> = Jtor * gm9NOTE: Jpar ≂̸ JparB
JparB = Jpar * B0IMAS.vloop — Function
vloop(cp1d::IMAS.core_profiles__profiles_1d{T}) where {T<:Real}Vloop = 2π * η * <J_oh⋅B> / (F * <R⁻²>): method emphasizes the resistive nature of the plasma
vloop(eq::IMAS.equilibrium{T}, n::Int=1; time0::Float64=global_time(eq)) where {T<:Real}Vloop = dψ/dt: method emphasizes the inductive nature of the loop voltage. Assumes COCOS 11.
n is the number of timeslices priors use to take the difference in psiboundary
vloop(eq::IMAS.equilibrium{T}, n::Int; time0::Float64=global_time(eq)) where {T<:Real}Vloop = dψ/dt: method emphasizes the inductive nature of the loop voltage. Assumes COCOS 11.
δt is the time difference used to take the difference in psi_boundary
vloop(ct::IMAS.controllers{T}; time0::Float64=global_time(ct)) where {T<:Real}Returns vloop at time0 from controller named ip
IMAS.vloop_time — Function
vloop_time(ct::IMAS.controllers{T}) where {T<:Real}Returns named tuple with time and data with vloop from controller named ip
vloop_time(eq::IMAS.equilibrium{T}) where {T<:Real}Returns named tuple with time and data with vloop from equilibrium
vloop_time(cp::IMAS.core_profiles{T}, eq::IMAS.equilibrium) where {T<:Real}Returns named tuple with time and data with vloop from core_profiles
IMAS.Ip_non_inductive — Function
Ip_non_inductive(cp1d::IMAS.core_profiles__profiles_1d{T}, eqt::IMAS.equilibrium__time_slice{T}) where {T<:Real}Integrated toroidal non-inductive current
IMAS.Ip_bootstrap — Function
Ip_bootstrap(cp1d::IMAS.core_profiles__profiles_1d{T}, eqt::IMAS.equilibrium__time_slice{T}) where {T<:Real}Integrated toroidal bootstrap current
IMAS.Ip_ohmic — Function
Ip_ohmic(cp1d::IMAS.core_profiles__profiles_1d{T}, eqt::IMAS.equilibrium__time_slice{T}) where {T<:Real}Integrated toroidal ohmic current
IMAS.plasma_lumped_resistance — Function
plasma_lumped_resistance(dd::IMAS.dd)Returns equivalent plasma lumped resistance in ohms
Physics diagnostics
IMAS.extract_subsystem_from_channel_name — Function
extract_subsystem_from_channel_name(name::AbstractString)Extract subsystem identifier from channel name. Returns the subsystem string, or the full name if no pattern is recognized.
Handles various naming conventions:
- Position-based: "TScorer+038" → "TScore"
- Position-based: "TStangentialr+125" → "TStangential"
- Simple numbered: "TSCORE01" → "TS_CORE"
- Generic pattern: "DIVERTOR_05" → "DIVERTOR"
Physics equilibrium
IMAS.ψ_interpolant — Function
ψ_interpolant(eqt2d::IMAS.equilibrium__time_slice___profiles_2d)Returns r, z, and ψ interpolant named tuple
ψ_interpolant(r::AbstractRange{T1}, z::AbstractRange{T1}, psi::Matrix{T2}) where {T1<:Real,T2<:Real}ψ_interpolant(eqt2dv::IDSvector{<:IMAS.equilibrium__time_slice___profiles_2d})ψ_interpolant(eqt::IMAS.equilibrium__time_slice)IMAS.ρ_interpolant — Function
ρ_interpolant(eqt2d::IMAS.equilibrium__time_slice___profiles_2d{T}, phi_norm::T) where {T<:Real}Returns r, z, and ρ interpolant named tuple
ρ_interpolant(r::AbstractRange{T}, z::AbstractRange{T}, rho::Matrix{T}) where {T<:Real}ρ_interpolant(eqt2dv::IDSvector{<:IMAS.equilibrium__time_slice___profiles_2d{T}}, phi_norm::T) where {T<:Real}ρ_interpolant(eqt::IMAS.equilibrium__time_slice)IMAS.calc_pprime_ffprim_f — Function
calc_pprime_ffprim_f(
psi::T,
R::T,
one_R::T,
one_R2::T,
R0::Real,
B0::Real;
pressure::Union{Missing,T}=missing,
dpressure_dpsi::Union{Missing,T}=missing,
j_tor::Union{Missing,T}=missing,
j_tor_over_R::Union{Missing,T}=missing,
f::Union{Missing,T}=missing,
f_df_dpsi::Union{Missing,T}=missing,
) where {T<:AbstractVector{<:Real}}This method returns the P' and FF' given P or P' and J or J/R based on the current equilibrium fluxsurfaces geometry
Arguments:
psi: poloidal fluxR: <R>one_R: <1/R>one_R2: <1/R²>R0: R at which B0 is definedB0: vacuum magnetic fieldpressure: Pdpressure_dpsi: P'j_tor: toroidal currentj_tor_over_R: flux surface averaged toroidal current density over major radiusf: Ff_df_dpsi: FF'
returns (P', FF', F)
IMAS.p_jtor_2_pprime_ffprim_f — Function
p_jtor_2_pprime_ffprim_f(eqt1::IMAS.equilibrium__time_slice___profiles_1d, R0::Real, B0::Real)Calculate P', FF' and F from pressure and j_tor in equilibrium
IMAS.symmetrize_equilibrium! — Function
symmetrize_equilibrium!(eqt::IMAS.equilibrium__time_slice)Update equilibrium time slice in place to be symmetric with respect to its magnetic axis.
This is done by averaging the upper and lower parts of the equilibrium.
Flux surfaces should re-traced after this operation.
NOTE: Use with care! This operation will change the flux surfaces (LCFS included) and as such quantities may change
IMAS.B0_geo — Function
B0_geo(eqt::IMAS.equilibrium__time_slice{T})::T where{T<:Real}Returns vacuum B0 at the plasma geometric center, which may or may not be eqt.global_quantities.vacuum_toroidal_field.r0
IMAS.elongation_limit — Function
elongation_limit(R0_over_a::Real)Returns elongation limit due to control limit from simple aspect ratio scaling
elongation_limit(eqt::IMAS.equilibrium__time_slice)IMAS.optimal_kappa_delta — Function
optimal_kappa_delta(li::T1, βp::T1, ϵ::T1, γτw::T2, ∆o::T2) where {T1<:Real,T2<:Real}An analytic scaling relation for the maximum tokamak elongation against n=0 MHD resistive wall modes Jungpyo Lee, Jeffrey P. Freidberg, Antoine J. Cerfon, Martin Greenwald https://doi.org/10.1088%2F1741-4326%2Faa6877
NOTE:
γτwis the feedback capability parameter and represents how fast a instability is controllable (𝛾is the instability growth rate andτwis the wall diffusion time) (typicallyγτw < 10is assumed for controllability, seeVacuumFields.normalized_growth_rate())∆ois the outer gap (NOTE: assumes∆o = ∆i = 1/3 * ∆v) detemines the relation betweenκandδof the plasma boundary and theκw=(κ+3∆o)(1+∆o)andδw=δ(1+∆o)of the wall boundary
optimal_kappa_delta(eqt::IMAS.equilibrium__time_slice, γτw::T, ∆o::T) where {T<:Real}Physics fast
IMAS.estrada_I_integrals — Function
estrada_I_integrals(ne::Real, Te::Real, ni::AbstractVector{<:Real}, Ti::AbstractVector{<:Real}, mi::AbstractVector{<:Real}, Zi::AbstractVector{Int}, Ef::Real, mf::Real, Zf::Real)Returns solution to i2 and i4 integrals from [Estrada et al., Phys of Plasm. 13, 112303 (2006)] Eq. 9 & 10
ne: electron density [m^-3]Te: electron temperature [eV]ni: list of ion densities [m^-3]Ti: list of ion temperatures [eV]mi: list of ion masses [amu]Zi: list of ion chargesEf: fast ion energy [eV]mf: mass of fast ion [amu]Zf: fast ion charge
estrada_I_integrals(ne::Real, Te::Real, ni::AbstractVector{<:Real}, Ti::AbstractVector{<:Real}, mi::AbstractVector{<:Real}, Zi::AbstractVector{Int}, Ef::Real, mf::Real, Zf::Real)Returns solution to i2 and i4 integrals from [Estrada et al., Phys of Plasm. 13, 112303 (2006)] Eq. 9 & 10
Ef: fast ion energy [eV]Ec: critical energy [eV]
IMAS.α_slowing_down_time — Function
α_slowing_down_time(cp1d::IMAS.core_profiles__profiles_1d, irho::Int=1)Returns the slowing down time in seconds of α particles evaluated at irho (by default on axis)
IMAS.critical_energy — Function
critical_energy(cp1d::IMAS.core_profiles__profiles_1d, irho::Int, mf::Real, Zf::Real; approximate::Bool=true)Returns Ec the critical energy by finding the root of the difference between the electron and ion drag
mf: mass of fast ion [amu]Zf: fast ion chargeapproximate: calculate critical energy assuminglnΛ_fe == lnΛ_fi. For DIII-D this results in a correction factor of (lnΛfi/lnΛfe)^(2/3) ≈ 1.2.
IMAS.thermalization_time — Function
thermalization_time(v_f, v_c, tau_s)Calculate thermalization time in seconds
v_f: fast ion velocityv_c: critical velocitytau_s: slowing down time
thermalization_time(cp1d::IMAS.core_profiles__profiles_1d, irho::Int, Ef::Real, mf::Real, Zf::Real)Calculate thermalization time in seconds of a fast ion with energy Ef and Ti*me/mi < 10Zi^2 eV < Te
Ef: fast ion energy [eV]mf: mass of fast ion [amu]Zf: fast ion charge
IMAS.α_thermalization_time — Function
α_thermalization_time(cp1d::IMAS.core_profiles__profiles_1d, irho::Int)Returns the thermalization time in seconds of α particles evaluated on axis
IMAS.fast_ion_thermalization_time — Function
fast_ion_thermalization_time(cp1d::IMAS.core_profiles__profiles_1d, irho::Int, ion::Union{IMAS.nbi__unit___species,IMAS.core_profiles__profiles_1d___ion___element}, ion_energy::Real)Returns the fast ion thermalization time in seconds evaluated on axis
IMAS.fast_particles_profiles! — Function
fast_particles_profiles!(cs::IMAS.core_sources, cp1d::IMAS.core_profiles__profiles_1d; verbose::Bool=false)Fills the core_profiles fast ion densities and pressures that result from fast ion sources (eg. fusion and nbi)
This calculation is done based on the slowing_down_time and thermalization_time of the fast ion species.
IMAS.sivukhin_fraction — Function
sivukhin_fraction(cp1d::IMAS.core_profiles__profiles_1d, particle_energy::Real, particle_mass::Real)Compute a low-accuracy but fast approximation of the ion to electron heating fraction for fast particles (like alpha particles and beam particles)
IMAS.smooth_beam_power — Function
smooth_beam_power(time::AbstractVector{Float64}, power::AbstractVector{T}, taus::Float64) where {T<:Real}Smooths out the beam power history based on a given thermalization constant taus
Such smoothing mimics the delayed contribution from the instantaneous source, as well as the gradual decay of the previous fast ion population over time.
smooth_beam_power(time::AbstractVector{Float64}, power::AbstractVector{T}, time0::Float64, taus::Float64) where {T<:Real}Return smoothed beam power at time0
Physics fields
IMAS.Br_Bz — Function
Br_Bz(eqt2d::IMAS.equilibrium__time_slice___profiles_2d)Returns Br and Bz named tuple evaluated at r and z starting from ψ interpolant
Br_Bz(eqt2dv::IDSvector{<:IMAS.equilibrium__time_slice___profiles_2d})Br_Bz(PSI_interpolant::Interpolations.AbstractInterpolation, r::T, z::T) where {T<:Real}Br_Bz(PSI_interpolant::Interpolations.AbstractInterpolation, r::AbstractArray{T}, z::AbstractArray{T}) where {T<:Real}IMAS.Br_Bphi_Bz — Function
Br_Bphi_Bz(PSI_interpolant::Interpolations.AbstractInterpolation, B0::T, R0::T, r::T, z::T) where {T<:Real}Returns Br, Bphi, Bz named tuple evaluated at r and z starting from ψ interpolant and B0 and R0
IMAS.Bx_By_Bz — Function
Bx_By_Bz(PSI_interpolant::Interpolations.AbstractInterpolation, B0::T, R0::T, r::T, z::T) where {T<:Real}Returns Bx, By, Bz named tuple evaluated at x,y,z starting from ψ interpolant and B0 and R0
IMAS.Bp — Function
Bp(eqt2d::IMAS.equilibrium__time_slice___profiles_2d)Returns Bp evaluated at r and z starting from ψ interpolant
Bp(PSI_interpolant::Interpolations.AbstractInterpolation, r::T, z::T) where {T<:Real}Bp(eqt2dv::IDSvector{<:IMAS.equilibrium__time_slice___profiles_2d})IMAS.trace_field_line — Function
trace_field_line(vector_field, start_point, step_size, stop_condition)Compute a trajectory along a vector field using the implicit midpoint method.
Arguments
vector_field: Vector field functionstart_point: Initial point of the trajectorystep_size: Size of each integration stepstop_condition: Function determining when to stop integration. Takes a ImplicitMidpointState struct as an argument.
Returns
Trajectory represented as a vector of points
trace_field_line(eqt::IMAS.equilibrium__time_slice, r, z; phi=zero(r), step_size=0.01, max_turns=1, stop_condition=nothing)Compute a field line trajectory for a specific equilibrium time slice.
Arguments
eqt: Equilibrium time slicer: Radial coordinatez: Vertical coordinatephi: Toroidal angle (default: 0)step_size: Size of each integration step (default: 0.01)max_turns: Maximum number of toroidal turns used in default stop condition (default: 1)stop_condition: User defined stop condition that takes a ImplicitMidpointState struct and returns a Boolean (default: nothing)
Returns
Named tuple with x, y, z, r, and phi coordinates of the trajectory
Physics flux-surfaces
IMAS.find_psi_boundary — Function
find_psi_boundary(
dimR::Union{AbstractVector{T1},AbstractRange{T1}},
dimZ::Union{AbstractVector{T1},AbstractRange{T1}},
PSI::Matrix{T2},
psi_axis::Real,
original_psi_boundary::Real,
RA::T3,
ZA::T3,
fw_r::AbstractVector{T4},
fw_z::AbstractVector{T5};
PSI_interpolant=IMAS.ψ_interpolant(dimR, dimZ, PSI).PSI_interpolant,
precision::Float64=flux_surfaces_precision,
raise_error_on_not_open::Bool,
raise_error_on_not_closed::Bool,
verbose::Bool=false
) where {T1<:Real,T2<:Real,T3<:Real,T4<:Real,T5<:Real}find_psi_boundary(
dimR::Union{AbstractVector{T1},AbstractRange{T1}},
dimZ::Union{AbstractVector{T1},AbstractRange{T1}},
PSI::Matrix{T2},
psi_axis::Real,
axis2bnd::Symbol,
RA::T3,
ZA::T3,
fw_r::AbstractVector{T4}=T1[],
fw_z::AbstractVector{T5}=T1[],
r_cache::AbstractVector{T1}=T1[],
z_cache::AbstractVector{T1}=T1[];
PSI_interpolant=IMAS.ψ_interpolant(dimR, dimZ, PSI).PSI_interpolant,
precision::Float64=flux_surfaces_precision,
raise_error_on_not_open::Bool,
raise_error_on_not_closed::Bool,
verbose::Bool=false
) where {T1<:Real,T2<:Real,T3<:Real,T4<:Real,T5<:Real}find_psi_boundary(
dimR::Union{AbstractVector{T1},AbstractRange{T1}},
dimZ::Union{AbstractVector{T1},AbstractRange{T1}},
PSI::Matrix{T2},
psirange_init::AbstractVector,
RA::T3,
ZA::T3,
fw_r::AbstractVector{T4},
fw_z::AbstractVector{T5},
r_cache::AbstractVector{T1}=T1[],
z_cache::AbstractVector{T1}=T1[];
PSI_interpolant=IMAS.ψ_interpolant(dimR, dimZ, PSI).PSI_interpolant,
precision::Float64=flux_surfaces_precision,
raise_error_on_not_open::Bool,
raise_error_on_not_closed::Bool,
verbose::Bool=false
)IMAS.find_psi_separatrix — Function
find_psi_separatrix(eqt::IMAS.equilibrium__time_slice{T}; precision::Float64=flux_surfaces_precision) where {T<:Real}Returns psi of the first magentic separatrix
Note: The first separatrix is the LCFS only in diverted plasmas
IMAS.find_psi_2nd_separatrix — Function
find_psi_2nd_separatrix(eqt::IMAS.equilibrium__time_slice{T}; precision::Float64=flux_surfaces_precision) where {T<:Real}Returns psi of the second magentic separatrix
returns diverted and not_diverted surface in a NamedTuple
divertedis when the surface starts and finishes on the same side of the midplanenot_divertedis when the surface starts and finishes on opposite sides of the midplane
IMAS.find_psi_last_diverted — Function
find_psi_last_diverted(
eqt::IMAS.equilibrium__time_slice,
wall_r::AbstractVector{<:Real},
wall_z::AbstractVector{<:Real},
PSI_interpolant::Interpolations.AbstractInterpolation;
precision::Float64=flux_surfaces_precision)Returns psi_last_lfs,psifirstlfsfar, andnullwithin_wall`
psi_first_lfs_farwill be the first surface insideOFL[:lfs_far]psi_last_lfswill be the last surface insideOFL[:lfs]
Precision between the two is defined on the poloidal crossection area at the OMP (Psol*precision = power flowing between psi_first_lfs_far and psi_last_lfs ~ 0)
IMAS.find_psi_tangent_omp — Function
find_psi_tangent_omp(
eqt::IMAS.equilibrium__time_slice,
wall_r::AbstractVector{<:Real},
wall_z::AbstractVector{<:Real},
PSI_interpolant::Interpolations.AbstractInterpolation;
precision::Float64=flux_surfaces_precision)Returns the psi of the magnetic surface in the SOL which is tangent to the wall near the outer midplane
IMAS.find_psi_max — Function
find_psi_max(eqt::IMAS.equilibrium__time_slice{T}; precision::Float64=1e-2) where {T<:Real}Returns the max psi useful for an ofl in the SOL with no wall.
IMAS.find_psi_wall_omp — Function
find_psi_wall_omp(eqt::IMAS.equilibrium__time_slice, wall_r::AbstractVector{<:Real}, wall_z::AbstractVector{<:Real})Returns the psi of the magnetic surface in the SOL which intersects the wall at the outer midplane
find_psi_wall_omp(
PSI_interpolant::Interpolations.AbstractInterpolation,
RA::T1,
ZA::T1,
wall_r::AbstractVector{T2},
wall_z::AbstractVector{T2},
psi_max::T1,
psi_sign::T1
) where {T1<:Real,T2<:Real}IMAS.interp_rmid_at_psi — Function
interp_rmid_at_psi(eqt::IMAS.equilibrium__time_slice, PSI_interpolant::Interpolations.AbstractInterpolation, R::AbstractVector{<:Real})Returns the interpolant r_mid(ψ) to compute the r at the midplane of the flux surface identified by ψ
The vector R defines the sampling of interest for thie interpolation
IMAS.SimpleSurface — Type
struct SimpleSurface{T<:Real} <: AbstractFluxSurface{T}
psi::T
r::Vector{T}
z::Vector{T}
ll::Vector{T}
fluxexpansion::Vector{T}
int_fluxexpansion_dl::T
endA simplified version of FluxSurface that only has the contour points and what is needed to compute flux surface averages
IMAS.FluxSurface — Type
struct FluxSurface{T<:Real} <: AbstractFluxSurface{T}
psi::T
r::Vector{T}
z::Vector{T}
r_at_max_z::T
max_z::T
r_at_min_z::T
min_z::T
z_at_max_r::T
max_r::T
z_at_min_r::T
min_r::T
Br::Vector{T}
Bz::Vector{T}
Bp::Vector{T}
Btot::Vector{T}
ll::Vector{T}
fluxexpansion::Vector{T}
int_fluxexpansion_dl::T
endIMAS.trace_simple_surfaces — Function
trace_simple_surfaces(eqt::IMAS.equilibrium__time_slice{T}, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}) where {T<:Real}Trace flux surfaces and returns vector of SimpleSurface structures. The result contains only the contours and what is needed to perform flux-surface averaging.
trace_simple_surfaces(
psi::AbstractVector{T},
r::AbstractVector{T},
z::AbstractVector{T},
PSI::Matrix{T},
PSI_interpolant::Interpolations.AbstractInterpolation,
RA::T,
ZA::T,
wall_r::AbstractVector{T},
wall_z::AbstractVector{T}
) where {T<:Real}Trace flux surfaces and returns vector of SimpleSurface structures. The result contains only the contours and what is needed to perform flux-surface averaging.
IMAS.trace_simple_surfaces! — Function
trace_simple_surfaces!(
surfaces::Vector{SimpleSurface{T}},
psi::AbstractVector{T},
r::AbstractVector{T},
z::AbstractVector{T},
PSI::Matrix{T},
PSI_interpolant::Interpolations.AbstractInterpolation,
RA::T,
ZA::T,
wall_r::AbstractVector{T},
wall_z::AbstractVector{T},
r_cache::AbstractVector{T}=T[],
z_cache::AbstractVector{T}=T[])
) where {T<:Real}Trace flux surfaces and store in surfaces vector of SimpleSurface structures. The result contains only the contours and what is needed to perform flux-surface averaging.
IMAS.trace_surfaces — Function
trace_surfaces(eqt::IMAS.equilibrium__time_slice{T}, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}; refine_extrema::Bool=true) where {T<:Real}Trace flux surfaces and returns vector of FluxSurface structures
trace_surfaces(
psi::AbstractVector{T},
f::AbstractVector{T},
r::AbstractVector{T},
z::AbstractVector{T},
PSI::Matrix{T},
BR::Matrix{T},
BZ::Matrix{T},
PSI_interpolant::Interpolations.AbstractInterpolation,
RA::T,
ZA::T,
wall_r::AbstractVector{T},
wall_z::AbstractVector{T};
refine_extrema::Bool=true
) where {T<:Real}IMAS.flux_surfaces — Function
flux_surfaces(eq::equilibrium{T}, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}) where {T<:Real}Update flux surface averaged and geometric quantities for all time_slices in the equilibrium IDS
flux_surfaces(eqt::equilibrium__time_slice{T1}, wall_r::AbstractVector{T2}, wall_z::AbstractVector{T2}) where {T1<:Real,T2<:Real}Update flux surface averaged and geometric quantities for a given equilibrum IDS time slice.
flux_surfaces(eqt::equilibrium__time_slice, wall::IMAS.wall)IMAS.flux_surface — Function
flux_surface(eqt::equilibrium__time_slice{T}, psi_level::Real, type::Symbol, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}) where {T<:Real}Returns a vector with the (r,z) coordiates of flux surface at given psi_level
The type parameter:
:any: return all contours:closed: all closed flux-surface that encircle the magnetic axis and do not cross the wall:open: all open flux-surfaces (considerning open even closed flux surfaces that hit the first wall):open_no_wall: all open flux-surfaces independently of wall:encircling: open flux-surfaces encircling the magnetic axis
flux_surface(
dim1::AbstractVector{T1},
dim2::AbstractVector{T1},
PSI::AbstractArray{T2},
RA::T3,
ZA::T3,
fw_r::AbstractVector{T4},
fw_z::AbstractVector{T4},
psi_level::T5,
type::Symbol
) where {T1<:Real,T2<:Real,T3<:Real,T4<:Real,T5<:Real}flux_surface(
dim1::Union{AbstractVector{T},AbstractRange{T}},
dim2::Union{AbstractVector{T},AbstractRange{T}},
cl::Contour.ContourLevel,
RA::T,
ZA::T,
fw_r::AbstractVector{T},
fw_z::AbstractVector{T},
type::Symbol) where {T<:Real}IMAS.flux_surface_avg — Function
flux_surface_avg(quantity::AbstractVector{T}, surface::FluxSurface{T}) where {T<:Real}Flux surface averaging of a quantity
flux_surface_avg(f::F1, surface::FluxSurface{T}) where {F1<:Function, T<:Real}Flux surface averaging of a function
IMAS.volume_integrate — Function
volume_integrate(eqt::IMAS.equilibrium__time_slice, what::AbstractVector{T})::AbstractVector{T} where {T<:Real}Integrate quantity over volume
IMAS.cumlul_volume_integrate — Function
volume_integrate(eqt::IMAS.equilibrium__time_slice, what::AbstractVector{T})::AbstractVector{T} where {T<:Real}Cumulative integrate quantity over volume
IMAS.surface_integrate — Function
surface_integrate(eqt::IMAS.equilibrium__time_slice, what::AbstractVector{T})::AbstractVector{T} where {T<:Real}Integrate quantity over surface
IMAS.cumlul_surface_integrate — Function
surface_integrate(eqt::IMAS.equilibrium__time_slice, what::AbstractVector{T})::AbstractVector{T} where {T<:Real}Cumulative integrate quantity over surface
IMAS.find_x_point! — Function
find_x_point!(eqt::IMAS.equilibrium__time_slice{T}, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}) where {T<:Real}Find the n X-points that are closest to the separatrix
IMAS.x_points_inside_wall — Function
x_points_inside_wall(x_points::IDSvector{<:IMAS.equilibrium__time_slice___boundary__x_point}, wall::IMAS.wall)Returns vector of x_points that are inside of the first wall
IMAS.miller_R_a_κ_δ_ζ — Function
miller_R_a_κ_δ_ζ(pr, pz, r_at_max_z, max_z, r_at_min_z, min_z, z_at_max_r, max_r, z_at_min_r, min_r)Returns named tuple with R0, a, κ, δu, δl, ζou, ζol, ζil, ζiu of a contour
miller_R_a_κ_δ_ζ(pr::Vector{T}, pz::Vector{T}) where {T<:Real}IMAS.fluxsurface_extrema — Function
fluxsurface_extrema(pr::Vector{T}, pz::Vector{T}) where {T<:Real}Returns extrema indexes and values of R,Z flux surfaces vectors:
imaxr, iminr,
imaxz, iminz,
r_at_max_z, max_z,
r_at_min_z, min_z,
z_at_max_r, max_r,
z_at_min_r, min_rIMAS.luce_squareness — Function
luce_squareness(pr::AbstractVector{T}, pz::AbstractVector{T}, r_at_max_z::T, max_z::T, r_at_min_z::T, min_z::T, z_at_max_r::T, max_r::T, z_at_min_r::T, min_r::T) where {T<:Real}Squareness from: "An analytic functional form for characterization and generation of axisymmetric plasma boundaries" T.C. Luce, Plasma Phys. Control. Fusion 55 (2013) http://dx.doi.org/10.1088/0741-3335/55/9/095009
Returns: zetaou, zetaol, zetail, zetaiu
IMAS.areal_elongation — Function
areal_elongation(eqt::IMAS.equilibrium__time_slice)A measure of the plasma elongation based on the averaged cross-sectional area of the plasma, most notably used in the H98y2 scaling
See: https://iopscience.iop.org/article/10.1088/0029-5515/48/9/099801/pdf
Physics interferometer
IMAS.line_average — Function
line_average(
rho_interp,
lcfs_r::AbstractVector{<:Real},
lcfs_z::AbstractVector{<:Real},
q::AbstractVector{<:Real},
rho_tor_norm::AbstractVector{<:Real},
r1::T1, z1::T1, ϕ1::T1,
r2::T1, z2::T1, ϕ2::T1;
n_points::Int=100
) where {T1<:Real}Low-level function that takes pre-computed ρ interpolant and LCFS boundary arrays.
This is more efficient when computing line averages for multiple channels or time-slices with the same equilibrium.
Returns named tuple with (:lineintegral, :lineaverage, :path_length)
line_average(
eqt::IMAS.equilibrium__time_slice,
q::AbstractVector{<:Real},
rho_tor_norm::AbstractVector{<:Real},
line_of_sight;
n_points::Int=100
)High-level convenience function for computing line averages
Returns named tuple with (:lineintegral, :lineaverage, :path_length)
IMAS.ne_line — Function
ne_line(::Nothing, cp1d::IMAS.core_profiles__profiles_1d)Calculates the averaged density along rhotornorm (to be used when equilibrium information is not available)
ne_line(eqt::IMAS.equilibrium__time_slice, ne_profile::AbstractVector{<:Real}, rho_ne::AbstractVector{<:Real})ne_line(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Calculates the line averaged density from the equilibrium midplane horizantal line
ne_line(dd::IMAS.dd; time0::Float64=dd.global_time)ne_line(ps::IMAS.pulse_schedule; time0=global_time(ps))returns neline from pulseschedule looking first in `pulseschedule.densitycontrol.neline.referenceand thenpulseschedule.densitycontrol.greenwald_fraction.reference`
ne_line(
equilibrium::IMAS.equilibrium,
core_profiles::IMAS.core_profiles,
interferometer::IMAS.interferometer;
times::Vector{Float64}=equilibrium.time
)Calculate time-dependent line average electron density for all interferometer channels.
Returns new interferometer IDS with synthetic data.
IMAS.interferometer! — Function
interferometer!(
intf::IMAS.interferometer,
eqt::IMAS.equilibrium__time_slice,
cp1d::IMAS.core_profiles__profiles_1d;
time0::Float64=eqt.time,
n_points::Int=100
)Populate interferometer IDS with synthetic line-averaged and line-integrated electron density measurements for a single time slice.
interferometer!(
intf::IMAS.interferometer,
eq::IMAS.equilibrium,
cp::IMAS.core_profiles;
n_points::Int=100
)Populate interferometer IDS with synthetic line-averaged and line-integrated electron density measurements for all time slices.
Physics magnetics
Physics neoclassical
IMAS.spitzer_conductivity — Function
spitzer_conductivity(ne, Te, Zeff)Calculates the Spitzer conductivity in [1/(Ω*m)]
IMAS.collision_frequencies — Function
collision_frequencies(cp1d::IMAS.core_profiles__profiles_1d{T}) where {T<:Real}Returns named tuple with nue, nui, nu_exch in [1/s]
NOTE: transpiled from the TGYRO collision_rates subroutine
IMAS.Sauter_neo2021_bootstrap — Function
Sauter_neo2021_bootstrap(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d; neo_2021::Bool=true, same_ne_ni::Bool=false)Calculates bootstrap current
neo_2021: A.Redl, et al., Phys. Plasma 28, 022502 (2021) instead of O Sauter, et al., Phys. Plasmas 9, 5140 (2002); doi:10.1063/1.1517052 (https://crppwww.epfl.ch/~sauter/neoclassical)
sameneni: assume same inverse scale length for electrons and ions
IMAS.collisionless_bootstrap_coefficient — Function
collisionless_bootstrap_coefficient(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Returns the collisional bootstrap coefficient Cbs defines as jbootfract = Cbs * sqrt(ϵ) * βp
See: Gi et al., Fus. Eng. Design 89 2709 (2014)
See: Wilson et al., Nucl. Fusion 32 257 (1992)
IMAS.nuestar — Function
nuestar(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Calculate the electron collisionality, ν_*e, as a dimensionless measure of the frequency of electron collisions relative to their characteristic transit frequency.
IMAS.nuistar — Function
nuistar(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Calculate the ion collisionality, ν_*i, as a dimensionless measure of the frequency of ion collisions relative to their characteristic transit frequency.
IMAS.neo_conductivity — Function
neo_conductivity(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Calculates the neo-classical conductivity in 1/(Ohm*meter) based on the NEO 2021 modifcation
IMAS.sawtooth_conductivity — Function
sawtooth_conductivity(conductivity::AbstractVector{T}, q::AbstractVector{T}; q_sawtooth::Float64=1.0) where {T<:Real}Jardin's model for stationary sawteeth to raise q>1
IMAS.ωtor2sonic — Function
ωtor2sonic(cp1d::IMAS.core_profiles__profiles_1d; ind::Int=0)Populates ExB rotation based on ion omega_tor assuming ωpol=0.0
Returns the updated cp1d with rotation_frequency_tor_sonic populated and all ion rotation frequencies updated.
ind is the ion species index to use as reference (0 = auto-select first non-zero rotation)
IMAS.ωtor2sonic! — Function
ωtor2sonic!(cp1d::IMAS.core_profiles__profiles_1d; ind::Int=0)Updates cp1d.rotation_frequency_tor_sonic in place based on the specified ion species toroidal rotation, assuming poloidal rotation ωpol=0.0.
IMAS.sonic2ωtor — Function
sonic2ωtor(cp1d::IMAS.core_profiles__profiles_1d, ion::IMAS.core_profiles__profiles_1d___ion)Returns ion ωtor based on ExB rotation assuming ωpol=0.0
Uses the existing rotation_frequency_tor_sonic field to compute individual ion toroidal rotation frequencies
IMAS.sonic2ωtor! — Function
sonic2ωtor!(cp1d::IMAS.core_profiles__profiles_1d, ion::IMAS.core_profiles__profiles_1d___ion)In-place conversion of ExB sonic rotation to ion toroidal rotation frequencies.
Updates individual ion rotation_frequency_tor fields based on the existing ExB sonic rotation, assuming poloidal rotation ωpol=0.0.
sonic2ωtor!(cp1d::IMAS.core_profiles__profiles_1d)All ions: Updates rotation frequencies for all ion species in the profile
Physics nuclear
IMAS.reactivity — Function
reactivity(Ti::AbstractVector{<:Real}, model::String; polarized_fuel_fraction::Real=0.0)Fusion reactivity coming from H.-S. Bosch and G.M. Hale, Nucl. Fusion 32 (1992) 611. Model can be ["D+T→He4", "D+He3→He4", "D+D→T", "D+D→He3"]")
IMAS.D_T_to_He4_reactions — Function
D_T_to_He4_reactions(cp1d::IMAS.core_profiles__profiles_1d)Calculates the number of D-T thermal fusion reactions to He4 in [reactions/m³/s]
IMAS.D_D_to_He3_reactions — Function
D_D_to_He3_reactions(dd::IMAS.dd)Calculates the number of D-D thermal fusion reactions to He3 in [reactions/m³/s]
IMAS.D_D_to_T_reactions — Function
D_D_to_T_reactions(dd::IMAS.dd)Calculates the number of D-D thermal fusion reactions to T in [reactions/m³/s]
IMAS.fusion_reaction_source — Function
fusion_reaction_source(
s1d::IMAS.core_sources__source___profiles_1d,
reactivity::Vector{<:Real},
in1::Symbol,
in2::Symbol,
out::Symbol,
eV::Float64
)Add a fusion reaction source for two isotopes coming in and one coming out
IMAS.D_T_to_He4_source! — Function
D_T_to_He4_source!(cs::IMAS.core_sources, cp::IMAS.core_profiles; combine_DT::Bool)Calculates DT fusion heating with an estimation of the alpha slowing down to the ions and electrons, modifies dd.core_sources
IMAS.D_D_to_He3_source! — Function
D_D_to_He3_source!(cs::IMAS.core_sources, cp::IMAS.core_profiles)Calculates the He-3 heating source from D-D fusion reactions, estimates energy transfer to ions and electrons, modifies dd.core_sources
IMAS.D_D_to_T_source! — Function
D_D_to_T_source!(cs::IMAS.core_sources, cp::IMAS.core_profiles)Calculates the T heating source from D-D fusion reactions, estimates energy transfer to ions and electrons, modifies dd.core_sources
IMAS.D_T_to_He4_heating — Function
D_T_to_He4_power(cp1d::IMAS.core_profiles__profiles_1d; polarized_fuel_fraction::Real=0.0)Volumetric heating source of He4 particles coming from D-T reactions [W m⁻³]
IMAS.D_D_to_He3_heating — Function
D_D_to_He3_power(cp1d::IMAS.core_profiles__profiles_1d; polarized_fuel_fraction::Real=0.0)Volumetric heating source of He3 particles coming from D-D reactions [W m⁻³]
IMAS.D_D_to_T_heating — Function
D_D_to_T_power(cp1d::IMAS.core_profiles__profiles_1d; polarized_fuel_fraction::Real=0.0)Volumetric heating source of T and H particles coming from D-D reactions [W m⁻³]
IMAS.D_T_to_He4_plasma_power — Function
D_T_to_He4_plasma_power(cp1d::IMAS.core_profiles__profiles_1d; polarized_fuel_fraction::Real=0.0)Total power in He4 from D-T reaction [W]
IMAS.D_D_to_He3_plasma_power — Function
D_D_to_He3_plasma_power(cp1d::IMAS.core_profiles__profiles_1d)Total power in He3 from D-D reaction [W]
IMAS.D_D_to_T_plasma_power — Function
D_D_to_T_plasma_power(cp1d::IMAS.core_profiles__profiles_1d)Total power in T from D-D reaction [W]
IMAS.fusion_plasma_power — Function
fusion_plasma_power(cp1d::IMAS.core_profiles__profiles_1d)Calculates the total fusion power in the plasma in [W]
If D+T plasma, then D+D is neglected
fusion_plasma_power(dd::IMAS.dd)IMAS.fusion_power — Function
fusion_power(cp1d::IMAS.core_profiles__profiles_1d{T}; DD_fusion::Bool=false) where {T<:Real}Calculates the total fusion power in [W]
If D+T plasma, then D+D is neglected
fusion_power(dd::IMAS.dd; DD_fusion::Bool=false)IMAS.fusion_energy — Function
fusion_energy(dd::IMAS.dd{T}; DD_fusion::Bool=false) where {T<:Real}Fusion energy in [J]: fusion power integrated over simulation time
IMAS.fusion_neutron_power — Function
fusion_neutron_power(cp1d::IMAS.core_profiles__profiles_1d)Calculates the total fusion power in the neutrons [W]
If D+T plasma, then D+D is neglected
fusion_neutron_power(dd::IMAS.dd)Physics particles
IMAS.Particle — Type
x::T
y::T
z::T
δvx::T
δvy::T
δvz::TCartesian coordinate system centered in (R=0, Z=0); R = sqrt(X^2 + Y^2) Z = Z
IMAS.Rcoord — Function
Rcoord(p::Particle)Return R coordinate of a Particle
IMAS.define_particles — Function
define_particles(eqt::IMAS.equilibrium__time_slice, psi_norm::Vector{T}, source_1d::Vector{T}, N::Int; random_seed::Int=0) where {T<:Real}Creates a vector of particles from a 1D source (psinorm, source1d) launching N particles.
Returns also a scalar (Ipertrace) which is the intensity per trace.
IMAS.find_flux — Function
find_flux(particles::Vector{Particle{T}}, I_per_trace::T, rwall::Vector{T}, zwall::Vector{T}, dr::T, dz::T; ns::Int=10, debug::Bool=false) where {T<:Real}Returns the flux at the wall of the quantity brought by particles, together with a vector of the surface elements on the wall (wall_s)
I_per_trace is the intensity per trace
(rwall, zwall) are the coordinates of the wall
dr, dz are the grid sizes used for the 2d source generation
ns is the window size
IMAS.toroidal_intersection — Function
toroidal_intersection(r1::Real, z1::Real, r2::Real, z2::Real, px::Real, py::Real, pz::Real, vx::Real, vy::Real, vz::Real, v2::Real, vz2::Real)Compute the time of intersection between a moving particle and a toroidal surface defined by two points (r1, z1) and (r2, z2).
Returns the smallest positive time at which the particle intersects the surface segment. Returns NaN if no valid intersection occurs.
r1,z1: Coordinates of the first endpoint of the toroidal surface segment.r2,z2: Coordinates of the second endpoint of the toroidal surface segment.px,py,pz: Current position of the particle in Cartesian coordinates.vx,vy,vz: Velocity components of the particle.v2: Squared radial velocity component,vx^2 + vy^2.vz2: Squared z-velocity component,vz^2.
toroidal_intersections(wallr::Vector{T}, wallz::vector{T}, px::Real, py::Real, pz::Real, vx::Real, vy::Real, vz::Real) where {T<:Real}Returns the first time of intersection between a moving particle and a wall in toroidal geometry
wallr: Vector with r coordinates of the wall (must be closed)wallz: Vector with z coordinates of the wall (must be closed)px,py,pz: Current position of the particle in Cartesian coordinates.vx,vy,vz: Velocity components of the particle
toroidal_intersections(wallr::Vector{T}, wallz::vector{T}, px::Real, py::Real, pz::Real, vx::Real, vy::Real, vz::Real) where {T<:Real}Returns first time of intersection between a moving particle and a wall in toroidal geometry
wallr: Vector with r coordinates of the wall (must be closed)wallz: Vector with z coordinates of the wall (must be closed)px,py,pz: Current position of the particle in Cartesian coordinates.pol_angle,tor_angle: poloidal and toroidal angles of injection
IMAS.toroidal_intersections — Function
toroidal_intersections(wallr::Vector{T}, wallz::vector{T}, px::Real, py::Real, pz::Real, vx::Real, vy::Real, vz::Real) where {T<:Real}Returns all times of intersection between a moving particle and a wall in toroidal geometry
wallr: Vector with r coordinates of the wall (must be closed)wallz: Vector with z coordinates of the wall (must be closed)px,py,pz: Current position of the particle in Cartesian coordinates.vx,vy,vz: Velocity components of the particle
toroidal_intersections(wallr::Vector{T}, wallz::vector{T}, px::Real, py::Real, pz::Real, vx::Real, vy::Real, vz::Real) where {T<:Real}Returns all times of intersection between a moving particle and a wall in toroidal geometry
wallr: Vector with r coordinates of the wall (must be closed)wallz: Vector with z coordinates of the wall (must be closed)px,py,pz: Current position of the particle in Cartesian coordinates.pol_angle,tor_angle: poloidal and toroidal angles of injection
IMAS.pol_tor_angles_2_vector — Function
pol_tor_angles_2_vector(pol_angle::T, tor_angle::T) where {T<:Real}Conversion from toroidal and poloidal angles expressed in radians to unit vector function poltorangles2vector(polangle::T, torangle::T) where {T<:Real} Uses ITER convention for EC launch angles and tor_angle
IMAS.pencil_beam — Function
pencil_beam(starting_position::Vector{T}, velocity_vector::Vector{T}, time::AbstractVector{Float64})returns named tuple with (x, y, z, r) of a beam injected at given (px, py, pz) starting_position with velocity (vx, vy, vz)
pencil_beam(starting_position::Vector{T}, pol_angle::T, tor_angle::T, time::AbstractVector{Float64}) where {T<:Real}returns named tuple with (x, y, z, r) of a beam injected at given (px, py, pz) starting_position with given poloidal and toroidal angles
Physics pedestal
IMAS.blend_core_edge_Hmode — Function
blend_core_edge_Hmode(
profile::AbstractVector{<:Real},
rho::AbstractVector{<:Real},
ped_height::Real,
ped_width::Real,
tr_bound0::Real,
tr_bound1::Real;
method::Symbol=:shift)Blends the core profiles to the pedestal for H-mode profiles, making sure the Z's at trbound0 and trbound1 match the Z's from the original profile
IMAS.blend_core_edge_EPED — Function
blend_core_edge_EPED(
profile::AbstractVector{<:Real},
rho::AbstractVector{<:Real},
ped_height::Real,
ped_width::Real,
nml_bound::Real,
ped_bound::Real;
expin::Real,
expout::Real;
method::Symbol=:shift)Blends the core and pedestal for given profile to match pedheight, pedwidth using nmlbound and pedbound as blending boundaries
IMAS.blend_core_edge_Lmode — Function
blend_core_edge_Lmode(
profile::AbstractVector{<:Real},
rho::AbstractVector{<:Real},
ped_height::Real,
ped_width::Real,
tr_bound0::Real,
tr_bound1::Real)Blends the core profiles to the pedestal for L-mode profiles, making sure the Z's at tr_bound1 matches the Z's from the original profile
NOTE: pedwidth, trbound0, tr_bound1 are not utilized
IMAS.pedestal_finder — Function
pedestal_finder(profile::Vector{T}, psi_norm::Vector{T}; do_plot::Bool=false) where {T<:Real}Finds the pedetal height and width using the EPED1 definition.
NOTE: The width is limited to be between 0.01 and 0.1. If the width is at the 0.1 boundary it is likely an indication that the profile is not a typical H-mode profile. The height is the value of the profile evaluated at (1.0 - width)
IMAS.ped_height_at_09 — Function
ped_height_at_09(rho::AbstractVector{T}, profile::AbstractVector{T}, height09::T) where {T<:Real}Scale density profile so that value at rho=0.9 is height09
IMAS.pedestal_tanh_width_half_maximum — Function
pedestal_tanh_width_half_maximum(rho::AbstractVector{T}, profile::AbstractVector{T}; rho_pedestal_full_height::Float64=0.9, tanh_width_to_09_factor::Float64=0.85)Estimates the tanh-width at half-maximum of a pedestal profile, based on a hyperbolic tangent fit
rho_pedestal_full_height: The normalized flux coordinate at which the pedestal is considered to reach full height.tanh_width_to_09_factor: Factor to convert the pedestal width at full height (ρ = 0.9) to a hyperbolic tangent profile width.
Physics pf_active
IMAS.is_ohmic_coil — Function
is_ohmic_coil(coil::IMAS.pf_active__coil)Returns true/false if coil is part of the OH
IMAS.set_coils_function — Function
set_coils_function(coils::IDSvector{<:IMAS.pf_active__coil}, R0::Float64; force::Bool=false)Setup the pf_active.coil[:].function
Physics profiles
IMAS.pressure_thermal — Function
pressure_thermal(cp1d::IMAS.core_profiles__profiles_1d)core_profiles thermal pressure
pressure_thermal(cp1de::IMAS.core_profiles__profiles_1d___electrons)electrons thermal pressure
pressure_thermal(ion::IMAS.core_profiles__profiles_1d___ion)ion thermal pressure
pressure_thermal(cp1di::IMAS.IDSvector{IMAS.core_profiles__profiles_1d___ion{T}}) where {T<:Real}thermal pressure for all ions
IMAS.beta_tor_thermal_norm — Function
beta_tor_thermal_norm(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Normalized toroidal beta from thermal pressure only, defined as 100 * betatorthermal * a[m] * B0 [T] / ip [MA]
IMAS.beta_tor_norm — Function
beta_tor_norm(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Normalized toroidal beta from total pressure, defined as 100 * beta_tor * a[m] * B0 [T] / ip [MA]
IMAS.beta_tor — Function
beta_tor(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d; norm::Bool, thermal::Bool)Toroidal beta, defined as the volume-averaged total perpendicular pressure divided by (B0^2/(2*mu0)), i.e. beta_toroidal = 2 mu0 int(p dV) / V / B0^2
beta_tor(pressure_average::Real, Bt::Real)Calculates Beta_tor from pressure and Bt
IMAS.list_ions — Function
list_ions(ids1::IDS, idss::Vararg{T}; time0::Float64) where {T<:IDS}List of ions mentioned in multiple IDSs at a given time
IMAS.ion_element! — Function
ion_element!(
ion::Union{IMAS.core_profiles__profiles_1d___ion,IMAS.core_sources__source___profiles_1d___ion},
ion_symbol::Symbol;
fast::Bool=false)Fills the ion.element structure with the a and z_n information, also updates the ion.label
IMAS.ion_properties — Function
ion_properties( ion_symbol::Symbol; fast::Bool=false)Returns named tuple with z_n, a, and label information of a given ion
IMAS.binding_energy — Function
binding_energy(Z::Int, N::Int)Return the estimate binding energy in MeV
IMAS.atomic_mass — Function
atomic_mass(Z::Int, N::Int)Returns the estimated nucleus mass including the estimated effect of binding energy
IMAS.energy_thermal — Function
energy_thermal(cp1d::IMAS.core_profiles__profiles_1d)Calculates the thermal stored energy
IMAS.energy_thermal_ped — Function
energy_thermal_ped(cp1d::IMAS.core_profiles__profiles_1d, su::IMAS.summary)Calculates the pedestal contribution to the thermal stored energy by integrating over entire domain but with pedestal pressure in the core
IMAS.tau_e_thermal — Function
tau_e_thermal(cp1d::IMAS.core_profiles__profiles_1d, sources::IMAS.core_sources; , include_radiation::Bool=true, include_time_derivative::Bool=trueEvaluate thermal energy confinement time
NOTE: This can go to infinity if there's more power coming out of the plasma than there is going in
IMAS.tau_e_h98 — Function
tau_e_h98(dd::IMAS.dd; time0::Float64=dd.global_time, include_radiation::Bool=true, include_time_derivative::Bool=true)H98y2 ITER elmy H-mode confinement time scaling
NOTE: H98y2 uses aereal elongation
NOTE: (IPBChap2.pdf, pg. 72) ...for practical reasons, the power lost by radiation inside the separatrix of the existing devices has been neglected when deriving the scalings. However, for ITER, such radiation is subtracted from the loss power when calculating the projected energy confinement time.
See Table 5 in https://iopscience.iop.org/article/10.1088/0029-5515/39/12/302/pdf and https://iopscience.iop.org/article/10.1088/0029-5515/48/9/099801/pdf for additional correction with plasma_volume
IMAS.tau_e_ds03 — Function
tau_e_ds03(dd::IMAS.dd; time0::Float64=dd.global_time, include_radiation::Bool=true, include_time_derivative::Bool=true)Petty's 2003 confinement time scaling
NOTE: Petty uses elongation at the separatrix and makes no distinction between volume and line-average density
IMAS.greenwald_density — Function
greenwald_density(eqt::IMAS.equilibrium__time_slice)Simple greenwald line-averaged density limit
greenwald_density(ip::T, minor_radius::T) where {T<:Real}greenwald_density(dd::IMAS.dd)greenwald_density(ps::IMAS.pulse_schedule; time0=global_time(ps))IMAS.greenwald_fraction — Function
greenwald_fraction(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d)Greewald fraction
greenwald_fraction(dd::IMAS.dd)IMAS.ne_vol_avg — Function
ne_vol_avg(cp1d::IMAS.core_profiles__profiles_1d)Volume averaged electron density
IMAS.beta_n — Function
beta_n(beta_tor::Real, minor_radius::Real, Bt::Real, Ip::Real)Calculates BetaN from beta_tor
IMAS.pressure_avg_from_beta_n — Function
pressure_avg_from_beta_n(beta_n::Real, minor_radius::Real, Bt::Real, Ip::Real)Calculates average pressure from BetaN
IMAS.Hmode_profiles — Function
Hmode_profiles(edge::Real, ped::Real, core::Real, ngrid::Int, expin::Real, expout::Real, width::Real; offset::Real=0.0)Generate H-mode density and temperature profiles evenly spaced in the radial coordinate
edge: separatrix valueped: pedestal valuecore: on-axis valuengrid: number of radial grid pointsexpin: inner core exponent for H-mode pedestal profileexpout: outer core exponent for H-mode pedestal profilewidth: full width of pedestal (from the separatrix to the pedestal itself)offset: offset of the pedestal center
Hmode_profiles(edge::T, ped::T, ngrid::Int, expin::T, expout::T, width::T) where {T<:Real}NOTE: The core value is allowed to float
IMAS.Lmode_profiles — Function
Lmode_profiles(edge::Real, ped::Real, core::Real, ngrid::Int, expin::Real, expout::Real, width::Real)Generate L-mode density and temperature profiles evenly spaced in the radial coordinate
edge: separatrix valueped: pedestal valuengrid: number of radial grid pointsexpin: inner core exponent for H-mode pedestal profileexpout: outer core exponent for H-mode pedestal profilewidth: width of pedestal
IMAS.A_effective — Function
A_effective(cp1d::IMAS.core_profiles__profiles_1d{T}) where {T<:Real}A_effective towards L to H scaling see: G. Birkenmeier et al 2022 Nucl. Fusion 62 086005
IMAS.scaling_L_to_H_power — Function
scaling_L_to_H_power(cp1d::IMAS.core_profiles__profiles_1d, eqt::IMAS.equilibrium__time_slice; metallic_wall:Bool=true)L to H transition power threshold for metal walls and isotope effect according to: G. Birkenmeier et al 2022 Nucl. Fusion 62 086005
See also: Reducing the L-H Transition Power Threshold in DIII-D ITER-Similar-Shape Hydrogen Plasmas, L. Schmitz IAEA FEC 2020
NOTE: This scaling does not reflect a well documented nonmonotonic density dependence.
LHthreshold doubles when ion ∇B drift is away from X-point
returns Infinity if the plasma is diverted or in negative triangularity
scaling_L_to_H_power(dd::IMAS.dd; time0::Float64=dd.global_time)IMAS.L_H_threshold — Function
L_H_threshold(cs::IMAS.core_sources, cp1d::IMAS.core_profiles__profiles_1d, eqt::IMAS.equilibrium__time_slice; time0::Float64=dd.global_time)Returns ratio of Psol to Plh
L_H_threshold(dd::IMAS.dd)IMAS.satisfies_h_mode_conditions — Function
satisfies_h_mode_conditions(dd::IMAS.dd; threshold_multiplier::Float64=1.0)Returns true if the plasma is diverted, has positive triangularity, and Psol > Plh * threshold_multiplier
IMAS.ITB_profile — Function
ITB_profile(rho::AbstractVector{T}, input_profile::AbstractVector{T}, rho0::T, width::T, height_ratio::T) where {T<:Real}Add ITB to existing profile. The ITB is centered at rho0, with a full width of width, and a height expressed as a ratio of the input_profile evaluated on axis
IMAS.species — Function
species(cp1d::IMAS.core_profiles__profiles_1d; only_electrons_ions::Symbol=:all, only_thermal_fast::Symbol=:all, return_zero_densities::Bool=false)Returns species index and names (followed by "fast" if densityfast is present), for example:
(0, :electrons)
(0, :electrons_fast)
(1, :DT)
(2, :Kr83)
(3, :He4)
(1, :DT_fast)
(3, :He4_fast)IMAS.is_quasi_neutral — Function
is_quasi_neutral(cp1d::IMAS.core_profiles__profiles_1d; rtol::Float64=0.001)Returns true if quasi neutrality is satisfied within a relative tolerance
is_quasi_neutral(dd::IMAS.dd)IMAS.enforce_quasi_neutrality! — Function
enforce_quasi_neutrality!(cp1d::IMAS.core_profiles__profiles_1d, species::Symbol)If species is :electrons then updates electrons.density_thermal to meet quasi neutrality condtion.
If species is a ion species, it evaluates the difference in number of charges needed to reach quasi-neutrality, and assigns positive difference to target ion densitythermal species and negative difference to electrons densitythermal
Also, sets density to the original expression
enforce_quasi_neutrality!(dd::IMAS.dd, species::Symbol)IMAS.lump_ions_as_bulk_and_impurity — Function
lump_ions_as_bulk_and_impurity(ions::IMAS.IDSvector{<:IMAS.core_profiles__profiles_1d___ion}, rho_tor_norm::Vector{<:Real})Changes core_profiles.ion to 2 species, one bulk species (H, D, T) and one combined impurity species
IMAS.unbundle_DT! — Function
unbundle_DT!(cp1d::IMAS.core_profiles__profiles_1d; dt_fraction::Real=0.5)Split DT ion species into separate D and T species.
Finds ion species labeled "DT" and splits them into separate "D" and "T" species.
The densities, pressures, and other quantities are distributed according to dt_fraction.
IMAS.bundle_DT! — Function
bundle_DT!(cp1d::IMAS.core_profiles__profiles_1d; d_label::String="D", t_label::String="T")Combine separate D and T ion species into a single DT species.
Finds "D" and "T" ion species and combines them into a single "DT" species.
Densities and pressures are summed, temperatures are density-weighted averaged.
IMAS.new_impurity_fraction! — Function
new_impurity_fraction!(cp1d::IMAS.core_profiles__profiles_1d, impurity_name::Symbol, impurity_ne_fraction::Real)Add thermal impurity to core_profiles given impurity fraction.
The new impurity will have the same density profile shape as the electron density profile and the same ion temperature as the main ion specie.
The main ion specie will be adjusted to have quasineutrality.
IMAS.new_impurity_radiation! — Function
new_impurity_radiation!(dd::IMAS.dd, impurity_name::Symbol, total_radiated_power::Real)Add thermal impurity to core_profiles to match a total radiated power.
The new impurity will have the same density profile shape as the electron density profile and the same ion temperature as the main ion specie.
The main ion specie will be adjusted to have quasineutrality.
new_impurity_radiation!(dd::IMAS.dd, impurity_name::Symbol, time_total_radiated_power::AbstractVector{Float64}, value_total_radiated_power::AbstractVector{<:Real})IMAS.t_i_average — Function
t_i_average(cp1d::IMAS.core_profiles__profiles_1d)::Vector{<:Real}Returns the average ion temperature weighted by each species density over the total number of ion particles
Missing docstring for IMAS.edge_profile. Check Documenter's build log for details.
IMAS.core_edge_energy — Function
core_edge_energy(cp1d::IMAS.core_profiles__profiles_1d, rho_ped::Real; thermal::Bool=true)Evaluate stored energy in the core (< rhoped) or the pedestal (> rhoped)
core_edge_energy(eqt::IMAS.equilibrium__time_slice, rho_ped::Real)Evaluate stored energy in the core (< rhoped) or the pedestal (> rhoped)
IMAS.scale_ion_densities_to_target_zeff — Function
scale_ion_densities_to_target_zeff(cp1d::IMAS.core_profiles__profiles_1d{T}, rho_scale::Real, target_zeff::Real) where {T<:Real}Returns scale coefficients for main ions and impurity density profiles needed to achieve a target zeff at a specific radial location
scale_ion_densities_to_target_zeff(cp1d::IMAS.core_profiles__profiles_1d{T}, target_zeff::Vector{T}) where {T<:Real}Returns scale coefficients for main ions and impurity density profiles needed to achieve a target zeff profile
IMAS.scale_ion_densities_to_target_zeff! — Function
scale_ion_densities_to_target_zeff!(cp1d::IMAS.core_profiles__profiles_1d{T}, rho_scale::Real, target_zeff::Real) where {T<:Real}Scale main ions and impurity density profiles in place to achieve a target zeff at a specific radial location
scale_ion_densities_to_target_zeff!(cp1d::IMAS.core_profiles__profiles_1d{T}, target_zeff::Vector{T}) where {T<:Real}Scale main ions and impurity density profiles in place to achieve a target zeff profile
NOTE: this function also enforces quasi-neutrality
Physics radiation
IMAS.radiation_losses — Function
radiation_losses(sources::IMAS.core_sources)Evaluate total plasma radiation losses [W] due to bremsstrahlung, synchrotron, and line radiation
These are energy losses, so they have a negative sign.
IMAS.bremsstrahlung_source! — Function
bremsstrahlung_source!(dd::IMAS.dd)Calculates approximate NRL Bremsstrahlung radiation source and modifies dd.core_sources
IMAS.rad_sync — Function
rad_sync(ϵ::T, a::T, B0::T, ne::T, Te::T; wall_reflection_coefficient) where {T<:Real}Synchrotron radiation from Trubnikov, JETP Lett. 16 (1972) 25.0
Transpiled from gacode/tgyro/src/tgyro_rad.f90
See also: Study of heat and synchrotron radiation transport in fusion tokamak plasmas (C. Villar 1997)
IMAS.synchrotron_source! — Function
synchrotron_source!(dd::IMAS.dd; wall_reflection_coefficient=0.0)Calculates synchrotron radiation source and modifies dd.core_sources
IMAS.line_radiation_source! — Function
line_radiation_source!(dd::IMAS.dd)Calculates line radiation sources and modifies dd.core_sources
IMAS.adas21 — Function
adas21(Te, name)NOTE: Te in [keV] and output is in [erg / cm^3 / s]
12-term Chebyshev polynomial fits to ADAS data
ln[ Lz(x) ] = sum c_n T_n(x)where
Lz = cooling rate in erg/cm^3/s
T_n(x) = cos[n*arccos(x)] (Chebyshev polynomials)
T = electron temperature
ln(T/T_min)
x = -1 + 2 ---------------
ln(T_max/T_min)
c_n = tabulated polynomial coefficients for each ionAcknowledgements:
- F. Sciortino for providing access to ADAS data via Aurora
- T. Pütterich for up-to-data ADAS data
- T. Odstrčil for help/checking of 2025 updates
References:
- Open ADAS: https://open.adas.ac.uk
- T. Pütterich et al 2019 Nucl. Fusion 59 056013
Notes:
- Lz = Lzline + Lzcontinuum
- Aurora follows the radiation nomenclature of ADAS (as described here), separating "line" and "continuum" radiation. Line radiation basically comes from ADF11 PLT files and continuum radiation comes from ADF11 PRB files. Bremsstrahlung is included in the continuum term.
- For generation of fit coefficients, see tgyro/tools/radiation # Min and max values of Te
Physics sol
IMAS.OpenFieldLine — Type
r::Vector{Float64}
z::Vector{Float64}
Br::Vector{Float64}
Bz::Vector{Float64}
Bp::Vector{Float64}
Bt::Vector{Float64}
pitch::Vector{Float64}
s::Vector{Float64}
midplane_index::Int
strike_angles::Vector{Float64} # Angle in radiants between flux surface and the wall; poloidal angle
pitch_angles::Vector{Float64} # Angle in radiants between B and Btoroidal; atan(Bp/Bt)
grazing_angles::Vector{Float64} # Angle in radiants between B and the wall; grazing angle
total_flux_expansion::Vector{Float64} # Total flux expansion
poloidal_flux_expansion::Vector{Float64} # Poloidal flux expansion
wall_index::Vector{Int} # index in dd.wall where strike points intersectIMAS.sol — Function
sol(eqt::IMAS.equilibrium__time_slice, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}; levels::Union{Int,AbstractVector}=20, use_wall::Bool=true) where {T<:Real}Returns vectors of hfs and lfs OpenFieldLine
If levels is a vector, it has the values of psi from 0 to max psiwallmidplane. The function will modify levels of psi to introduce relevant sol surfaces
sol(eqt::IMAS.equilibrium__time_slice, wall::IMAS.wall; levels::Union{Int,AbstractVector}=20, use_wall::Bool=true)sol(dd::IMAS.dd; levels::Union{Int,AbstractVector}=20, use_wall::Bool=true)IMAS.find_levels_from_P — Function
find_levels_from_P(
eqt::IMAS.equilibrium__time_slice,
wall_r::AbstractVector{<:Real},
wall_z::AbstractVector{<:Real},
PSI_interpolant::Interpolations.AbstractInterpolation,
r::Vector{<:Real},
q::Vector{<:Real},
levels::Int
)Function for the discretization of the poloidal flux ψ on the SOL, based on an hypotesis of OMP radial transport through arbitrary q(r) returns vector with level of ψ, vector with matching rmidplane and q. Discretization with even steps of P = integralsep^wal 2πrq(r)dr (same power in each flux tube)
find_levels_from_P(
eqt::IMAS.equilibrium__time_slice,
wall::IMAS.wall,
PSI_interpolant::Interpolations.AbstractInterpolation,
r::Vector{<:Real},
q::Vector{<:Real},
levels::Int
)find_levels_from_P(dd::IMAS.dd, r::Vector{<:Real}, q::Vector{<:Real}, levels::Int)IMAS.find_levels_from_wall — Function
find_levels_from_wall(
eqt::IMAS.equilibrium__time_slice,
wall_r::AbstractVector{<:Real},
wall_z::AbstractVector{<:Real},
PSI_interpolant::Interpolations.AbstractInterpolation
)Function for that computes the value of psi at the points of the wall mesh in dd
find_levels_from_wall(eqt::IMAS.equilibrium__time_slice, wall::IMAS.wall, PSI_interpolant::Interpolations.AbstractInterpolation)find_levels_from_wall(dd::IMAS.dd)IMAS.line_wall_2_wall — Function
line_wall_2_wall(r::AbstractVector{T}, z::AbstractVector{T}, wall_r::AbstractVector{T}, wall_z::AbstractVector{T}, RA::Real, ZA::Real) where {T<:Real}Returns r, z coordinates of open field line contained within wall, as well as angles of incidence at the strike locations
RA and ZA are the coordinate of the magnetic axis
IMAS.identify_strike_surface — Function
identify_strike_surface(ofl::OpenFieldLine, divertors::IMAS.divertors)Returns vector of two tuples with three integers each, identifying the indexes of the divertor/target/tile that the field line intersections
When a field line does not intersect a divertor target, then the tuple returned is (0, 0, 0)
IMAS.divertor_totals_from_targets! — Function
divertor_totals_from_targets!(divertor::IMAS.divertors__divertor)Sets total values of properties by summing over the corresponding value on the targets
- powerblackbody
- power_conducted
- power_convected
- power_currents
- power_incident
- power_neutrals
- power_radiated
- powerrecombinationneutrals
- powerrecombinationplasma
IMAS.Bpol_omp — Function
Bpol_omp(eqt::IMAS.equilibrium__time_slice)Poloidal magnetic field magnitude evaluated at the outer midplane
IMAS.power_sol — Function
power_sol(core_sources::IMAS.core_sources, cp1d::IMAS.core_profiles__profiles_1d; time0::Float64=global_time(cp1d))Total power coming out of the SOL [W]
NOTE: This function returns 1.0 [W] if power is less than that so that SOL quantities remain finite
power_sol(dd::IMAS.dd; time0::Float64=dd.global_time)IMAS.widthSOL_loarte — Function
widthSOL_loarte(B0::T, q95::T, Psol::T) where {T<:Real}Returns midplane power decay length λ_q in meters
widthSOL_loarte(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d, core_sources::IMAS.core_sources)widthSOL_loarte(dd::IMAS.dd)IMAS.widthSOL_sieglin — Function
widthSOL_sieglin(R0::T, a::T, Bpol_omp::T, Psol::T, ne_ped::T) where {T<:Real}Returns integral power decay length λ_int in meters Eich scaling(NF 53 093031) & B. Sieglin PPCF 55 (2013) 124039
widthSOL_sieglin(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d, core_sources::IMAS.core_sources)widthSOL_sieglin(dd::IMAS.dd)IMAS.widthSOL_eich — Function
widthSOL_eich(R0::T, a::T, Bpol_omp::T, Psol::T) where {T<:Real}Returns midplane power decay length λ_q in meters
Eich scaling (NF 53 093031)
widthSOL_eich(eqt::IMAS.equilibrium__time_slice, Psol::Real)widthSOL_eich(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d, core_sources::IMAS.core_sources)widthSOL_eich(dd::IMAS.dd)IMAS.q_pol_omp_eich — Function
q_pol_omp_eich(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d, core_sources::IMAS.core_sources)Poloidal heat flux [W/m²] at the outer midplane based on Eich λ_q
q_pol_omp_eich(dd::IMAS.dd)IMAS.q_par_omp_eich — Function
q_par_omp_eich(eqt::IMAS.equilibrium__time_slice, cp1d::IMAS.core_profiles__profiles_1d, core_sources::IMAS.core_sources)Parallel heat flux [W/m²] at the outer midplane based on Eich λ_q
q_par_omp_eich(dd::IMAS.dd)Missing docstring for IMAS.q_par_omp_eich. Check Documenter's build log for details.
IMAS.find_strike_points — Function
find_strike_points(pr::AbstractVector{T1}, pz::AbstractVector{T1}, wall_r::AbstractVector{T2}, wall_z::AbstractVector{T2}) where {T1<:Real, T2<:Real}Real}Finds strike points and angles of incidence between two paths
find_strike_points(
eqt::IMAS.equilibrium__time_slice{T1},
wall_r::AbstractVector{T2},
wall_z::AbstractVector{T2},
psi_last_closed::Real,
psi_first_open::Real;
strike_surfaces_r::AbstractVector{T3}=wall_r,
strike_surfaces_z::AbstractVector{T3}=wall_z,
private_flux_regions::Bool=true
) where {T1<:Real,T2<:Real,T3<:Real}Finds equilibrium strike points, angle of incidence between wall and strike leg, and the minimum distance between the surface where a strike point is located (both private and encircling) and the last closed flux surface (encircling)
find_strike_points(
eqt::IMAS.equilibrium__time_slice{T1},
wall_r::AbstractVector{T2},
wall_z::AbstractVector{T2},
psi_last_closed::Real,
psi_first_open::Real,
dv::IMAS.divertors{T1};
private_flux_regions::Bool=true
) where {T1<:Real,T2<:Real}Return strike points location in the divertors
IMAS.find_strike_points! — Function
find_strike_points!(
eqt::IMAS.equilibrium__time_slice{T1},
wall_r::AbstractVector{T2},
wall_z::AbstractVector{T2},
psi_last_closed::Real,
psi_first_open::Real,
dv::IMAS.divertors{T1};
private_flux_regions::Bool=true,
in_place::Bool=true
) where {T1<:Real,T2<:Real}Adds strike points location to equilibrium IDS
find_strike_points!(
eqt::IMAS.equilibrium__time_slice{T1},
wall_r::AbstractVector{T2},
wall_z::AbstractVector{T2},
psi_last_closed::Real,
psi_first_open::Real
) where {T1<:Real,T2<:Real}Physics sources
IMAS.fusion_source! — Function
fusion_source!(cs::IMAS.core_sources, cp::IMAS.core_profiles; DD_fusion::Bool=false)Calculates fusion source from D-T and D-D reactions and adds them to dd.core_sources
If D+T plasma, then D+D is neglected
If D+D plasma fusion is included depending on DD_fusion switch
NOTE: if found, it removes :fusion source that were not generated by FUSE, to avoid double counting
fusion_source!(dd::IMAS.dd; DD_fusion::Bool=false)IMAS.collisional_exchange_source! — Function
collisional_exchange_source!(dd::IMAS.dd)Calculates collisional exchange source and adds it to dd.core_sources
IMAS.ohmic_source! — Function
ohmic_source!(dd::IMAS.dd)Calculates the ohmic source from data in dd.core_profiles and adds it to dd.core_sources
IMAS.bootstrap_source! — Function
bootstrap_source!(dd::IMAS.dd)Calculates the bootsrap current from profiles in dd.core_profiles and adds it both as source in dd.core_sources and cp1d.j_bootstrap
IMAS.radiation_source! — Function
radiation_source!(dd::IMAS.dd)Calculates the total radiation by calling:
- bremsstrahlung_source!()
- lineradiationsource!()
- synchrotron_source!()
NOTE: if found, it removes total :radiation source to avoid double counting radiation
IMAS.sources! — Function
sources!(dd::IMAS.dd; bootstrap::Bool=true, DD_fusion::Bool=false)Calculates intrisic sources and sinks, and adds them to dd.core_sources
IMAS.time_derivative_source! — Function
time_derivative_source!(dd::IMAS.dd, cp1d_old::IMAS.core_profiles__profiles_1d, Δt::Float64; zero_out::Bool)Calculates time dependent sources and sinks, and adds them to dd.core_sources
These are the ∂/∂t term in the transport equations
time_derivative_source!(dd::IMAS.dd; zero_out::Bool=false)IMAS.total_mass_density — Function
total_mass_density(cp1d::IMAS.core_profiles__profiles_1d)Finds the total mass density [kg/m^-3]
IMAS.total_power_inside — Function
total_power_inside(
core_sources::IMAS.core_sources,
cp1d::IMAS.core_profiles__profiles_1d;
time0::Float64=global_time(cp1d),
include_radiation::Bool=true,
include_time_derivative::Bool=true
)Returns total power inside of the separatrix
IMAS.total_power_source — Function
total_power_source(source::IMAS.core_sources__source___profiles_1d)Returns the total power (electron + ion) for a single source
IMAS.total_power_time — Function
total_power_time(core_sources::IMAS.core_sources, include_indexes::Vector{<:Integer})Returns tuple of vectors with the total thermal power and time_array for given set of sources selected by identifier.index
IMAS.retain_source — Function
retain_source(source::IMAS.core_sources__source, all_indexes::Vector{Int}, include_indexes::Vector{Int}, exclude_indexes::Vector{Int})Function that decides whether a source should be kept or ignored when totaling sources
Missing docstring for IMAS.retain_source. Check Documenter's build log for details.
IMAS.total_radiation_sources — Function
total_radiation_sources(dd::IMAS.dd; time0::Float64=dd.global_time, kw...)IMAS.sawteeth_source! — Function
sawteeth_source!(dd::IMAS.dd; qmin_desired::Float64=1.0)Model sawteeth by flattening all sources where abs(q) drops below qmin_desired
sawteeth_source!(dd::IMAS.dd{T}, rho0::T) where {T<:Real}Model sawteeth by flattening all sources within the inversion radius rho0
IMAS.new_source — Function
new_source(
source::IMAS.core_sources__source,
index::Int,
name::String,
rho::Union{AbstractVector,AbstractRange},
volume::Union{AbstractVector,AbstractRange},
area::Union{AbstractVector,AbstractRange};
electrons_energy::Union{AbstractVector,Missing}=missing,
electrons_power_inside::Union{AbstractVector,Missing}=missing,
total_ion_energy::Union{AbstractVector,Missing}=missing,
total_ion_power_inside::Union{AbstractVector,Missing}=missing,
electrons_particles::Union{AbstractVector,Missing}=missing,
electrons_particles_inside::Union{AbstractVector,Missing}=missing,
j_parallel::Union{AbstractVector,Missing}=missing,
current_parallel_inside::Union{AbstractVector,Missing}=missing,
momentum_tor::Union{AbstractVector,Missing}=missing,
torque_tor_inside::Union{AbstractVector,Missing}=missing
)Populates the IMAS.core_sources__source with given heating, particle, current, momentun profiles
IMAS.total_power — Function
total_power(
ps::Union{IMAS.pulse_schedule,IMAS.pulse_schedule__ec,IMAS.pulse_schedule__ic,IMAS.pulse_schedule__lh,IMAS.pulse_schedule__nbi},
times::AbstractVector{Float64};
time_smooth::Float64)Total injected power interpolated on a given time basis
Physics technology
IMAS.coil_technology — Function
coil_technology(technology::Symbol, coil_type::Symbol)Return coil parameters from technology and coil type [:oh, :tf, :pf_active]"
IMAS.fraction_conductor — Function
fraction_conductor(coil_tech::Union{IMAS.build__pf_active__technology,IMAS.build__oh__technology,IMAS.build__tf__technology})returns the fraction of (super)conductor in a coil technology
IMAS.GAMBL_blanket — Function
GAMBL_blanket(bm::IMAS.blanket__module)Define layers for the dd.blanket.module for a GAMBL type blanket technology
Physics tf
IMAS.tf_ripple — Function
tf_ripple(r, R_tf::Real, N_tf::Integer)Evaluate fraction of toroidal magnetic field ripple at r [m] generated from N_tf toroidal field coils with outer leg at R_tf [m]
IMAS.R_tf_ripple — Function
R_tf_ripple(r, ripple::Real, N_tf::Integer)Evaluate location of toroidal field coils outer leg R_tf[m] at whichN_tftoroidal field coils generate a given fraction of toroidal magnetic field ripple atr` [m]
IMAS.top_view_outline — Function
top_view_outline(tf::IMAS.build__tf, n::Int=1; cutouts::Bool=false)returns x,y outline (top view) of the n'th TF coil
Physics thermal loads
IMAS.particle_heat_flux — Function
particle_heat_flux(
SOL::OrderedCollections.OrderedDict{Symbol,Vector{OpenFieldLine}},
rmesh::AbstractVector{<:Real},
zmesh::AbstractVector{<:Real},
r::Vector{<:Real},
q::Vector{<:Real})Computes the heat flux on the wall due to the influx of charged particles using the Scrape Off-Layer, the mesh, and an hypothesis of the decay of the parallel heat flux at the OMP.
Assumption: Points at the extrema of OpenFieldLines in SOL are in the mesh (see mesherheatflux).
Physics transport
IMAS.profile_from_z_transport — Function
profile_from_z_transport(
profile_old::AbstractVector{<:Real},
rho::AbstractVector{<:Real},
transport_grid::AbstractVector{<:Real},
z_transport_grid::AbstractVector{<:Real},
rho_ped::Real=0.0)Updates profileold with the scale lengths given by ztransport_grid
If rho_ped > transport_grid[end] then scale-length is linearly interpolated between transport_grid[end] and rho_ped
if rho_ped < transport_grid[end] then scale-length then boundary condition is at transport_grid[end]
IMAS.profile_from_rotation_shear_transport — Function
profile_from_rotation_shear_transport(
profile_old::AbstractVector{<:Real},
rho::AbstractVector{<:Real},
transport_grid::AbstractVector{<:Real},
rotation_shear_transport_grid::AbstractVector{<:Real},
rho_ped::Real=0.0)Updates rotation profile using the rotation shear given by rotationsheartransport_grid
We integrate dω/dr to get ω(r).
If rho_ped > transport_grid[end] then rotation shear is linearly interpolated between transport_grid[end] and rho_ped
If rho_ped < transport_grid[end] then boundary condition is at transport_grid[end]
IMAS.total_fluxes — Function
total_fluxes(
core_transport::IMAS.core_transport{T},
cp1d::IMAS.core_profiles__profiles_1d,
rho_total_fluxes::AbstractVector{<:Real};
time0::Float64) where {T<:Real}Sums up all the fluxes and returns it as a core_transport.model IDS
total_fluxes(dd::IMAS.dd, rho_total_fluxes::AbstractVector{<:Real}=dd.core_profiles.profiles_1d[].grid.rho_tor_norm; time0::Float64=dd.global_time)Plot
Real
Base.fill! — Function
fill!(ids_new::IDS{<:T1}, ids::IDS{<:T2}, field::Symbol) where {T1<:Measurement{<:Real},T2<:Real}Function used to map fields in IDS{T} to IDS{Measurement{T}}
fill!(ids_new::IDS{<:T1}, ids::IDS{<:T2}, field::Symbol) where {T1<:Real,T2<:Measurement{<:Real}}Function used to map fields in IDS{Measurement{T}} to IDS{T}
fill!(@nospecialize(ids_new::IDS{<:T1}), @nospecialize(ids::IDS{<:T2}), field::Symbol) where {T1<:Float64,T2<:ForwardDiff.Dual{Float64,Float64,0}}Function used to map fields in IDS{Dual{T,T,0}} to IDS{T}
Signal
IMAS.heaviside — Function
heaviside(t::Float64)Heaviside step triggered at t=0
heaviside(t::Float64, t_start::Float64)Heaviside step triggered at t=t_start
IMAS.pulse — Function
pulse(t::Float64)Unitary pulse with width of 1, starting at t=0
pulse(t::Float64, t_start::Float64, Δt::Float64)Unitary pulse with given width Δt, starting at t=t_start
IMAS.ramp — Function
ramp(t::Float64)Unitary ramp from t=0 to t=1
ramp(t::Float64, ramp_fraction::Float64)Unitary ramp
The ramp_fraction defines the fraction of ramp with respect to 1.0 and must be between [0.0,1.0]
NOTE: This function is designed as is to be able to switch between ramp(t, ramp_fraction) and trap(t, ramp_fraction).
ramp(t::Float64, t_start::Float64, Δt::Float64)Unitary ramp of duration Δt, starting at t=t_start
IMAS.trap — Function
trap(t::Float64, ramp_fraction::Float64)Unitary trapezoid
The ramp_fraction defines the fraction of ramp with respect to flattop and must be between [0.0,0.5]
trap(t::Float64, t_start::Float64, Δt::Float64, ramp_fraction::Float64)Unitary trapezoid of duration Δt, starting at t=t_start
The ramp_fraction defines the fraction of ramp with respect to flattop and must be between [0.0,0.5]
IMAS.beta — Function
beta(t::Float64, mode::Float64)Unitary beta distribution
The mode [-1.0, 1.0] defines how skewed the distribution is
beta(t::Float64, t_start::Float64, Δt::Float64, mode::Float64)Unitary beta distribution of duration Δt, starting at t=t_start
The mode [-1.0, 1.0] defines how skewed the distribution is
IMAS.sequence — Function
sequence(t::Float64, t_y_sequence::Vector{Tuple{Float64,Float64}}; scheme::Symbol=:linear)returns interpolated data given a sequence (tuple) of time/value points
IMAS.moving_average — Function
moving_average(data::Vector{<:Real}, window_size::Int)Calculate the moving average of a data vector using a specified window size. The window size is always rounded up to the closest odd number to maintain symmetry around each data point.
IMAS.highpassfilter — Function
highpassfilter(signals, fs, cutoff, order=4)Apply butterworth high pass filter of given order to signals sampled with frequency fs
IMAS.lowpassfilter — Function
lowpassfilter(signals, fs, cutoff, order=4)Apply butterworth low pass filter of given order to signals sampled with frequency fs
get from
IMAS.get_from — Function
get_from(dd::IMAS.dd, what::Symbol, from_where::Symbol; time0::Float64=dd.global_time)IMAS stores the same physical quantities in different IDSs, and get_from() abstracts away the details of which IDS to access, depending on the requested quantity (what) and the specified source (from_where). This is generally useful when coupling different codes/modules/actors.
Supported quantities for what:
:ip- Plasma current [A]- Possible sources (
from_where)::equilibrium,:core_profiles,:pulse_schedule
- Possible sources (
:vacuum_r0_b0- Vacuum magnetic field parameters (major radiusr0[m], toroidal fieldb0[T])- Possible sources (
from_where)::equilibrium,:pulse_schedule
- Possible sources (
:vloop- Loop voltage [V]- Possible sources (
from_where)::equilibrium,:core_profiles,:pulse_schedule,:controllers__ip
- Possible sources (
:βn- Normalized beta [-]- Possible sources (
from_where)::equilibrium,:core_profiles
- Possible sources (
:ne_ped- Electron density at the pedestal [m^-3]- Possible sources (
from_where)::core_profiles,:summary,:pulse_schedule
- Possible sources (
:zeff_ped- Effective charge at the pedestal [-]- Possible sources (
from_where)::core_profiles,:summary,:pulse_schedule
- Possible sources (
time0 defines the time point at which to retrieve the value, default is dd.global_time.
Returns the requested physical quantity from the specified location in the IMAS data structure.