GGDUtils.jl

Installation

GGDUtils is registered with public repository FuseRegistry. For installation:

using Pkg
Pkg.Registry.add(RegistrySpec(url="https://github.com/ProjectTorreyPines/FuseRegistry.jl.git"))
Pkg.Registry.add("General")
Pkg.add("GGDUtils)

Interpolations

Several interpolation functions are available to create interpolaiton functions for data present in a GGD represented over a particular grid subset:

GGDUtils.interpFunction
interp(
    prop_values::Vector{T},
    kdtree::KDTree;
    use_nearest_n::Int=4,
    weighing::Function=(d) -> 1 / d,
) where {T <: Real}

Lowest level interpolation function. It takes a vector of property values and a KDTree defined over a 2D space with the same number of nodes as the property values. It returns a function that can be used to interpolate the property values at any point in the space.

source
interp(
    y::Vector{T},
    TPS_mats::Tuple{Matrix{U}, Matrix{U}, Matrix{U}, Vector{Tuple{U, U}}},
) where {T <: Real, U <: Real}

Lowest level function for Thin Plate Spline method

source
interp(y::Vector{T}, x::Vector{Tuple{U, U}}) where {T <: Real, U <: Real}

Thin plate smoothing interpolation function for a 2d space scalar function. The algorithm has been adopted from:

http://www.geometrictools.com/Documentation/ThinPlateSplines.pdf

This is an implementation of Euler-Lagrange equation for minimizing bending energy of a surface.

source
interp(
    prop_values::Vector{T},
    space::IMASdd.edge_profiles__grid_ggd___space
) where {T <: Real}

If the whole space is provided instead of a kdtree, calculate the kdtree for whole space. Again, here it is assumed that the property values are porvided for each node of the space.

source
interp(
    prop_values::Vector{Real},
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset
)

If a subset of the space is provided, calculate the kdtree for the subset. In this case it is assumed that the property values are provided for each element of the subset.

source
interp(
    prop::edge_profiles__prop_on_subset,
    grid_ggd::IMASdd.edge_profiles__grid_ggd,
    value_field::Symbol=:values
)

Example:

grid_ggd = dd.edge_profiles.grid_ggd[1]
get_electron_density = interp(dd.edge_profiles.ggd[1].electrons.density[1], grid_ggd)
get_e_field_par = interp(dd.edge_profiles.ggd[1].e_field[1], grid_ggd, :parallel)
source
interp(
    prop_arr::AbstractVector{T},
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    value_field::Symbol=:values
) where {T <: edge_profiles__prop_on_subset}

Example:

sol = get_grid_subset_with_index(dd.edge_profiles.grid_ggd[1], 23)
get_electron_density = interp(dd.edge_profiles.ggd[1].electrons.density, space, sol)
source
interp(
    prop_arr::AbstractVector{T},
    grid_ggd::IMASdd.edge_profiles__grid_ggd,
    grid_subset_index::Int,
    value_field::Symbol=:values
) where {T <: edge_profiles__prop_on_subset}

Example:

get_n_e_sep = interp(dd.edge_profiles.ggd[1].electrons.density, grid_ggd, 16)
source
interp(
    prop_arr::AbstractVector{T},
    TPS_mats::Tuple{Matrix{U}, Matrix{U}, Matrix{U}, Vector{Tuple{U, U}}},
    grid_subset_index::Int,
    value_field::Val{V}=Val(:values),
) where {T <: edge_profiles__prop_on_subset, U <: Real, V}

Same use case as above but allows one to reuse previously calculated TPS matrices.

Example:

TPS_mat = get_TPS_mats(dd.edge_profiles.grid_ggd[1], 5);

for it ∈ eachindex(dd.edge_profiles.ggd)
    get_n_e = interp(dd.edge_profiles.ggd[it].electrons.density, TPS_mat_sep, 5)
    println("This time step has n_e at (0, 0) = ", get_n_e(0, 0))
end

This will run faster as heavy matrix calculations will happen only once.

source
interp(eqt::IMASdd.equilibrium__time_slice)

For a given equilibrium time slice, return a function that can be used to interpolate from (r, z) space to rho (normalized toroidal flux coordinate)space.

Example:

rz2rho = interp(dd.equilibrium.time_slice[1])
rho = rz2rho.([(r, z) for r ∈ 3.0:0.01:9.0, z ∈ -5.0:0.01:5.0])
source
interp(
    prop::Vector{T},
    prof::IMASdd.core_profiles__profiles_1d,
) where {T <: Real}

Returns an inteprolation function for the core profile property values defined on normalized toroidal flux coordinate rho.

Example:

core_profile_n_e = dd.core_profiles.profiles_1d[1].electrons.density
get_n_e = interp(core_profile_n_e, dd.core_profiles.profiles_1d[1])
get_n_e(1) # Returns electron density at rho = 1 (separatix)
source
interp(
    prop::Vector{T},
    prof::IMASdd.core_profiles__profiles_1d,
    rz2rho::Function,
)

Returns an inteprolation function in (R, Z) domain for the core profile property values defined on normalized toroidal flux coordinate rho and with a provided function to convert (R,Z) to rho.

Example:

rz2rho = interp(dd.equilibrium.time_slice[1])
core_profile_n_e = dd.core_profiles.profiles_1d[1].electrons.density
get_n_e = interp(core_profile_n_e, dd.core_profiles.profiles_1d[1], rz2rho)
get_n_e(5.0, 3.5) # Returns electron density at (R, Z) = (5.0, 3.5)
source
interp(
    prop::Vector{T},
    prof::IMASdd.core_profiles__profiles_1d,
    eqt::IMASdd.equilibrium__time_slice,
) where {T <: Real}

Returns an inteprolation function in (R, Z) domain for the core profile property values defined on normalized toroidal flux coordinate rho and with a provided equilibrium time slice to get (R, Z) to rho conversion.

Example:

eqt = dd.equilibrium.time_slice[1]
core_profile_n_e = dd.core_profiles.profiles_1d[1].electrons.density
get_n_e = interp(core_profile_n_e, dd.core_profiles.profiles_1d[1], eqt)
get_n_e(5.0, 3.5) # Returns electron density at (R, Z) = (5.0, 3.5)
source
GGDUtils.get_TPS_matsFunction
get_TPS_mats(x::Vector{Tuple{U, U}}) where {U <: Real}

We can separate matrix operations of thin plate spline method to utilize the calculation for itnerpolating over same space. Similar to the idea of reusing KDTree like above.

source
GGDUtils.get_kdtreeFunction
get_kdtree(space::IMASdd.edge_profiles__grid_ggd___space)

Get a KDTree for all the cells in the space for search for nearest neighbours.

source
get_kdtree(
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
)

Get a KDTree for a subset of the space for search for nearest neighbours.

source

Subset Tools

GGDUtils.add_subset_element!Function
add_subset_element!(
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    sn::Int,
    dim::Int,
    index::Int,
    in_subset=(x...) -> true;
    kwargs...,
)

Adds a new element to girdsubset with properties space number (sn), dimension (dim), and index (index). The element is added only if the function insubset returns true.

source
add_subset_element!(
    subset,
    sn,
    dim,
    index::Vector{Int},
    in_subset=(x...) -> true;
    kwargs...,
)

Overloaded to work differently (faster) with list of indices to be added.

source
GGDUtils.get_subset_spaceFunction
get_subset_space(
    space::IMASdd.edge_profiles__grid_ggd___space,
    elements::AbstractVector{<:IMASdd.edge_profiles__grid_ggd___grid_subset___element}
)

Returns an array of space object indices corresponding to the correct objectsperdimension (nodes, edges or cells) for the subset elements.

source
GGDUtils.get_grid_subsetFunction
get_grid_subset(
    grid_ggd::IMASdd.edge_profiles__grid_ggd,
    grid_subset_index::Int,
)

Returns the gridsubset in a gridggd with the matching gridsubsetindex

source
get_grid_subset(
    grid_ggd::IMASdd.edge_profiles__grid_ggd,
    grid_subset_name::String,
)

Returns the gridsubset in a gridggd with the matching gridsubsetname

source
GGDUtils.get_subset_boundary_indsFunction
get_subset_boundary_inds(
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
)

Returns an array of space object indices corresponding to the boundary of the subset. That means, it returns indices of nodes that are at the end of open edge subset or it returns the indices of edges that are the the boundary of a cell subset. Returns an empty array if the subset is 1D (nodes).

source
GGDUtils.get_subset_boundaryFunction
get_subset_boundary(
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
)

Returns an array of elements of grid_subset generated from the boundary of the subset provided. The dimension of these elments is reduced by 1.

source
GGDUtils.subset_doFunction
subset_do(
    set_operator,
    itrs::Vararg{
        AbstractVector{<:IMASdd.edge_profiles__grid_ggd___grid_subset___element},
    };
    space::IMASdd.edge_profiles__grid_ggd___space=IMASdd.edge_profiles__grid_ggd___space(),
    use_nodes=false
)

Function to perform any set operation (intersect, union, setdiff etc.) on subset.element to generate a list of elements to go to subset object. If usenodes is true, the set operation will be applied on the set of nodes from subset.element, space argument is required for this. Note: that the arguments are subset.element (not the subset itself). Similarly, the return object is a list of IMASdd.edgeprofiles__gridggd___gridsubset___element.

source
GGDUtils.get_subset_centersFunction
get_subset_centers(
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
)

Returns an array of tuples corresponding to (r,z) coordinates of the center of cells or the center of edges in the subset space.

source
GGDUtils.project_prop_on_subset!Function
project_prop_on_subset!(
    prop_arr::AbstractVector{T},
    from_subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    to_subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    space::IMASdd.edge_profiles__grid_ggd___space,
    value_field::Symbol=:values;
    TPS_mats::Union{
        Nothing,
        Tuple{Matrix{U}, Matrix{U}, Matrix{U}, Vector{Tuple{U, U}}},
    }=nothing,
) where {T <: edge_profiles__prop_on_subset, U <: Real}

This function can be used to add another instance on a property vector representing the value in a new subset that can be taken as a projection from an existing larger subset.

Input Arguments:

  • prop: A property like electrons.density that is a vector of objects with fields coefficients, gridindex, gridsubsetindex, and values. The different instances in the vector correspond to different gridsubset for which the property is provided.
  • fromsubset: gridsubset object which is already represented in the property instance. grid subset with index 5 is populated in electrons.density already if the values for all cells are present.
  • tosubset: gridsubset which is either a smaller part of fromsubset (core, sol, idr, odr) but has same dimensions as fromsubset OR is smaller in dimension that goes through the fromsubset (coreboundary, separatix etc.)
  • space: (optional) space object in gridggd is required only when fromsubset is higher dimensional than to_subset.

Returns: NOTE: This function ends in ! which means it updates prop argument in place. But for the additional utility, this function also returns a tuple (tosubsetcenters, topropvalues) when fromsubset dimension is greater than tosubset dimension OR (tosubseteleobjinds, topropvalues) when fromsubset dimension is same as tosubset dimension)

Descriptions: tosubsetcenters: center of cells or center of edges of the tosubset where property values are defined and stored tosubseteleobjinds: Indices of the elements of tosubset where property values are defined and stored topropvalues: The projected values of the properties added to prop object in a new instance

source
project_prop_on_subset!(
    prop_arr::AbstractVector{T},
    from_subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    to_subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
    value_field::Symbol=:values,
) where {T <: edge_profiles__prop_on_subset}

If the dimensions of fromsubset and tosubset are the same, this function can be used to add another instance on a property vector representing the value in tosubset without any interpolation or use of space object. The function returns a tuple of indices of elements of tosubset and the values of the property in to_subset.

source
GGDUtils.deepcopy_subsetFunction
deepcopy_subset(subset::IMASdd.edge_profiles__grid_ggd___grid_subset)

Faster deepcopy function for gridsubset object. This function is used to create a deep copy of a gridsubset object bypassing several checks performed by IMASdd.

source
Base.:∈Function
Base.:∈(
    point::Tuple{Real, Real},
    subset_of_space::Tuple{
        IMASdd.edge_profiles__grid_ggd___grid_subset,
        IMASdd.edge_profiles__grid_ggd___space,
    },
)

Overloading ∈ operator to check if a point is inside a subset of space.

If the subset is 1-dimensional, all points are searched. If the subset is 2-dimensional, it is checked if the point is within the enclosed area. It is assumed that a 2-dimensional subset used in such a context will form a closed area. If the subset is 3-dimensional, its boundary is calculated on the fly. If used multiple times, it is recommended to calculate the boundary once and store it in a variable.

Example:

if (5.5, 0.0) ∈ (subset_sol, space)
    println("Point (5.5, 0.0) is inside the SOL subset.")
else
    println("Point (5.5, 0.0) is outside the SOL subset.")
end
source
GGDUtils.get_prop_with_grid_subset_indexFunction
get_prop_with_grid_subset_index(
    prop::AbstractVector{T},
    grid_subset_index::Int,
) where {T <: edge_profiles__prop_on_subset}

Find the edgeprofiles property instance in an array of properties that corresponds to the gridsubset_index provided.

source

Types

GGDUtils.get_types_withFunction
get_types_with(parent::Type, field::Symbol)

A type creation utility meant for searching types in IMAS database. This function returns a list of types that are fields at any level below the parent data type which have a particular field present in it.

Example:

get_types_with(IMASdd.edge_profiles, :grid_subset_index)

returns all edgeprofiles types that have a subfield named gridsubset_index.

source

Plot recipes

Several plot recipes have been defined for easy visualization.

RecipesBase.apply_recipeMethod
plot(space::IMASdd.edge_profiles__grid_ggd___space)

Plot the grid_ggd space object. Defaults to size of [600, 900] and linecolor of :black, linewidth of 0.2, and no legend.

source
RecipesBase.apply_recipeMethod
plot(
    space::IMASdd.edge_profiles__grid_ggd___space,
    subset::IMASdd.edge_profiles__grid_ggd___grid_subset,
)

Plot the a subset of a space. Defaults to size of [600, 900] and linecolor of :black, linewidth of 0.2, and no legend.

source
RecipesBase.apply_recipeMethod
plot(
    grid_ggd::IMASdd.edge_profiles__grid_ggd,
    prop::IMASdd.IDSvectorElement,
)

Plot 2D heatmap of edgeprofilesggd property on a gridggd space object. Defaults to size of [635, 900], xaxis of "R / m", yaxis of "Z / m", and no legend. If :seriescolor is not provided, :inferno color scheme is used. If :colorbartitle is not provided, the property name is used. This function creates a plot with layout [a{0.95w} b] where a is the heatmap and b is the colorbar.

source
RecipesBase.apply_recipeMethod
plot(
    grid_ggd_arr::AbstractVector{<:IMASdd.edge_profiles__grid_ggd},
    prop::IMASdd.IDSvectorElement,
)

Plot 2D heatmap of edgeprofilesggd property on a gridggd space object. Defaults to size of [635, 900], xaxis of "R / m", yaxis of "Z / m", and no legend. If :seriescolor is not provided, :inferno color scheme is used. If :colorbartitle is not provided, the property name is used. This function creates a plot with layout [a{0.95w} b] where a is the heatmap and b is the colorbar.

source
RecipesBase.apply_recipeMethod
plot(
    ifo::IMASdd.interferometer,
)

Plot all the channels of interferometer object. Optional keywords:

  • :plottype: :los(default), :ne, or :neaverage. :los plots the line of sight of the channel in a 2D plot, :ne plots the integrated ne along the line of sight vs time, and :neaverage plots the average n_e vs time.
  • :mirror: true(default) or false.
  • :mirror_length: 0.5(default).
  • :mirror_thickness: 0.1(default).
source
RecipesBase.apply_recipeMethod
plot(
    ifo_ch::IMASdd.interferometer__channel,
)

Plot individual channel of interferometer. Optional keywords:

  • :plottype: :los(default), :ne, or :neaverage. :los plots the line of sight of the channel in a 2D plot, :ne plots the integrated ne along the line of sight vs time, and :neaverage plots the average n_e vs time.
  • :mirror: true(default) or false.
  • :mirror_length: 0.5(default).
  • :mirror_thickness: 0.1(default).
source
RecipesBase.apply_recipeMethod
plot(
    ifo_ch_los::IMASdd.interferometer__channel___line_of_sight,
)

Plot line of sight of a channel of interferometer.

Default plot settings:

  • subplot: 1
  • size: [600, 900]
  • xaxis: "R / m"
  • yaxis: "Z / m"

Optional keywords:

  • :mirror: true(default) or false.
  • :mirror_length: 0.5(default).
  • :mirror_thickness: 0.1(default).
source
RecipesBase.apply_recipeMethod
plot(
    ifo_ch_n_e_line::IMASdd.interferometer__channel___n_e_line,
)

Plot line integrated electron density of a channel of interferometer.

Default plot settings:

  • subplot: 1
  • xaxis: "time / s"
  • yaxis: "Integerated n_e / m^-2"

Optional keywords:

  • :average: true or false(default). If true, plot the average n_e vs time.
source
RecipesBase.apply_recipeMethod
plot(
    ifo_ch_n_e_line_average::IMASdd.interferometer__channel___n_e_line_average,
)

Plot average electron density of a channel of interferometer.

Default plot settings:

  • subplot: 1
  • xaxis: "time / s"
  • yaxis: "Average n_e / m^-3"
source