APIs

Batsrus

Batsrus.HeadType

BATSRUS output standalone header information.

Batsrus.convertIDLtoVTKMethod
convertIDLtoVTK(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.convertTECtoVTUFunction
convertTECtoVTU(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_pvdMethod
create_pvd(filepattern)

Generate PVD file for a time series collection of VTK data.

Example

create_pvd("*.vtu)
Batsrus.cutdataMethod
cutdata(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!Method
fillCellNeighbors!(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.find_grid_blockMethod
find_grid_block(batl, xyz_D)

Return processor local block index that contains a point. Input location should be given in Cartesian coordinates.

Batsrus.find_tree_nodeMethod
find_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.getfileheadMethod
getfilehead(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.getvarMethod
getvars(bd::BATLData, var::AbstractString) -> Array

Return variable data from string var. This is also supported via direct indexing,

Examples

bd["rho"]

See also: getvars.

Batsrus.getvarsMethod
getvars(bd::BATLData, names::Vector) -> Dict

Return variables' data as a dictionary from vector of names. See also: getvar.

Batsrus.ibitsMethod

Logical shifts as the Fortran instrinsic function.

Batsrus.interp1dMethod
interp1d(bd::BATLData, var::AbstractString, loc::AbstractVector{<:AbstractFloat})

Interpolate var at spatial point loc in bd.

Batsrus.loadMethod
load(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.meshgridMethod

Generating consistent 2D arrays for passing to plotting functions.

Batsrus.order_children!Method
order_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_treeMethod
order_tree(batl)

Return maximum AMR level in the used block and the Morton curve order. Set iNodeMorton_I indirect index arrays according to

  1. root node order
  2. Morton ordering for each root node
Batsrus.readheadMethod

Return header from info file. Currently only designed for 2D and 3D.

Batsrus.readtecdataMethod
readtecdata(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.setunitsMethod
setunits(filehead, 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.showheadFunction
showhead(file, filehead)

Displaying file header information.

Batsrus.showheadMethod
showhead(data)
showhead(io, data)

Display file information of data.

Batsrus.slice1dMethod
slice1d(bd::BATLData, var::AbstractString, point1::Vector, point2::Vecto)

Interpolate var along a line from point1 to point2 in bd.

Batsrus.slice2dMethod
slice2d(bd::BATLData, var::AbstractString, plotrange=[-Inf, Inf, -Inf, Inf],
   plotinterval=Inf; kwargs...)

Return 2D 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 (somehow faster) or NaturalNeighbours for scattered interpolation.
Batsrus.subsurfaceMethod
subsurface(x, y, data, limits)
subsurface(x, y, u, v, limits)

Extract subset of 2D surface dataset in ndgrid format. See also: subvolume.

Batsrus.subvolumeMethod
subvolume(x, y, z, data, limits)
subvolume(x, y, z, u, v, w, limits)

Extract subset of 3D dataset in ndgrid format. See also: subsurface.

UnitfulBatsrus

Batsrus.UnitfulBatsrus.@bu_strMacro
macro 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.HDFModule

Module for BATSRUS HDF5 file processing.

Batsrus.HDF.HDF5CommonType

BATSRUS 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.extract_fieldMethod
extract_field(file::BatsrusHDF5Uniform, var::String; kwargs...)

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 field.
Batsrus.HDF.prep_extractMethod
prep_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_blockMethod
prep_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.prepsliceMethod
prepslice(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 matching slc_new.start.
  • xu_new: adjusted lower corner coordinate matching slc_new.stop.
Batsrus.HDF.trimsliceMethod
trimslice(start, stop, step, stop_max)

Set slice's start to be nonnegative and start/stop to be within bound. Reverse slicing is not handled.