Internal
Public APIs
Vlasiator.MetaVLSV — TypeVLSV meta data.
Vlasiator.VarInfo — TypeVariable information from the VLSV footer.
Vlasiator.check_plasma_characteristics — Methodcheck_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.
Vlasiator.compute_flux_function — Methodcompute_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
Vlasiator.extractsat — Methodextractsat(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.
Vlasiator.find_reconnection_points — Methodfind_reconnection_points(ψ::Array{T,2}; retol::Float64=1e-4,
method::Int=1) -> indices_x, indices_oFind 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.
Vlasiator.getKLdivergence — MethodgetKLdivergence(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.
Vlasiator.getcell — Methodgetcell(meta::MetaVLSV, location:::AbstractVector{<:AbstractFloat}) -> IntReturn cell ID containing the given spatial location in meter, excluding domain boundaries. Only accept 3D location.
Vlasiator.getcellcoordinates — Methodgetcellcoordinates(meta::MetaVLSV, cid::Int) -> SVector{3,Float64}Return a given cell's spatial coordinates.
Vlasiator.getcellinline — Methodgetcellinline(meta, point1::Vector, point2::Vector) -> cellids, distances, coordsReturns cell IDs, distances and coordinates for every cell in a line between two given points point1 and point2. TODO: preallocation?
Vlasiator.getchildren — Methodgetchildren(meta::MetaVLSV, cid::Int) -> Vector{Int}Return direct children of cid.
Vlasiator.getdensity — Methodgetdensity(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.
Vlasiator.getheatfluxvector — Methodgetheatfluxvector(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.
Vlasiator.getlevel — Methodgetlevel(meta::MetaVLSV, cid::Int) -> IntReturn the AMR level of a given cell ID cid.
Vlasiator.getmaxwellianity — Methodgetmaxwellianity(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.
Vlasiator.getnearestcellwithvdf — Functiongetnearestcellwithvdf(meta, id::Int, species::String="proton") -> IntFind the nearest spatial cell with VDF saved for species of a given cell id associated with meta.
Vlasiator.getparent — Methodgetparent(meta::MetaVLSV, cid::Int) -> IntReturn the parent cell ID of given child cid.
Vlasiator.getpressure — Methodgetpressure(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.
Vlasiator.getsiblings — Methodgetsiblings(meta::MetaVLSV, cid::Int) -> Vector{Int}Return sibling cells of a given cid, including itself.
Vlasiator.getslicecell — Methodgetslicecell(meta, sliceoffset, dir, minCoord, maxCoord) -> idlist, indexlistFind 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.
Vlasiator.getvcellcoordinates — Methodgetvcellcoordinates(meta::MetaVLSV, vcellids::Vector; species="proton")Return velocity cells' coordinates of species and vcellids.
Vlasiator.getvelocity — Methodgetvelocity(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.
Vlasiator.hasname — MethodCheck if the XML nodes ns contain a node of name.
Vlasiator.hasparameter — Methodhasparameter(meta::MetaVLSV, param::String) -> BoolCheck if the VLSV file contains a certain parameter param.
Vlasiator.hasvariable — Methodhasvariable(meta::MetaVLSV, var::String) -> BoolCheck if the VLSV file associated with meta contains a variable var.
Vlasiator.isparent — Methodisparent(meta::MetaVLSV, cid::Int) -> BoolCheck if cid is a parent cell.
Vlasiator.issame — Functionissame(file1, file2, tol=1e-4; strict=true, verbose=false) -> BoolCheck 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%.
Vlasiator.load — Methodload(file::AbstractString)) -> MetaVLSVGenerate a MetaVLSV object from file of VLSV format.
Vlasiator.readlog — Methodreadlog(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.
Vlasiator.readparameter — Methodreadparameter(meta::MetaVLSV, param::String)Return the parameter value from the VLSV file associated with meta.
Vlasiator.readvariable — Functionreadvariable(meta::MetaVLSV, var::String, sorted::Bool=true, usemmap::Bool=false) -> ArrayReturn 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.
Vlasiator.readvariable — Methodreadvariable(meta::MetaVLSV, var::String, cid::Int) -> ArrayRead variable var in cell cid associated with meta.
Vlasiator.readvariable — Methodreadvariable(meta::MetaVLSV, var::String, ids::Vector{Int}) -> ArrayRead variable var in a collection of cells ids associated with meta. if ids is empty, return the whole sorted array of var.
Vlasiator.readvariablemeta — Methodreadvariablemeta(meta, var) -> VarInfoReturn VarInfo of var in the VLSV file associated with meta.
Vlasiator.readvcells — Methodreadvcells(meta::MetaVLSV, cid::Int; species="proton") -> vcellids, vcellfRead 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.
Vlasiator.refineslice — Methodrefineslice(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.
Vlasiator.vdfvolume — Functionvdfvolume(meta, location; species="proton", unit=SI, flimit=-1.0, verbose=false)Meshscatter plot of VDFs in 3D.
Vlasiator.viz — Functionviz(meta, var; args)Visualize Vlasiator output var in meta with various options:
Keyword arguments
axisunit- unit of axis of typeAxisUnitcolorscale- scale of colormap of typeColorScalenormal- slice normal directionvmin- minimum color valuevmax- maximum color valuecomp- 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.
Vlasiator.viz! — Functionviz!(meta, var; args)Visualize Meshes.jl object in an existing scene with options forwarded to viz.
Vlasiator.vlslices — Functionvlslices(meta::MetaVLSV, var; axisunit=SI, comp=0, origin=[0.0, 0.0, 0.0])Three orthogonal slices of var from meta.
Vlasiator.write_vlsv — Methodwrite_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.
Vlasiator.write_vtk — Methodwrite_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.
PyPlot Wrapper APIs
VlasiatorPyPlot has been moved to a standalone package in a subdirectory.
Private APIs
Vlasiator.AxisUnit — TypeAxis unit type. Currently supported: SI, EARTH.
Vlasiator.ColorScale — TypeColor scales type for 2D plots. Currently supported: Log, Linear, SymLog.
Vlasiator.PlotArgs — TypePlotting arguments.
Base.isopen — MethodDetermine whether meta is not yet closed.
Base.ndims — Methodndims(meta::MetaVLSV) -> IntReturn the simulation spatial dimension of VLSV data.
Base.size — MethodFile size in bytes.
Vlasiator.curl — Methodcurl(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.
Vlasiator.divergence — Methoddivergence(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.
Vlasiator.divergence2d — Methoddivergence2d(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.
Vlasiator.downsample_fg — Methoddownsample_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.
Vlasiator.downsample_fg_cell — MethodReturn a field solver grid subarray contained inside spatial cell cid.
Vlasiator.fillmesh — Methodfillmesh(meta::MetaVLSV, vars::Vector{String};
skipghosttype=true, maxamronly=false, verbose=false) -> celldata, vtkGhostTypeFill 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!).
Vlasiator.findindex — MethodGet the original vcell index without blocks from raw vcell index i (0-based).
Vlasiator.get1stcell — MethodReturn the first cell ID on mylevel given ncells on this level.
Vlasiator.getDomainDecomposition — MethodgetDomainDecomposition(globalsize, nprocs)Obtain decomposition of this grid over the given number of processors. Reference: fsgrid.hpp
Vlasiator.getObjInfo — MethodGeneral inquiry of element tag with tagname and attr.
Vlasiator.getRotationMatrix — MethodgetRotationMatrix(e1::AbtractMatrix) -> AngleAxisObtain the rotation matrix from the Cartesian coordinate to orthgonal base vectors vbase, where each column represents a base vector.
Example
e1 = [0.0 1.0 0.0; 1.0 0.0 0.0; 0.0 0.0 1.0]
R = getRotationMatrix(e1)Vlasiator.get_axis — Methodget_axis(axisunit::AxisUnit, plotrange::NTuple{4,<:Real}, sizes::NTuple{2,<:Integer})
get_axis(pArgs::PlotArgs)Return x and y ranges for 2D.
Vlasiator.get_fg_array_cell — Methodget_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.
Vlasiator.get_fg_indices — MethodReturn the field solver grid cell indexes containing coords (low-inclusive).
Vlasiator.get_fg_indices_cell — MethodReturns a slice tuple of fsgrid indices that are contained in the spatial cell cid.
Vlasiator.get_fg_indices_subvolume — Functionget_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.
Vlasiator.getdata2d — Methodgetdata2d(meta::MetaVLSV, var::String)Return 2d scalar/vector data. Nonpublic since it won't work with DCCRG AMR.
Vlasiator.getfooter — MethodReturn the XML footer of opened VLSV file.
Vlasiator.getindexes — MethodCompute x, y and z index of cell id on a given refinement level ilevel(0-based).
Vlasiator.getindexes — MethodCompute x, y and z indexes of all cell ids on the given refinement level (0-based).
Vlasiator.getmaxrefinement — MethodGet maximum refinement level, assuming 3D grid and dichotomy method.
Vlasiator.gradient — Methodgradient(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.
Vlasiator.init_cellswithVDF! — MethodLazily initialize arrays of cell ID with VDF and number of vblock per cell.
Vlasiator.isextrema — MethodCheck if the center point in a 3x3 matrix ψ is an extrema point.
Vlasiator.issaddle — MethodCheck if the center point in a 3x3 matrix ψ is a saddle point.
Vlasiator.prep1d — Methodprep1d(meta::MetaVLSV, var::String; i1::Int=0, i2::Int=0, comp::Int=0) -> ArrayObtain 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.
Vlasiator.prep2d — Functionprep2d(meta::MetaVLSV, var::String, comp::Union{Symbol, Int}=0) -> ArrayObtain data from meta of var for 2D plotting. Use comp to select vector components.
Vlasiator.prep2dslice — Methodprep2dslice(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.
Vlasiator.prep_vdf — Methodprep_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 inSI,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 whencenteris 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:particleand:flux.flimit: minimum VDF threshold for plotting.verbose: display the selection process.
Vlasiator.read_variable_as_fg — Methodread_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.
Vlasiator.readmesh — MethodReturn spatial mesh information.
Vlasiator.readvector — MethodReturn vector of name from the VLSV file associated with stream fid.
Vlasiator.readvmesh — MethodReturn velocity mesh information.
Vlasiator.reconstruct — Methodreconstruct(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.
Vlasiator.reorder — Methodreorder(vmesh::VMeshInfo, vcellids::Vector{Int32}) -> vcellids_originReorder vblock-organized VDF indexes into x–>y–>z indexes. vcellids are raw indices of nonzero VDFs ordered by blocks.
Vlasiator.rotateTensorToVectorZ — MethodrotateTensorToVectorZ(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
Vlasiator.save_image — Functionsave_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.
Vlasiator.set_args — Methodset_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.
Vlasiator.set_lim — Methodset_lim(vmin, vmax, data, colorscale=Linear)Set colormap limits vmin, vmax for data under scale colorscale.
Vlasiator.upsample_fsgrid_subarray! — MethodSet the elements of the fsgrid array to the value of corresponding cell ID cid.