Internal

Public APIs

Vlasiator.check_plasma_characteristicsMethod
check_plasma_characteristics(n, v, T, B)

Display characteristic plasma parameters given density n, bulk velocity v, temperature T, and magnetic field strength B in SI units.

source
Vlasiator.compute_flux_functionMethod
compute_flux_function(b::AbstractArray{T,N}, Δ::Vector{T}, nG::Int=2) where {T,N}

Calculate the 2D magnetic flux function ψ from the magnetic field b and discrete steps Δ. nG is the number of ghost cells along each dimension in the vector field. ψ is defined as $\psi = \int B_x dz = - \int B_z dx$ from $\mathbf{B} = \hat{y}\times\nabla\psi$ in the X-Z plane and Y is the out-of-plane direction. This is strictly true if B is divergence-free and the guide field By is constant. However, numerically there will be errors. The current implementation calculates ψ by integrating along -z boundary first, and then going along z. Reference: Flux function

source
Vlasiator.extractsatMethod
extractsat(files::AbstractVector{String}, var::String, cid::Int)

Multi-threaded extraction of var at a fixed cell ID cid from files. This assumes that files come from the same grid structure.

source
Vlasiator.find_reconnection_pointsMethod
find_reconnection_points(ψ::Array{T,2}; retol::Float64=1e-4,
   method::Int=1) -> indices_x, indices_o

Find X-point and O-point indices in 2D magnetic field topology from flux function ψ. The current implementation does not work for the 3 layers near the boundary.

Keywords

  • retol=1e-1: determines the relative tolerance of the ratio w.r.t. |∇ψ|² to accept a gradient as 0.
  • method=1: method 1 compute the cell-centered 1st and 2nd order derivatives and check the Hessian matrix; method 2 check the flux function at each point against its 8 neighbors, which is more deterministic.
source
Vlasiator.getKLdivergenceMethod
getKLdivergence(meta, VDF; species="proton")
getKLdivergence(meta, vcellids, vcellf; species="proton")

Obtain the KL-divergence ∫ f*log(f/g)dv, where f is the VDF from Vlasiator and g is the analytical Maxwellian distribution that generates the same density as f. The value ranges from [0, +∞], with 0 meaning perfect Maxwellian. Usually the values are quite small. Alternatively, one can pass original vcellids and vcellf directly.

source
Vlasiator.getcellMethod
getcell(meta::MetaVLSV, location:::AbstractVector{<:AbstractFloat}) -> Int

Return cell ID containing the given spatial location in meter, excluding domain boundaries. Only accept 3D location.

source
Vlasiator.getcellinlineMethod
getcellinline(meta, point1::Vector, point2::Vector) -> cellids, distances, coords

Returns cell IDs, distances and coordinates for every cell in a line between two given points point1 and point2. TODO: preallocation?

source
Vlasiator.getdensityMethod
getdensity(meta, VDF; species="proton")
getdensity(meta, vcellf; species="proton")
getdensity(vmesh::VMeshInfo, vcellf)

Get density from VDF of species associated with meta, n = ∫ f(r,v) dV. Alternatively, one can directly pass vcellids as original indices of nonzero VDFs and vcellf as their corresponding values.

source
Vlasiator.getheatfluxvectorMethod
getheatfluxvector(meta, VDF; species="proton")
getheatfluxvector(meta, vcellids, vcellf; species="proton")
getheatfluxvector(vmesh::VMeshInfo, vcellids, vcellf)

Get heat flux vector (3 components) of species from VDF associated with meta, qᵢ = m/2 * ∫ (v - u)²(v - u)ᵢ * f(r,v) dV. Alternatively, one can directly pass vcellids, vcellf, as in getdensity.

source
Vlasiator.getlevelMethod
getlevel(meta::MetaVLSV, cid::Int) -> Int

Return the AMR level of a given cell ID cid.

Warning

This function does not check if the VLSV file of meta actually contains cid; it may be shadowed by refined children.

source
Vlasiator.getmaxwellianityMethod
getmaxwellianity(meta, VDF; species="proton")
getmaxwellianity(meta, vcellids, vcellf; species="proton")

Obtain the Maxwellian similarity factor -log(1/(2n) * ∫ |f - g| dv), where f is the VDF from Vlasiator and g is the analytical Maxwellian distribution that generates the same density as f. The value ranges from [0, +∞], with 0 meaning not Maxwellian-distributed at all, and +∞ a perfect Maxwellian distribution. Alternatively, one can pass original vcellids and vcellf directly.

source
Vlasiator.getnearestcellwithvdfFunction
getnearestcellwithvdf(meta, id::Int, species::String="proton") -> Int

Find the nearest spatial cell with VDF saved for species of a given cell id associated with meta.

source
Vlasiator.getpressureMethod
getpressure(meta, VDF; species="proton")
getpressure(meta, vcellids, vcellf; species="proton")
getpressure(vmesh::VMeshInfo, vcellids, vcellf)

Get pressure tensor (6 components: Pxx, Pyy, Pzz, Pyz, Pzx, Pxy) of species from VDF associated with meta, pᵢⱼ = m * ∫ (v - u)ᵢ(v - u)ⱼ * f(r,v) dV. Alternatively, one can directly pass vcellids, vcellf, as in getdensity.

source
Vlasiator.getsiblingsMethod
getsiblings(meta::MetaVLSV, cid::Int) -> Vector{Int}

Return sibling cells of a given cid, including itself.

source
Vlasiator.getslicecellMethod
getslicecell(meta, sliceoffset, dir, minCoord, maxCoord) -> idlist, indexlist

Find the cell IDs idlist which are needed to plot a 2d cut through of a 3d mesh, in a direction dir at sliceoffset, and the indexlist, which is a mapping from original order to the cut plane and can be used to select data onto the plane.

source
Vlasiator.getvcellcoordinatesMethod
getvcellcoordinates(meta::MetaVLSV, vcellids::Vector; species="proton")

Return velocity cells' coordinates of species and vcellids.

source
Vlasiator.getvelocityMethod
getvelocity(meta, VDF; species="proton")
getvelocity(meta, vcellids, vcellf; species="proton")
getvelocity(vmesh::VMeshInfo, vcellids, vcellf)

Get bulk velocity from VDF of species, u = ∫ v * f(r,v) dV / n. Alternatively, one can directly pass vcellids, vcellf, as in getdensity.

source
Vlasiator.hasparameterMethod
hasparameter(meta::MetaVLSV, param::String) -> Bool

Check if the VLSV file contains a certain parameter param.

source
Vlasiator.hasvariableMethod
hasvariable(meta::MetaVLSV, var::String) -> Bool

Check if the VLSV file associated with meta contains a variable var.

source
Vlasiator.issameFunction
issame(file1, file2, tol=1e-4; strict=true, verbose=false) -> Bool

Check if two VLSV files file1 and file2 are approximately identical, under relative tolerance tol. If strict=true, the file size difference should be less than 1%.

source
Vlasiator.loadMethod
load(file::AbstractString)) -> MetaVLSV

Generate a MetaVLSV object from file of VLSV format.

source
Vlasiator.readlogMethod
readlog(file::AbstractString)

Read the run log file, check the iteration status and return the timestamps (exluding the last) as well as the model running speed in physical seconds per simulated seconds.

source
Vlasiator.readparameterMethod
readparameter(meta::MetaVLSV, param::String)

Return the parameter value from the VLSV file associated with meta.

source
Vlasiator.readvariableFunction
readvariable(meta::MetaVLSV, var::String, sorted::Bool=true, usemmap::Bool=false) -> Array

Return variable value of var from the VLSV file associated with meta. By default for DCCRG variables are sorted by cell ID. usemmap decides whether to use memory-mapped IO, which is especially useful for large arrays.

source
Vlasiator.readvariableMethod
readvariable(meta::MetaVLSV, var::String, cid::Int) -> Array

Read variable var in cell cid associated with meta.

source
Vlasiator.readvariableMethod
readvariable(meta::MetaVLSV, var::String, ids::Vector{Int}) -> Array

Read variable var in a collection of cells ids associated with meta. if ids is empty, return the whole sorted array of var.

source
Vlasiator.readvcellsMethod
readvcells(meta::MetaVLSV, cid::Int; species="proton") -> vcellids, vcellf

Read velocity cells of species from a spatial cell of ID cid associated with meta, and return a map of velocity cell ids vcellids and corresponding value vcellf.

source
Vlasiator.refinesliceMethod
refineslice(meta, idlist::Vector{Int}, data::AbstractVector, normal::Symbol) -> Vector
refineslice(meta, idlist::Vector{Int}, data::AbstractMatrix, normal::Symbol) -> Matrix
refineslice(meta, idlist::Vector{Int}, data::AbstractArray, dir::Int)

Generate data on the finest refinement level given cellids idlist and variable data on the slice perpendicular to normal. If data is 2D, then the first dimension is treated as vector components.

source
Vlasiator.vdfvolumeFunction
vdfvolume(meta, location; species="proton", unit=SI, flimit=-1.0, verbose=false)

Meshscatter plot of VDFs in 3D.

source
Vlasiator.vizFunction
viz(meta, var; args)

Visualize Vlasiator output var in meta with various options:

Keyword arguments

  • axisunit - unit of axis of type AxisUnit
  • colorscale - scale of colormap of type ColorScale
  • normal - slice normal direction
  • vmin - minimum color value
  • vmax - maximum color value
  • comp - selection of vector components

Notes

  • This function will only work in the presence of a Makie.jl backend via package extensions in Julia v1.9 or later versions of the language.
source
Vlasiator.viz!Function
viz!(meta, var; args)

Visualize Meshes.jl object in an existing scene with options forwarded to viz.

source
Vlasiator.vlslicesFunction
vlslices(meta::MetaVLSV, var; axisunit=SI, comp=0, origin=[0.0, 0.0, 0.0])

Three orthogonal slices of var from meta.

source
Vlasiator.write_vlsvMethod
write_vlsv(filein, fileout, newvars::Vector{Tuple{Vector, String, VarInfo}};
   force=false)

Generate a new VLSV fileout based on filein, with newvars added. force=true overwrites the existing fileout.

source
Vlasiator.write_vtkMethod
write_vtk(meta::MetaVLSV; kwargs...)
write_vtk(file::AbstractString; kwargs...)

Convert VLSV file to VTK format.

Keywords

  • vars::Vector{String}=[""]: select which variables to convert.
  • ascii::Bool=false: output stored in ASCII or compressed binary format.
  • maxamronly::Bool=false: generate image files on the highest refinement level only.
  • box::Vector: selected box range in 3D.
  • outdir::String="": output directory.
  • verbose::Bool=false: display logs during conversion.
source

PyPlot Wrapper APIs

VlasiatorPyPlot has been moved to a standalone package in a subdirectory.

Private APIs

Base.ndimsMethod
ndims(meta::MetaVLSV) -> Int

Return the simulation spatial dimension of VLSV data.

source
Vlasiator.curlMethod
curl(A::AbstractArray{T,N}, dx:Union{Vector{U}, NTuple{3,U}}) where {T,N,U}

Calculate 2nd order cell-centered ∇×A where A is a 4D array of size (3, nx, ny, nz) and dx is a vector of grid intervals in each dimension.

source
Vlasiator.divergenceMethod
divergence(A::AbstractArray{T,N}, dx::Vector{U}=ones(T, 3)) where {T,N,U}

Calculate 2nd order cell-centered ∇⋅A where A is a 4D array of size (3, nx, ny, nz) and dx is a vector of grid intervals in each dimension.

source
Vlasiator.divergence2dMethod
divergence2d(A::AbstractArray{T,N}, dx::AbstractFloat) where {T,N}

Calculate 2nd order cell-centered ∇⋅A where A is a 4D array of size (3, nx, ny, nz) and dx is the uniform grid interval in each dimension. One of the dimension must be 1.

source
Vlasiator.downsample_fgMethod
downsample_fg(meta::MetaVLSV, v_fg::Array)
downsample_fg(meta::MetaVLSV, var::String)

Downsample a field solver array v_fg to the spatial grid associated with meta.

source
Vlasiator.fillmeshMethod
fillmesh(meta::MetaVLSV, vars::Vector{String};
   skipghosttype=true, maxamronly=false, verbose=false) -> celldata, vtkGhostType

Fill the DCCRG mesh with quantity of vars on all refinement levels.

Return arguments

  • celldata::Vector{Vector{Array}}: data for each variable on each AMR level.
  • vtkGhostType::Array{UInt8}: cell status (to be completed!).
source
Vlasiator.getRotationMatrixMethod
getRotationMatrix(e1::Matrix, e2::Matrix) -> SMatrix{3,3}

Obtain the rotation matrix from orthgonal base vectors e1 to e2, such that a vector $\mathbf{u}_1$ in e1 can be expressed as $\mathbf{u}_1 = M\cdot \mathbf{u}_2$, where $M$ is the rotation matrix and $\mathbf{u}_2$ is the same vector in e2.

Example

e1 = [1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0]
e2 = [0.0 1.0 0.0; 1.0 0.0 0.0; 0.0 0.0 1.0]
R = getRotationMatrix(e1, e2)
source
Vlasiator.getRotationMatrixMethod
getRotationMatrix(axis::AbstractVector, angle::Real) --> SMatrix{3,3}

Create a rotation matrix for rotating a 3D vector around a unit axis by an angle in radians. Reference: Rotation matrix from axis and angle

Example

using LinearAlgebra
v = [-0.5, 1.0, 1.0]
v̂ = normalize(v)
θ = deg2rad(-74)
R = getRotationMatrix(v̂, θ)
source
Vlasiator.get_axisMethod
get_axis(axisunit::AxisUnit, plotrange::NTuple{4,<:Real}, sizes::NTuple{2,<:Integer})
get_axis(pArgs::PlotArgs)

Return x and y ranges for 2D.

source
Vlasiator.get_fg_array_cellMethod
get_fg_array_cell(meta::MetaVLSV, v_fg::Array, cid::Int)

Return a subarray of the field solver grid array, corresponding to the fsgrid covered by the spatial cell ID cid.

source
Vlasiator.get_fg_indices_subvolumeFunction
get_fg_indices_subvolume(meta::MetaVLSV, lower, upper, tol::Float64=1e-3)

Get indices for subarray of fsgrid variables, in a cuboid defined by lower and upper vertices. This is used for mapping a set of fsgrid cells to a given DCCRG cell. Shift the corners (lower, upper) inward by a distance controlled by tol. If direct low-inclusive behaviour is required, tol shall be set to 0.

source
Vlasiator.getdata2dMethod
getdata2d(meta::MetaVLSV, var::String)

Return 2d scalar/vector data. Nonpublic since it won't work with DCCRG AMR.

source
Vlasiator.gradientMethod
gradient(A::AbstractArray{T,N}, dx::Vector{U}) where {T,N,U}

Calculate 2nd order cell-centered ∇A where A is a scalar array and dx is a vector of grid intervals in each dimension.

Warning

The current implementation has issues at the boundary if gradient is taken multiple times.

source
Vlasiator.prep1dMethod
prep1d(meta::MetaVLSV, var::String; i1::Int=0, i2::Int=0, comp::Int=0) -> Array

Obtain a 1D slice from 2D meta of var. i1 is the index for the first dimension, and i2 is the index for the second dimension; only one of them is needed. Use comp to select vector components.

source
Vlasiator.prep2dFunction
prep2d(meta::MetaVLSV, var::String, comp::Union{Symbol, Int}=0) -> Array

Obtain data from meta of var for 2D plotting. Use comp to select vector components.

source
Vlasiator.prep2dsliceMethod
prep2dslice(meta::MetaVLSV, var::String, normal, comp, pArgs::PlotArgs)

Return data of var on a uniform 2D mesh on the finest AMR level. Use normal to select the plane orientation, and comp to select the component of a vector.

source
Vlasiator.prep_vdfMethod
prep_vdf(meta::MetaVLSV, location::AbstractVector; kwargs...)

Return the cell velocities v1, v2, bin ranges r1, r2, cell VDF values fweight, and strings of labels and titles for VDF plots.

Keywords

  • unit::AxisUnit: location unit in SI, EARTH.
  • unitv::String: velocity unit in ("km/s", "m/s").
  • limits::Vector{Real}: velocity space range given in [xmin, xmax, ymin, ymax].
  • slicetype: symbol for choosing the slice type from :xy, :xz, :yz, :bperp, :bpar1, :bpar2.
  • center: symbol for setting the reference frame from :bulk, :peak.
  • vslicethick: setting the velocity space slice thickness in the normal direction. If set to 0, the whole distribution along the normal direction is projected onto a plane. Currently this is only meaningful when center is set such that a range near the bulk/peak normal velocity is selected.
  • weight::Symbol: choosing distribution weights from phase space density or particle flux between :particle and :flux.
  • flimit: minimum VDF threshold for plotting.
  • verbose: display the selection process.
source
Vlasiator.read_variable_as_fgMethod
read_variable_as_fg(meta::MetaVLSV, var::String)

Interpolate DCCRG variable var to field solver grid size. This is an alternative method to fillmesh, but not optimized for performance.

source
Vlasiator.reconstructMethod
reconstruct(vmesh::VMeshInfo, vcellids::Vector{Int32}, vcellf::Vector{<:AbstractFloat})

Reconstruct the full VDFs in 3D. vcellids are raw indices of nonzero VDFs ordered by blocks, and vcellf are the corresponding values in each cell.

source
Vlasiator.reorderMethod
reorder(vmesh::VMeshInfo, vcellids::Vector{Int32}) -> vcellids_origin

Reorder vblock-organized VDF indexes into x–>y–>z indexes. vcellids are raw indices of nonzero VDFs ordered by blocks.

source
Vlasiator.rotateTensorToVectorZMethod
rotateTensorToVectorZ(tensor::AbstractMatrix, vector::AbstractVector) -> SMatrix{3,3}

Rotate tensor with a rotation matrix that aligns the 3rd direction with vector, which is equivalent to change the basis from (i,j,k) to (i′,j′,k′) where k′ ∥ vector. Reference: Tensor rotation

source
Vlasiator.save_imageFunction
save_image(meta::MetaVLSV, file, vars, data, vtkGhostType, level,
   ascii=false, append=true, box=[-Inf, Inf, -Inf, Inf, -Inf, Inf])

Save data of name vars at AMR level into VTK image file of name file.

Arguments

  • file::String: output file name.
  • vars::Vector{String}: variable names to be saved.
  • data::Vector{Vector{Array}}: data for all the variables on each refinement level.
  • vtkGhostType::Array{UInt8}: array for visibility control.
  • level::Int: refinement level (0-based).
  • ascii::Bool=false: save output in ASCII or binary format.
  • append::Bool=true: determines whether to append data at the end of file or do in-block writing.
  • box::Vector: selected box range in 3D.
source
Vlasiator.set_argsMethod
set_args(meta::MetaVLSV, var::String, axisunit::AxisUnit;
   normal::Symbol=:none, origin=0.0)

Set plot-related arguments of var in axisunit. normal and origin are used for 2D slices of 3D data.

source
Vlasiator.set_limMethod
set_lim(vmin, vmax, data, colorscale=Linear)

Set colormap limits vmin, vmax for data under scale colorscale.

source