APIs
Batsrus
Batsrus.BATS
— TypeBatsrus data container, with Dim
being the dimension of output.
Batsrus.Batl
— TypeBATSRUS output high-level struct.
Batsrus.BatsHead
— TypeBatsrus file head information.
Batsrus.FileList
— TypeType for Batsrus file information.
Batsrus.Head
— TypeBATSRUS output standalone header information.
Batsrus.adjust_plotrange!
— MethodAdjust 2D plot ranges.
Batsrus.allocateBuffer
— MethodCreate buffer for x and w.
Batsrus.convertIDLtoVTK
— MethodconvertIDLtoVTK(filename; gridType=1, verbose=false)
Convert 3D BATSRUS *.out to VTK. If gridType==1
, it converts to the rectilinear grid; if gridType==2
, it converts to the structured grid. If filename
does not end with "out", it tries to find the ".info" and ".tree" file with the same name tag and generates 3D unstructured VTU file.
Batsrus.convertTECtoVTU
— FunctionconvertTECtoVTU(file::AbstractString, outname="out")
convertTECtoVTU(head, data, connectivity, outname="out")
Convert unstructured Tecplot data to VTK. Note that if using voxel type data in VTK, the connectivity sequence is different from Tecplot: the 3D connectivity sequence in Tecplot is the same as the hexahedron
type in VTK, but different with the voxel
type. The 2D connectivity sequence is the same as the quad
type, but different with the pixel
type. For example, in 3D the index conversion is:
# PLT to VTK voxel index_ = [1 2 4 3 5 6 8 7]
for i in 1:2
connectivity = swaprows!(connectivity, 4*i-1, 4*i)
end
Batsrus.create_pvd
— Methodcreate_pvd(filepattern)
Generate PVD file for a time series collection of VTK data.
Example
create_pvd("*.vtu)
Batsrus.cutdata
— Methodcutdata(data, var; plotrange=[-Inf,Inf,-Inf,Inf], dir="x", sequence=1)
Get 2D plane cut in orientation dir
for var
out of 3D box data
within plotrange
. The returned 2D data lies in the sequence
plane from - to + in dir
.
Batsrus.fillCellNeighbors!
— MethodfillCellNeighbors!(batl, iCell_G, DiLevelNei_III, iNodeNei_III, nBlock_P)
Fill neighbor cell indexes for the given block. The faces, edges, and vertices are ordered from left (-) to right (+) in x-y-z sequentially.
Vertices: Edges: (10,11 ignored)
7 ––- 8 . –4– .
- . - . 7 . 8 .
5 ––- 6 . . –3– . 12 . . . . . . . . . 3 ––- 4 9 . –2– . . - . - . 5 . 6 1 ––- 2 . –1– .
Only tested for 3D.
Batsrus.fill_vector_from_scalars
— Methodfill_vector_from_scalars(bd::BATS, var)
Construct vector of var
from its scalar components. Alternatively, check get_vectors
for returning vector components as saparate arrays.
Batsrus.find_grid_block
— Methodfind_grid_block(batl, xyz_D)
Return processor local block index that contains a point. Input location should be given in Cartesian coordinates.
Batsrus.find_neighbor_for_anynode
— MethodFind neighbors for any node in the tree. Only for Cartesian now.
Batsrus.find_tree_node
— Methodfind_tree_node(batl, Coord_D)
Find the node that contains a point. The point coordinates should be given in generalized coordinates normalized by the domain size.
Batsrus.findindex
— MethodFind variable index in the BATSRUS data.
Batsrus.getConnectivity
— MethodGet cell connectivity list.
Batsrus.getRotationMatrix
— MethodgetRotationMatrix(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̂, θ)
Batsrus.getSibling
— MethodReturn sibling index (1-8) for the given block node matrix.
Batsrus.get_anisotropy
— Methodget_anisotropy(bd::BATS, species=0)
Calculate the pressure anisotropy for species
, indexing from 0.
Batsrus.get_convection_E
— MethodReturn the convection electric field from PIC outputs.
Batsrus.get_hall_E
— MethodReturn the Hall electric field from PIC outputs.
Batsrus.get_magnitude
— Functionget_magnitude(bd::BATS, var)
Calculate the magnitude of vector var
. See get_vectors
for the options.
Batsrus.get_magnitude2
— Functionget_magnitude2(bd::BATS, var)
Calculate the magnitude square of vector var
. See get_vectors
for the options.
Batsrus.get_var_range
— MethodReturn value range of var
in bd
.
Batsrus.get_vectors
— Methodget_vectors(bd::BATS, var)
Return a tuple of vectors of var
. var
can be :B
, :E
, :U
, :U0
, or :U1
.
Batsrus.getascii!
— MethodRead ascii format coordinates and data values.
Batsrus.getbinary!
— MethodRead binary format coordinates and data values.
Batsrus.getfilehead
— Methodgetfilehead(fileID::IoStream, filelist::FileList) -> NameTuple
Obtain the header information from BATSRUS output file of type
linked to fileID
.
Input arguments
fileID::IOStream
: file identifier.filelist::FileList
: file information.
Batsrus.getfilesize
— Methodgetfilesize(fileID::IOStream, lenstr::Int32, ::Val{FileType})
Return the size in bytes for one snapshot.
Batsrus.getfiletype
— MethodObtain file type.
Batsrus.getvar
— Methodgetvar(bd::BATS, var::AbstractString) -> Array
Return variable data from string var
. This is also supported via direct indexing,
Examples
bd["rho"]
Batsrus.getview
— MethodReturn view of variable var
in bd
.
Batsrus.ibits
— MethodLogical shifts as the Fortran instrinsic function.
Batsrus.interp1d
— Methodinterp1d(bd::BATS, var::AbstractString, loc::AbstractVector{<:AbstractFloat})
Interpolate var
at spatial point loc
in bd
.
Batsrus.interp1d
— Methodinterp1d(bd::BATS, var::AbstractString, point1::Vector, point2::Vector)
Interpolate var
along a line from point1
to point2
in bd
.
Batsrus.interp2d
— Methodinterp2d(bd::BATS, var::AbstractString, plotrange=[-Inf, Inf, -Inf, Inf],
plotinterval=Inf; kwargs...)
Return 2D interpolated slices of data var
from bd
. If plotrange
is not set, output data resolution is the same as the original.
Keyword Arguments
innermask=false
: Whether to mask the inner boundary with NaN.rbody=1.0
: Radius of the inner mask. Used when the rbody parameter is not found in the header.useMatplotlib=true
: Whether to Matplotlib (faster) or NaturalNeighbours for scattered
interpolation. If true, a linear interpolation is performed on a constructed triangle mesh.
Batsrus.interpolate2d_generalized_coords
— MethodPerform Triangle interpolation of 2D data W
on grid X
, Y
.
Batsrus.load
— Methodload(filename; npict=1, verbose=false)
Read BATSRUS output files. Stores the npict
snapshot from an ascii or binary data file into the arrays of coordinates x
and data w
.
Batsrus.meshgrid
— FunctionReturn the axis range for 2D outputs. See interp2d
.
Batsrus.meshgrid
— MethodGenerating consistent 2D arrays for passing to plotting functions.
Batsrus.nodeToGlobalBlock
— MethodReturn global block index for the node.
Batsrus.order_children!
— Methodorder_children!(batl::Batl, iNode, iMorton::Int, iNodeMorton_I::Vector{Int32})
Recursively apply Morton ordering for nodes below a root block. Store result into iNodeMortonI and iMortonNodeA using the iMorton index.
Batsrus.order_tree
— Methodorder_tree(batl)
Return maximum AMR level in the used block and the Morton curve order. Set iNodeMorton_I indirect index arrays according to
- root node order
- Morton ordering for each root node
Batsrus.readhead
— MethodReturn BATL header from info file. Currently only designed for 2D and 3D.
Batsrus.readlogdata
— MethodRead information from log file.
Batsrus.readtecdata
— Methodreadtecdata(file; verbose=false)
Return header, data and connectivity from BATSRUS Tecplot outputs. Both 2D and 3D binary and ASCII formats are supported.
Examples
file = "3d_ascii.dat"
head, data, connectivity = readtecdata(file)
Batsrus.readtree
— MethodReturn BATL tree structure.
Batsrus.rotateTensorToVectorZ
— MethodrotateTensorToVectorZ(tensor::AbstractMatrix, v::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
Batsrus.setunits
— Methodsetunits(head, type; distance=1.0, mp=1.0, me=1.0)
Set the units for the output files. If type is given as "SI", "CGS", "NORMALIZED", "PIC", "PLANETARY", "SOLAR", set typeunit = type
, otherwise try to guess from the fileheader. Based on typeunit
set units for distance [xSI], time [tSI], density [rhoSI], pressure [pSI], magnetic field [bSI] and current density [jSI] in SI units. Distance unit [rplanet | rstar], ion and electron mass in [amu] can be set with optional distance
, mp
and me
.
Also calculate convenient constants ti0, cs0 ... for typical formulas. This function is currently not used anywhere!
Batsrus.showhead
— Functionshowhead(file, head)
Displaying file header information.
Batsrus.showhead
— Methodshowhead(data)
showhead(io, data)
Display file information of data
.
Batsrus.slice1d
— Functionslice1d(bd, var, icut::Int=1, dir::Int=2)
Return view of variable var
in bd
along 1D slice. icut
is the index along axis dir
. dir == 1
means align with the 2nd (y) axis, dir == 2
means align with the 1st (x) axis.
Batsrus.squeeze
— MethodSqueeze singleton dimensions for an array A
.
Batsrus.subsurface
— Methodsubsurface(x, y, data, limits)
subsurface(x, y, u, v, limits)
Extract subset of 2D surface dataset in ndgrid format. See also: subvolume
.
Batsrus.subvolume
— Methodsubvolume(x, y, z, data, limits)
subvolume(x, y, z, u, v, w, limits)
Extract subset of 3D dataset in ndgrid format. See also: subsurface
.
Batsrus.swaprows!
— MethodReturn matrix X with swapped rows i and j.
UnitfulBatsrus
Batsrus.UnitfulBatsrus.@bu_str
— Macromacro bu_str(unit)
String macro to easily recall Batsrus units located in the UnitfulBatsrus
package. Although all unit symbols in that package are suffixed with _bu
, the suffix should not be used when using this macro. Note that what goes inside must be parsable as a valid Julia expression.
Example
julia> 1.0bu"u"
1.0 u
julia> 1.0bu"amucc"
1.0 amu/cc
HDF
Batsrus.HDF
— ModuleModule for BATSRUS HDF5 file processing.
Batsrus.HDF.BatsrusHDF5Uniform
— TypeBATSRUS HDF5 file with uniform Cartesian mesh.
Batsrus.HDF.HDF5Common
— TypeBATSRUS hdf5 file wrapper.
The data are stored in blocks, i.e., each field component is stored in a 4D array in the order (iblock, iz, iy, ix). This is a generic wrapper and does not assume grid type, i.e., uniform, stretched nonuniform, or AMR, etc. Classes to handle data with different grids can be derived from this class.
Batsrus.HDF.coord2index
— MethodReturn lower corner index.
Batsrus.HDF.extract_var
— Methodextract_var(file::BatsrusHDF5Uniform, var::String; kwargs...)
Extract variable var
from HDF5 file
.
Keywords
xmin
: minimum extracted coordinate in x.xmax
: maximum extracted coordinate in x.stepx
: extracted stride in x.ymin
: minimum extracted coordinate in y.ymax
: maximum extracted coordinate in y.stepy
: extracted stride in y.zmin
: minimum extracted coordinate in z.zmax
: maximum extracted coordinate in z.stepz
: extracted stride in z.verbose::Bool=true
: display type and size information of output variable.
Batsrus.HDF.global_slice_to_local_slice
— Methodglobal_slice_to_local_slice(file::BatsrusHDF5Uniform, dim, gslc, ib)
Convert global slice gslc
to local slice lslc
on a given block index ib
.
Batsrus.HDF.prep_extract
— Methodprep_extract(file::BatsrusHDF5Uniform, vmin=-Inf, vmax=Inf, step=1)
Get info for data extraction in 1D.
Keywords
vmin, vmax
: requested coordinate range (corner values).step
: stride.
Returns:
gslc
: global slice.vmin_new, vmax_new
: adjusted coordinate range (corner values).ibmin:ibmax
: block range.
Batsrus.HDF.prep_extract_per_block
— Methodprep_extract_per_block(file::BatsrusHDF5Uniform, dim, gslc, ib)
Get info for data extraction on a single block.
Arguments
gslc
: global slice from prep_extract.ib::Int
: block index.
Returns
lslc
: range to be used on the current block.ix0:ix1
: index range in the global array.
Batsrus.HDF.prepslice
— Methodprepslice(file::BatsrusHDF5Uniform; dim::Int, vmin, vmax, step=1)
Return range that covers [vmin
, vmax
) along dimension dim
.
Returns
slc_new
: trimed slice. If the object's Nx==1, then 1:1 will be returned.xl_new
: adjusted lower corner coordinate matchingslc_new.start
.xu_new
: adjusted lower corner coordinate matchingslc_new.stop
.
Batsrus.HDF.trimslice
— Methodtrimslice(start, stop, step, stop_max)
Set slice's start to be nonnegative and start/stop to be within bound. Reverse slicing is not handled.