API Reference

Exported Functions

GadgetIO.block_presentFunction
block_present(filename::String, blockname::String, blocks::Vector{String}=[""])

Checks if a given block is present in the snapshot file, or in the supplied blocks Vector.

source
GadgetIO.filter_cubeMethod
filter_cube(snap_file::String, corner_lowerleft::Array{<:Real}, corner_upperright::Array{<:Real}, 
                      parttype::Integer)

Reads positions from snap_file and returns the indices of particles contained in a box defined by corner_lowerleft and corner_upperright.

source
GadgetIO.filter_cylinderMethod
filter_cylinder(filename::String, pt1::Array{<:Real}, pt2::Array{<:Real}, r::Real)

Reads the positions contained in a file and returns the indices of particles contained in a cylinder defined by the endpoints pt1 and pt2 and radius r.

source
GadgetIO.filter_idsMethod
filter_ids(snap_file::String, selected_ids::Vector{<:Integer})

Reads IDs from snap_file and returns the indices of particles matching the selected_ids

source
GadgetIO.filter_subfindFunction
filter_subfind(sub_base::String, filter_function::Function, files=nothing)

Filters all entries in a subfind file that fulfill the 'filter_funcion' requirements and returns a Vector of HaloIDs.

Examples

# load packages
using GadgetIO, GadgetUnits

# define filter function
function find_mvir_gt_1e15(filename) 
    h = read_header(filename)
    GU = GadgetPhysical(h) # unit conversion

    # read Mvir and convert to solar masses 
    M = read_subfind(filename, "MVIR") .* GU.m_msun

    return findall(M .> 1.0e15)
end

# basename of subfind output (without .*)
sub_base = /path/to/groups_000/sub_000

# get relevant halos from first 10 files
halo_ids = filter_subfind(sub_base, find_mvir_gt_1e15, 0:9)
source
GadgetIO.find_most_massive_haloFunction
find_most_massive_halo(filebase::String [, nfiles::Int=1])

Reads the selected file and its subfiles and returns position, virial radius and a HaloID object that contains the subfile which contains the most massive halo and the position in the block.

source
GadgetIO.find_read_positionsMethod
find_read_positions( snap_base::String, geometry::AbstractGadgetGeometry;
                     parttype::Integer=0,
                     verbose::Bool=true )

Finds the number of particles and their storage location in a snapshot directory that are contained in the space defined by geometry.

source
GadgetIO.find_read_positionsMethod
find_read_positions( snap_base::String, filter_function::Function)

Finds the number of particles and their storage location in a snapshot directory that pass the filter_function.

source
GadgetIO.get_index_listMethod
get_index_list(list_to_find::Vector{<:Integer}, list_to_check::Vector{<:Integer})

Finds the indices at which list_to_check contains elements from list_to_find. If both either of the lists are not sorted it uses a Set lookup, otherwise it uses a Vector forward-search.

source
GadgetIO.get_total_particlesMethod
get_total_particles(h::AbstractGadgetHeader, parttype::Integer)

Calculates to total number of particles present in the simulation. Accounts for integer overflow.

source
GadgetIO.get_total_particlesMethod
get_total_particles(h::AbstractGadgetHeader, parttype::Integer)

Calculates to total number of particles present in the simulation. Accounts for integer overflow.

source
GadgetIO.global_idxs_to_halo_idMethod
global_idxs_to_halo_id(sub_base::String, offset::Integer, n_to_read::Integer;
                        parttype::Integer=0)

Converts a given number of indices defined by offset and n_to_read to HaloIDs.

source
GadgetIO.parse_balanceMethod
parse_balance(filename)

Reads the balance.txt log file and returns a tuple of Arrays of (step_number, timing, active).

source
GadgetIO.print_blocksMethod
print_blocks(filename::String; verbose::Bool=true)

Reads the block names of blocks in a snapshot and returns them in an array. Outputs them to console if verbose=true

source
GadgetIO.read_blockMethod
read_block(filename::String, blockname::String;
            parttype::Integer=0,
            block_position::Integer=-1,
            info::Union{Nothing,InfoLine}=nothing,
            h::Union{Nothing,SnapshotHeader}=nothing,
            offset=0, n_to_read=-1)

Reads a block in a snapshot with given name. Defaults to reading gas particles. Block Names are case sensitive.

Keyword Arguments

  • parttype: Which particle type to read (0-indexed)
    • 0: Gas (default)
    • 1: Dark Matter
    • 2: Boundary
    • 3: Bulge
    • 4: Stars
    • 5: Black Holes
    • -1: All particle types
  • block_position: Position of the block in the data file [bytes]. Can be used to speed up IO.
  • info: InfoLine for the given block. Can be used to speed up IO, or if no INFO block is present.
  • h: SnapshotHeader for given file. Can be used to speed up IO.
  • offset: Adds an offset to start reading the block at a later point. Can be used for sub-IO.
  • n_to_read: How many particles to read. Can be used for sub-IO.

Examples

In case an INFO block is present:

gas_pos = read_block(filename, "POS", parttype=0)

If not you need to supply your own InfoLine

pos_info = InfoLine("POS", Float32, 1, [1, 1, 1, 1, 1, 1])
gas_pos = read_block(filename, "POS", info=pos_info, parttype=0)
source
GadgetIO.read_blocks_filteredMethod
read_blocks_filtered( snap_base::String, blocks::Array{String};
                            filter_function::Union{Function, Nothing}=nothing, 
                            read_positions::Union{Dict, Nothing}=nothing, 
                            parttype::Integer=0, verbose::Bool=true )

Reads the specified blocks from all distributed files where particles pass the filter_function, or are given by a Dict of read_positions. For read_positions please see find_read_positions.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, block::AbstractString, i_global::AbstractVector{<:Integer}; verbose::Bool=true)

Get halo properties defined by the block for an Array of global indices. Returns an array with the requested block.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, block::AbstractString, haloids::AbstractVector{HaloID}; verbose::Bool=true)

Get halo properties defined by the block for an Array of HaloIDs. Returns an array with the requested block.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, haloid::HaloID, blockname::AbstractString; verbose::Bool=true)

Get halo property from block blockname by halo id. Returns an array or scalar depending on the block type.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, i_global::Integer, blockname::AbstractString; verbose::Bool=true)

Get halo property from block blockname by global halo index i_global (zero-based index).

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, blocks::AbstractVector{<:AbstractString}, i_global::AbstractVector{<:Integer}; verbose::Bool=true)

Get halo properties defined by an Array of blocks for an Array of global indices. Please note that this only works if all blocks are of the same halo type. Returns a dictionary with all requested blocks.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, blocks::AbstractVector{<:AbstractString}, haloids::AbstractVector{HaloID}; verbose::Bool=true)

Get halo properties defined by an Array of blocks for an Array of HaloIDs. Please note that this only works if all blocks are of the same halo type. Returns a dictionary with all requested blocks.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, blocks::AbstractVector{<:AbstractString}, haloid::HaloID; verbose::Bool=true)

Get halo properties defined by an Array of blocks for a HaloID. Please note that this only works if all blocks are of the same halo type. Returns a dictionary with all requested blocks.

source
GadgetIO.read_halo_propMethod
read_halo_prop(sub_base, blocks::AbstractVector{<:AbstractString}, i_global::Integer; verbose::Bool=true)

Get halo properties defined by an Array of blocks for a global index. Please note that this only works if all blocks are of the same halo type. Returns a dictionary with all requested blocks.

source
GadgetIO.read_halo_prop_and_idMethod
read_halo_prop_and_id(sub_base, blockname::AbstractString, i_global::Integer; verbose::Bool=true)

Get halo property and HaloID from block blockname by global halo index i_global (zero-based index). nfiles should generally be set to h.num_files, obtained from read_header. When nfiles is not passed, it is read automatically from the header.

source
GadgetIO.read_ids_in_haloMethod
read_ids_in_halo( sub_base::String, halo::HaloID; 
                  halo_type::Integer=1, verbose::Bool=true)

Reads the IDs of all particles contained in a halo.

source
GadgetIO.read_infoMethod
read_info(filename; verbose::Bool=false)

Reads the info block of a snapshot and returns the information in an array of InfoLine types. If verbose=true the blocknames are also printed to console.

source
GadgetIO.read_particles_by_idMethod
read_particles_by_id(snap_base::String, ids::Array{<:Integer}, 
                     blocks::Array{String}; 
                     parttype::Integer=0, verbose::Bool=true,
                     pos0::Array{<:Real}=[-1.234, -1.234, -1.234],
                     r0::Real=0.0)

Reads particles filtered by the provided IDs. Returns all requested blocks as entries in a Dict.

source
GadgetIO.read_particles_by_idMethod
read_particles_by_id(snap_base::String, ids::Array{<:Integer}, 
                     block::String; 
                     parttype::Integer=0, verbose::Bool=true,
                     pos0::Array{<:Real}=[-1.234, -1.234, -1.234],
                     r0::Real=0.0)

Reads particles filtered by the provided IDs. Returns the requested block as an Array.

source
GadgetIO.read_particles_in_boxMethod
read_particles_in_box(filename::String, block::String,
                      corner_lowerleft::Array{<:Real}, corner_upperright::Array{<:Real};
                      parttype::Integer=0, verbose::Bool=true,
                      use_keys::Bool=true)

Like read_particles_in_box but for a single block. Returns the block as an array.

See also: read_particles_in_box

source
GadgetIO.read_particles_in_boxMethod
read_particles_in_box(filename::String, blocks::Vector{String},
                    corner_lowerleft::Array{<:Real}, corner_upperright::Array{<:Real};
                    parttype::Integer=0, verbose::Bool=true,
                    use_keys::Bool=true)

Reads all particles within a box defined by a lower left and upper right corner for a given particle type. Returns a dictionary with all requested blocks. If use_keys=true it uses Peano-Hilbert keys to constrain the read-in, otherwise it uses a brute-force read-in with a filter function. Peano-Hilbert key based read-in is significantly faster.

source
GadgetIO.read_particles_in_geometryMethod
read_particles_in_geometry( filename::String, block::String,
                            geometry::AbstractGadgetGeometry;
                            parttype::Integer=0, verbose::Bool=true,
                            use_keys::Bool=true, do_shift_across_box_border::Bool=true)

Reads all particles within a space defined by an AbstractGeometry struct for a given particle type. Returns a dictionary with the requested block.

source
GadgetIO.read_particles_in_geometryMethod
read_particles_in_geometry( filename::String, blocks::Vector{String},
                            geometry::AbstractGadgetGeometry;
                            parttype::Integer=0, verbose::Bool=true,
                            use_keys::Bool=true, do_shift_across_box_border::Bool=true)

Reads all particles within a space defined by an AbstractGeometry struct for a given particle type. Returns a dictionary with all requested blocks. If shift_across_box_border is true, the particles are moved beyond the borders of a periodic box, if false, the periodicity of the box is still accounted for, but the particles' positions are not shifted.

source
GadgetIO.read_particles_in_haloMethod
read_particles_in_halo(snap_base::String, blocks::Array{String},
                            sub_base::String, halo::HaloID; 
                            radius::Union{Real,Nothing}=nothing,
                            rad_scale::Real=1.0, halo_type::Integer=1,
                            parttype::Integer=0, verbose::Bool=true,
                            use_keys::Bool=true)

Reads all particles of type parttype that are contained in a halo defined by its HaloID. If radius is given (in simulation units), particles are read within at least this radius times rad_scale. Otherwise, R200, RMEA, or RHMS times rad_scale is used depending on halo_type (1, 1, and 2, respectively). Returns a Dict with each of the blocks as entries.

source
GadgetIO.read_particles_in_haloMethod
read_particles_in_halo( snap_base::String, block::String,
                        sub_base::String, halo::HaloID; 
                        kwargs...)

Reads all particles of type parttype that are contained in a halo defined by its HaloID. If radius is given (in simulation units), particles are read within at least this radius times rad_scale. Otherwise, R200, RMEA, or RHMS times rad_scale is used depending on halo_type (1, 1, and 2, respectively). Returns an Array with the requested block.

source
GadgetIO.read_particles_in_volumeMethod
read_particles_in_box(filename::String, blocks::Vector{String},
                      center_pos, radius;
                      parttype::Integer=0, verbose::Bool=true,
                      use_keys::Bool=true)

Reads all particles within a box encapsulating a volume defined by center position and radius for a given particle type. Returns a dictionary with all requested blocks.

See also: read_particles_in_box

source
GadgetIO.read_snapFunction
read_snap(filename::String [, blockname::String="", parttype::Integer=-1] )

Wrapper function to read snapshot in various ways: filename only: returns the entire snapshot as a dictionary. blockname: Returns only that block. If parttype specified only for that particle type.

Examples

julia> gas_pos = read_snap(filename, "POS", 0)
source
GadgetIO.read_subfindMethod
read_subfind(filename::String, blockname::String, ids::AbstractVector{<:Integer}; return_haloid::Bool=false)

Reads the block at the given subfind indices (ids, 0-indexed).

If return_haloid is true, returns a tuple of the block array and the corresponding HaloIDs.

Example

# Read the virial masses of the first four halos in subfind:
mvir = read_subfind(filebase, "MVIR", [0, 1, 2, 3])

# or:
mvir, haloids = read_subfind(filebase, "MVIR", [0, 1, 2, 3]; return_haloid=true)
source
GadgetIO.read_subfindMethod
read_subfind(filename::String, blockname::String)

Example

To read e.g. the virial radius of halos use

R_vir = read_subfind(filename, "RVIR")
source
GadgetIO.snap_to_dictMethod
snap_to_dict(filename::String)

Reads whole snapshot into memory and returns a dictionary sorted by particle type and block names. Be cautious with large snapshots!

source
GadgetIO.write_blockFunction
write_block(f::IOStream, data,
            blockname::String="";
            snap_format::Integer=2)

Write data to a block to an opened file f.

source
GadgetIO.write_headerMethod
write_header(f::IOStream, h::SnapshotHeader; snap_format::Integer=2)

Writes the header block to an opened file f.

source
GadgetIO.write_info_blockMethod
write_info_block(f::IOStream, info::Vector{InfoLine}; snap_format::Integer=2)

Writes the info block to an opened file f.

source

Exported Types

GadgetIO.GadgetCubeMethod
GadgetCube(center::Vector{T}, side_length::T) where T

Define a GadgetCube by its center and side length.

source
GadgetIO.GadgetShellType
struct GadgetShell{T} <: AbstractGadgetGeometry
    center::Vector{T}
    radius::T
    width::T
end

Defines a shell by center, radius and width. Radius is defined to the middle of the shell thickness. To be used for read_particles_in_geometry

source
GadgetIO.HaloIDType
struct HaloID
    file::Int64
    id::Int64
end

Stores the subfile that contains the halo in file and the position in the block in id.

source
GadgetIO.InfoLineType
struct InfoLine([  block_name="", data_type=Float32, n_dim=Int32(0),
                is_present=zeros(Int32, 6) ])

Contains the data of a single entry in the INFO block of a Gadget snapshot.

Fields

NameMeaning
block_name::Stringname of the data block, e.g. "POS"
data_type::DataTypedatatype of the block, e.g. Float32 for single precision, Float64 for double
n_dim::Int32number of dimensions of the block, usually 1 or 3
is_present::Vector{Int32}array of flags for which particle type this block is present,
e.g. gas only: [ 1, 0, 0, 0, 0, 0 ],
or gas + BHs: [ 1, 0, 0, 0, 0, 1 ]
source
GadgetIO.SnapshotHeaderType
mutable struct SnapshotHeader <: AbstractGadgetHeader

Contains the data of the HEAD block of a Gadget snapshot.

Fields

NameMeaning
npart::Vector{Int32}an array of particle numbers per type in this snapshot
massarr::Vector{Float64}an array of particle masses per type in this snapshot - if zero: MASS block present
time::Float64time / scale factor of the simulation
z::Float64redshift of the simulation
flag_sfr::Int321 if simulation was run with star formation, else 0
flag_feedback::Int321 if simulation was run with stellar feedback, else 0
nall::Vector{UInt32}total number of particles in the simulation
flag_cooling::Int321 if simulation was run with cooling, else 0
num_files::Int32number of snapshots over which the simulation is distributed
omega_0::Float64Omega matter
boxsize::Float64total size of the simulation box
omega_l::Float64Omega dark enery
h0::Float64little h
flag_stellarage::Int321 if simulation was run with stellar age, else 0
flag_metals::Int321 if simulation was run with metals, else 0
npartTotalHighWord::Vector{UInt32}If Npart > 1584^3 (>2^32) this contains a high bit: ntotal = npartTotalHighWord*2^32 + nall
flag_entropy_instead_u::Int321 if snapshot U field contains entropy instead of internal energy, else 0
flag_doubleprecision::Int321 if snapshot is in double precision, else 0
flag_ic_info::Int321 if initial snapshot file contains an info block, else 0
lpt_scalingfactor::Float32factor to use second order ic generation
fill::Vector{Int32}the HEAD block needs to be filled with zeros to have a size of 256 bytes
source
GadgetIO.SubfindHeaderType
struct SubfindHeader

Contains the data of the HEAD block in the subfind output

Fields

NameMeaning
nhalos::Int32number of halos in the output file
nsubhalos::Int32number of subhalos in the output file
nfof::Int32number of particles in the FoF
ngroups::Int32number of large groups in the output file
time::Float64time / scale factor of the simulation
z::Float64redshift of the simulation
tothalos::UInt32total number of halos over all output files
totsubhalos::UInt32total number of subhalos over all output files
totfof::UInt32total number of particles in the FoF
totgroups::UInt321 if simulation was run with cooling, else 0
num_files::Int32number of files over which subfind data is distributed
boxsize::Float64total size of the simulation box
omega_0::Float64Omega matter
omega_l::Float64Omega dark enery
h0::Float64little h
flag_doubleprecision::Int321 if snapshot is in double precision, else 0
flag_ic_info::Int321 if initial snapshot file contains an info block, else 0
source

Private Functions

GadgetIO.allocate_data_dictMethod
allocate_data_dict( blocks::Array{String}, N_to_read::Integer, 
                    snap_info::Array{InfoLine}, no_mass_block::Bool )

Helper function to allocate the data Dict.

source
GadgetIO.check_blocksizeMethod
check_blocksize(f::IOStream, position_before::Integer, blocksize_before::Integer)

Checks for integer overflow in the size of the block.

source
GadgetIO.check_in_cylinderMethod
check_in_cylinder(x, pt1, pt2, r)

Checks if a 3D point x is in a cylinder defined by its endpoints pt1 and pt2 and radius r.

source
GadgetIO.check_infoMethod
check_info(filename::String, blockname::String)

Helper function to read INFO block or construct InfoLine for MASS block, if no INFO block is present. Returns a single InfoLine struct.

source
GadgetIO.construct_matched_dictMethod
construct_matched_dict(data_in::Dict{String, Union{Vector{T}, Array{T,2}}}, 
                            blocks::Array{String,1}, matched::Array{<:Integer,1}) where T

Write all matching particles to a new dictionary.

source
GadgetIO.filter_positionsMethod
filter_positions(snap_file::String, corner_lowerleft::Array{<:Real}, corner_upperright::Array{<:Real}, 
                      parttype::Integer)

Reads positions from snap_file and returns the indices of particles contained in a box defined by corner_lowerleft and corner_upperright.

source
GadgetIO.filter_sphereMethod
filter_sphere(filename::String, center::Array{<:Real}, r::Real; parttype::Integer=0)

Reads the positions contained in a file and returns the indices of particles contained in a cylinder defined by the endpoints pt1 and pt2 and radius r.

source
GadgetIO.find_files_for_keysMethod
find_files_for_keys(filebase::String, nfiles::Integer, keylist::Vector{<:Integer})

Selects the files in which the particles associated with the given Peano-Hilbert keys are stored.

source
GadgetIO.get_block_positionsMethod
get_block_positions(filename::String)

Returns a dictionary with the starting positions of all blocks in a snapshot in bits.

Example

block_positions = get_block_positions(filename)
block_positions["POS"]
source
GadgetIO.get_datatype_nameMethod
function get_datatype_name(dt::DataType, dim::Integer)

Convert the datatype into a string and then to a vector of 8bit integers to match name in Gadget.

source
GadgetIO.get_geometry_box_cornersMethod
get_geometry_box_corners(cylinder::GadgetCylinder)

Returns a tuple with the lower left and upper right corner of a box which contains the cylinder.

source
GadgetIO.get_geometry_maskMethod
get_geometry_mask(cube::GadgetCube, pos::Matrix{T}) where T

Returns the indices of all particles contained in the cube.

source
GadgetIO.get_geometry_maskMethod
get_geometry_mask(cylinder::GadgetCylinder, pos::Matrix{T}) where T

Returns the indices of all particles contained in the cylinder.

source
GadgetIO.get_geometry_maskMethod
get_geometry_mask(shell::GadgetShell, pos::Matrix{T}) where T

Returns the indices of all particles contained in the shell.

source
GadgetIO.get_geometry_maskMethod
get_geometry_mask(sphere::GadgetSphere, pos::Matrix{T}) where T

Returns the indices of all particles contained in the sphere.

source
GadgetIO.get_index_boundsMethod
get_index_bounds(ids::Vector{<:Integer}, low_bounds::Vector{<:Integer}, high_bounds::Vector{<:Integer})

Returns sorted Vector of indices i, for which low_bounds[i] ≤ ids[j] ≤ high_bounds[i] for any j. All parameters ids, low_bounds, and high_bounds have to already be sorted.

source
GadgetIO.get_index_list_arrMethod
get_index_list_arr(list_to_find::Vector{<:Integer}, list_to_check::Vector{<:Integer})

Get positions in list_to_check where list_to_check matches list_to_find. Uses forward-searching in sorted array. Both arrays have to be sorted.

source
GadgetIO.get_index_list_dictMethod
get_index_list_dict(list_to_find::Vector{<:Integer}, list_to_check::Vector{<:Integer})

Get positions in list_to_check where list_to_check matches list_to_find. Uses a Dict for lookup -> slower than the array search, but works on unsorted arrays.

source
GadgetIO.get_index_list_setMethod
get_index_list_set(list_to_find::Vector{<:Integer}, list_to_check::Vector{<:Integer})

Get positions in list_to_check where list_to_check matches list_to_find. Uses a Set for lookup -> slower than the array search get_index_list_arr, but works on unsorted arrays like get_index_list_dict, just faster.

source
GadgetIO.get_int_posMethod
get_int_pos(pos::Real, domain_corner::Real, domain_fac::Real )

Computes the integer position along the PH line.

source
GadgetIO.get_keylistMethod
get_keylist(h_key::KeyHeader, x0::Array{<:Real}, x1::Array{<:Real})

Get all Peano-Hilbert keys for domain defined by the corner points x0 and x1.

source
GadgetIO.get_requested_infoMethod
get_requested_info(snap_info::Array{InfoLine}, block::AbstractString)

Checks if the InfoLine is present for the requested block, or if the MASS block is requested, but not in the INFO block.

source
GadgetIO.peano_hilbert_keyMethod
peano_hilbert_key(bits::Integer, x::Integer, y::Integer, z::Integer)

Computes a Peano-Hilbert key for an integer triplet (x,y,z) with x,y,z typically in the range between 0 and 2^bits-1. Values outside this range can occur when reading across the periodic box borders.

source
GadgetIO.read_all_parttypesMethod
read_all_parttypes( filename::String, blockname::String;
                    block_position::Integer=-1,
                    info::Union{Nothing,InfoLine}=nothing,
                    h::Union{Nothing,SnapshotHeader}=nothing)

Reads the requested block for all particle types.

source
GadgetIO.read_block!Function
read_block!(a::AbstractArray, f::IOStream,
            offset::Integer,
            nread::Integer, n_to_read::Integer;
            parttype::Integer,
            block_position::Integer,
            info::InfoLine,
            h::SnapshotHeader)

Read part of a block from a given file stream into a pre-allocated array.

source
GadgetIO.read_block!Method
read_block!(a::AbstractArray, f::IOStream,
            offset::Array{<:Integer},
            nread::Integer, n_to_read::Array{<:Integer};
            parttype::Integer,
            block_position::Integer,
            info::InfoLine,
            h::SnapshotHeader)

Read part of a block from a given file stream into a pre-allocated array.

source
GadgetIO.read_key_blockMethod
read_key_block( filename_keyfile::String, key_info::Vector{InfoLine}, 
                fields::Vector{String}, key::String, parttype::Integer)

Returns key block for different key. Valid values of key: KEY, NKEY, OKEY

source
GadgetIO.read_particles_in_box_peanoMethod
read_particles_in_box_peano(filename::String, blocks::Vector{String},
                            corner_lowerleft::Array{<:Real}, corner_upperright::Array{<:Real};
                            parttype::Integer=0, verbose::Bool=true)

Reads all particles within a box defined by a lower left and upper right corner for a given particle type based on peano hilbert key reading. Returns a dictionary with all requested blocks.

source
GadgetIO.read_pidsMethod
read_pids(sub_base::String, offset::Integer, N_ids::Integer)

Reads the PID block in the subfind output.

source
GadgetIO.read_positions_from_PH_keysMethod
read_positions_from_PH_keys(filebase::String,
                            corner_lowerleft::Array{<:Real}, 
                            corner_upperright::Array{<:Real};
                            parttype::Integer, verbose::Bool)

Finds the Peano-Hilbert keys in a cube defined by corner_lowerleft and corner_upperright.

source
GadgetIO.read_positions_from_keys_filesMethod
read_positions_from_PH_keys(files::Vector{<:Integer}, filebase::String, 
                            blocks::AbstractVector{String}, parttype::Integer,
                            keylist::Vector{<:Integer}, key_info::Vector{InfoLine},
                            verbose::Bool)

Helper function to get positions and length of particle blocks in files.

source
GadgetIO.read_subfind_lengthMethod
read_subfind_length(filename::String, blockname::String)

Get number of entries for block blockname. Uses the header, so it is faster than reading the whole block.

source
GadgetIO.reduce_read_positionsMethod
reduce_read_positions(sel::Array{<:Integer})

Reduces the individual read positions by finding sub-ranges. Returns an array of read indices and an array with the number of entries per index.

source
GadgetIO.select_fileMethod
select_file(filebase::String, filenum::Integer)

Checks if the requested files exists and returns the path.

source
GadgetIO.shift_across_box_borderMethod
shift_across_box_border(x::Real, x_halo::Real, boxsize::Real, boxsize_half::Real)

Shift coordinate x across the box border if the zero coordinate x₀ is on the other side.

source

Private Types

GadgetIO.KeyHeaderType
struct KeyHeader

Helper struct to store header information of .key files.

Fields

NameMeaning
nkeys_fileNumber of keys in this .key file
domain_cornersCorners of the domains defined by these keys
domain_facFactor needed for reconstructung int positions
bitsSize of PH keys it bits
nkeys_totalTotal number of keys in all files
nkeys_total_highwordTotal number of keys in all files
source