Types
Abstract type hierarchies

- HydroMapsType <: DataMapsType
- PartMapsType <: DataMapsType
List of types
Mera.ArgumentsTypeMera.CheckOutputNumberTypeMera.ClumpDataTypeMera.CompilationInfoTypeMera.ContainMassDataSetTypeMera.DataMapsTypeMera.DataSetTypeMera.DescriptorTypeMera.GridInfoTypeMera.Histogram2DMapTypeMera.HydroDataTypeMera.HydroMapsTypeMera.HydroPartTypeMera.InfoTypeMera.MaskTypeMera.PartDataTypeMera.PartInfoTypeMera.PartMapsTypeMera.PhysicalUnitsTypeMera.PhysicalUnitsType001Mera.ScalesTypeMera.ScalesType001Mera.WStatType
Functions
Mera.amroverviewMera.amroverviewMera.amroverviewMera.average_velocityMera.bellMera.bulk_velocityMera.center_of_massMera.center_of_massMera.checkoutputsMera.checksimulationsMera.comMera.comMera.construct_datatypeMera.convertdataMera.createpathMera.createscales!Mera.dataoverviewMera.dataoverviewMera.dataoverviewMera.dataoverviewMera.export_vtkMera.export_vtkMera.getclumpsMera.getextentMera.getgravityMera.gethydroMera.getinfoMera.getmassMera.getparticlesMera.getpositionsMera.gettimeMera.getvarMera.getvelocitiesMera.humanizeMera.infodataMera.loaddataMera.makefileMera.msumMera.notifymeMera.patchfileMera.prep_cylindrical_shellrangesMera.prep_spherical_shellrangesMera.preprangesMera.preprangesMera.printtimeMera.projectionMera.projectionMera.savedataMera.shellregionMera.storageoverviewMera.subregionMera.timerfileMera.usedmemoryMera.viewallfieldsMera.viewdataMera.viewfieldsMera.viewmoduleMera.wstat
Macros
Mera.@apply — MacroFind examples in the Mera Documentation for: filter data with pipeline macros
Mera.@filter — MacroFind examples in the Mera Documentation for: filter data with macros
Mera.@where — MacroFind examples in the Mera Documentation for: filter data with pipeline macros
Documentation Types
Mera.ArgumentsType — TypeMutable Struct: Contains fields to use as arguments in functions
Mera.CheckOutputNumberType — TypeMutable Struct: Contains the existing simulation snapshots in a folder and a list of the empty output-folders
Mera.ClumpDataType — TypeMutable Struct: Contains clump data and information about the selected simulation
ClumpDataType <: ContainMassDataSetType
Mera.CompilationInfoType — TypeMutable Struct: Contains the collected information about the compilation of RAMSES
Mera.ContainMassDataSetType — TypeAbstract Supertype of all datasets that contain mass variables
HydroPartType <: ContainMassDataSetType <: DataSetType
Mera.DataMapsType — TypeAbstract Supertype of all the different dataset type maps HydroMapsType <: DataMapsType PartMapsType <: DataMapsType
Mera.DataSetType — TypeAbstract Supertype of all the different dataset types
HydroPartType <: ContainMassDataSetType <: DataSetType
Mera.DescriptorType — TypeMutable Struct: Contains the collected information about the descriptors
Mera.GridInfoType — TypeMutable Struct: Contains the collected information about grid
Mera.Histogram2DMapType — TypeMutable Struct: Contains the 2D histogram returned by the function: histogram2 and information about the selected simulation
Mera.HydroDataType — TypeMutable Struct: Contains hydro data and information about the selected simulation
HydroDataType <: HydroPartType
Mera.HydroMapsType — TypeMutable Struct: Contains the maps/units returned by the hydro-projection information about the selected simulation
Mera.HydroPartType — TypeAbstract Supertype of data-sets that contain hydro and particle data
HydroPartType <: ContainMassDataSetType <: DataSetType
Mera.InfoType — TypeMutable Struct: Collected information about the selected simulation output
Mera.PartDataType — TypeMutable Struct: Contains particle data and information about the selected simulation
PartDataType <: HydroPartType
Mera.PartInfoType — TypeMutable Struct: Contains the collected information about particles
Mera.PartMapsType — TypeMutable Struct: Contains the maps/units returned by the particles-projection information about the selected simulation
Mera.PhysicalUnitsType — TypeMutable Struct: Contains the physical constants in cgs units
Mera.PhysicalUnitsType001 — TypeMutable Struct: Contains the physical constants in cgs units
Mera.ScalesType — TypeMutable Struct: Contains the created scale factors from code to physical units
Mera.ScalesType001 — TypeMutable Struct: Contains the created scale factors from code to physical units
Mera.WStatType — TypeMutable Struct: Contains the output statistics returned by wstat
Mera.MaskType — TypeUnion Type: Mask-array that is of type Bool or BitArray MaskType = Union{Array{Bool,1},BitArray{1}}
Documentation Functions
Mera.amroverview — MethodGet the number of cells and/or the CPUs per level
function overview_amr(dataobject::GravDataType, verbose::Bool=true)
function overview_amr(dataobject::GravDataType; verbose::Bool=true)
return a JuliaDB tableMera.amroverview — MethodGet the number of cells and/or the CPUs per level
function overview_amr(dataobject::HydroDataType, verbose::Bool=true)
function overview_amr(dataobject::HydroDataType; verbose::Bool=true)
return a JuliaDB tableMera.amroverview — MethodGet the number of particles and/or the CPUs per level
function overview_amr(dataobject::PartDataType, verbose::Bool=true)
function overview_amr(dataobject::PartDataType; verbose::Bool=true)
return a JuliaDB tableMera.average_velocity — MethodCalculate the average velocity (w/o mass-weight) of any ContainMassDataSetType:
average_velocity(dataobject::ContainMassDataSetType; unit::Symbol=:standard, weighting::Symbol=:mass, mask::MaskType=[false])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "ContainMassDataSetType"
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units) :kms, :ms, :cm_s (of typye Symbol) ..etc. ; see for defined velocity-scales viewfields(info.scale)weighting: use different weightings: :mass (default), :volume (hydro), :nomask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Mera.bell — MethodGet a notification sound, e.g., when your calculations are finished.
This may not apply when working remotely on a server:
julia> bell()Mera.bulk_velocity — MethodCalculate the average velocity (w/o mass-weight) of any ContainMassDataSetType:
bulk_velocity(dataobject::ContainMassDataSetType; unit::Symbol=:standard, weighting::Symbol=:mass, mask::MaskType=[false])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "ContainMassDataSetType"
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units) :kms, :ms, :cm_s (of typye Symbol) ..etc. ; see for defined velocity-scales viewfields(info.scale)weighting: use different weightings: :mass (default), :volume (hydro), :nomask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Mera.center_of_mass — MethodCalculate the center-of-mass of any ContainMassDataSetType:
center_of_mass(dataobject::ContainMassDataSetType; unit::Symbol=:standard, mask::MaskType=[false])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "ContainMassDataSetType"
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)mask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Mera.center_of_mass — MethodCalculate the joint center-of-mass of any HydroPartType:
center_of_mass(dataobject::Array{HydroPartType,1}, unit::Symbol; mask::MaskArrayAbstractType=[[false],[false]])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "Array{HydroPartType,1}""
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)mask: needs to be of type MaskArrayAbstractType which contains two entries with supertype of Array{Bool,1} or BitArray{1} and the length of the database (rows)
Mera.checkoutputs — FunctionGet the existing simulation snapshots in a given folder
- returns field
outputswith Array{Int,1} containing the output-numbers of the existing simulations - returns field
misswith Array{Int,1} containing the output-numbers of empty simulation folders - returns field
pathas String
checkoutputs(path::String="./"; verbose::Bool=true)
return CheckOutputNumberTypeExamples
# Example 1:
# look in current folder
julia> N = checkoutputs();
julia> N.outputs
julia> N.miss
julia> N.path
# Example 2:
# look in given path
# without any keyword
julia>N = checkoutputs("simulation001");Mera.checksimulations — FunctionList the existing simulation snapshots in a given folder
- returns a Dictonary with existing simulations and their folder names with:
- field
outputswith Array{Int,1} containing the output-numbers of the existing simulations - field
misswith Array{Int,1} containing the output-numbers of empty simulation folders - field
pathas String
checksimulations(path::String="./"; verbose::Bool=true, filternames=String[])
return Dict with named Tuple: simulation name and N=CheckOutputNumberTypeMera.com — MethodCalculate the center-of-mass of any ContainMassDataSetType:
com(dataobject::ContainMassDataSetType; unit::Symbol=:standard, mask::MaskType=[false])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "ContainMassDataSetType"
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)mask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Mera.com — MethodCalculate the joint center-of-mass of any HydroPartType:
com(dataobject::Array{HydroPartType,1}, unit::Symbol; mask::MaskArrayAbstractType=[[false],[false]])
return Tuple{Float64, Float64, Float64,}Arguments
Required:
dataobject: needs to be of type: "Array{HydroPartType,1}""
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)mask: needs to be of type MaskArrayAbstractType which contains two entries with supertype of Array{Bool,1} or BitArray{1} and the length of the database (rows)
Mera.construct_datatype — MethodCreate a New DataSetType from a Filtered Data Table
function construct_datatype(data::IndexedTables.AbstractIndexedTable, dataobject::HydroDataType)
return HydroDataType
function construct_datatype(data::IndexedTables.AbstractIndexedTable, dataobject::PartDataType)
return PartDataType
function construct_datatype(data::IndexedTables.AbstractIndexedTable, dataobject::ClumpDataType)
return ClumpDataType
function construct_datatype(data::IndexedTables.AbstractIndexedTable, dataobject::GravDataType)
return GravDataTypeExample
# read simulation information
julia> info = getinfo(420)
julia> gas = gethydro(info)
# filter and create a new` data table
julia> density = 3. /gas.scale.Msol_pc3
julia> filtered_db = @filter gas.data :rho >= density
# construct a new HydroDataType
# (comparable to the object "gas" but only with filtered data)
julia> gas_new = construct_datatype(filtered_db, gas)Mera.convertdata — MethodConverts full simulation data into a compressed/uncompressed JLD2 format:
- all existing datatypes are sequentially loaded and stored
- supports :hydro, :particles, :hydro, :clumps
- running number is taken from original RAMSES folders
- use different compression methods
- select a certain data range, smallr, smallc, lmax
- add a string to describe the simulation
- the individual loading and storing processes are timed and stored into the file (and more statistics)
- toggle progressbar mode
- toggle verbose mode
function convertdata(output::Int; datatypes::Array{<:Any,1}=[missing], path::String="./", fpath::String="./",
fname = "output_",
compress::Any=nothing,
comments::Any=nothing,
lmax::Union{Int, Missing}=missing,
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
smallr::Real=0.,
smallc::Real=0.,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )
return statistics (dictionary)
Arguments
Required:
output: output number
Predefined/Optional Keywords:
datatypes: default -> all available (known) data is converted; pass an array with only selected datatypes, e.g.: datatypes=[:hydro, :particles]fname: default name of the files "output_" and the running number is added. Change the string to apply a user-defined name.compress: by default compression is activated. compress=false (deactivate).
If necessary, choose between different compression types: LZ4FrameCompressor() (default), Bzip2Compressor(), ZlibCompressor(). Load the required package to choose the compression type and to see their parameters: CodecZlib, CodecBzip2 or CodecLz4
comments: add a string that includes e.g. a description about your simulationlmax: the maximum level to be read from the datapath: path to the RAMSES folders; default is local path.fpath: path to the JLD23 file; default is local path.xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..smallr: set lower limit for density; zero means inactivesmallc: set lower limit for thermal pressure; zero means inactivemyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of xrange, yrange, zrange, center, range_unit, verboseverbose: print timestamp and further information on screen; default: true
Defined Methods - function defined for different arguments
- convertdata(output::Int64; ...)
- convertdata(output::Int64, datatypes::Vector{Symbol}; ...)
- convertdata(output::Int64, datatypes::Symbol; ...)
Mera.createpath — Method```julia createpath(output::Real, path::String; namelist::String="")
return FileNamesType ```
Mera.createscales! — MethodCreate an object with predefined scale factors from code to pysical units
function createscales!(dataobject::InfoType)
return ScalesType001Mera.dataoverview — MethodGet the extrema of each variable in the database
function dataoverview(dataobject::ClumpDataType)
return a JuliaDB tableMera.dataoverview — Methodtodo: needs to be tested
Get the total epot and min/max value of each variable in the database per level
function dataoverview(dataobject::GravDataType, verbose::Bool=true)
function dataoverview(dataobject::GravDataType; verbose::Bool=true)
return a JuliaDB tableMera.dataoverview — MethodGet the mass and min/max value of each variable in the database per level
function dataoverview(dataobject::HydroDataType, verbose::Bool=true)
function dataoverview(dataobject::HydroDataType; verbose::Bool=true)
return a JuliaDB tableMera.dataoverview — MethodGet the min/max value of each variable in the database per level
function dataoverview(dataobject::PartDataType, verbose::Bool=true)
function dataoverview(dataobject::PartDataType; verbose::Bool=true)
return a JuliaDB tableMera.export_vtk — MethodExport hydro data to VTK format for visualization in tools like ParaView.
- export data that is present in your database and can be processed by getvar() (done internally)
- select scalar(s) and their unit(s)
- select a vector and its unit (like velocity)
- export data in log10
- creates binary files with optional compression
- supports multi-threading
-> generating per-level VTU files for scalar and optionally vector data and creates corresponding VTM multiblock container files to reference these VTU files.
export_vtk(
dataobject::HydroDataType, outprefix::String;
scalars::Vector{Symbol} = [:rho],
scalars_unit::Vector{Symbol} = [:nH],
scalars_log10::Bool=false,
vector::Array{<:Any,1}=[missing, missing, missing],
vector_unit::Symbol = :km_s,
vector_name::String = "velocity",
vector_log10::Bool=false,
positions_unit::Symbol = :standard,
lmin::Int = dataobject.lmin,
lmax::Int = dataobject.lmax,
chunk_size::Int = 50000,
compress::Bool = true,
interpolate_higher_levels::Bool = true,
max_cells::Int = 100_000_000,
verbose::Bool = true,
myargs::ArgumentsType=ArgumentsType()
)Arguments
Required:
dataobject::HydroDataType: The AMR data structure from MERA.jl containing variables like level, position, and physical quantities.outprefix::String: The base path and prefix for output files (e.g., "output/data" will create files like "output/data_L0.vtu").
Predefined/Optional Keywords:
scalars: List of scalar variables to export (default is :rho); from the database or a predefined quantity (see field: info, function getvar(), dataobject.data)scalars_unit: Sets the unit for the list of scalars (default is hydrogen number density in cm^-3).scalars_log10: Apply log10 to the scalars (default false).vector: List of vector component variables to export (default is missing); exports vector data as separate VTU filesvector_unit: Sets the unit for the vector components (default is km/s).vector_name: The name of the vector field in the VTK file (default: "velocity").vector_log10: Apply log10 to the vector components (default: false).positions_unit: Sets the unit of the cell positions (default: code units); usefull in paraview to select regionslmin: Minimum AMR level to process (default: simulations lmin); smaller levels are excluded in exportlmax: Maximum AMR level to process (default: simulations lmax); existing higher levels are interpolated down if interpolatehigherlevels is true, otherwise excluded from exportchunk_size::Int = 50000: Size of data chunks for processing (currently unused but reserved for future optimizations).chunk_size::Int = 50000: Size of data chunks for processing (currently unused but reserved for future optimizations).compress: Iftrue(default), enable compression.interpolate_higher_levels: Iftrue, interpolate data from higher levels down to givenlmax.max_cells: Maximum number of cells to export per level (caps output if exceeded, prioritizing denser regions), (default: 100000000)verbose: Iftrue(default), print detailed progress and diagnostic messages.
Mera.export_vtk — MethodExport particle data to VTK format for visualization in tools like ParaView.
- export data that is present in your database and can be processed by getvar() (done internally)
- select scalar(s) and their unit(s)
- select a vector and its unit (like velocity)
- export data in log10
- creates binary files with optional compression
- supports multi-threading
-> generates VTU files; each particle is represented as a vertex point with associated scalar and vector data.
export_vtk(
dataobject::PartDataType, outprefix::String;
scalars::Vector{Symbol} = [:mass],
scalars_unit::Vector{Symbol} = [:Msol],
scalars_log10::Bool=false,
vector::Array{<:Any,1}=[missing, missing, missing],
vector_unit::Symbol = :km_s,
vector_name::String = "velocity",
vector_log10::Bool=false,
positions_unit::Symbol = :standard,
chunk_size::Int = 50000,
compress::Bool = false,
max_particles::Int = 100_000_000,
verbose::Bool = true,
myargs::ArgumentsType=ArgumentsType()
)Arguments
Required:
- **
dataobject::PartDataType:*** needs to be of type "PartDataType" outprefix: The base path and prefix for output file (e.g., "foldername/particles" will create "foldername/particles.vtu").
Predefined/Optional Keywords:
scalars: List of scalar variables to export (default is particle mass); from the database or a predefined quantity (see field: info, function getvar(), dataobject.data)scalars_unit: Sets the unit for the list of scalars (default is Msun).scalars_log10: Apply log10 to the scalars (default false).vector: List of vector component variables to export (default is missing).vector_unit: Sets the unit for the vector components (default is km/s).vector_name: The name of the vector field in the VTK file (default: "velocity").vector_log10: Apply log10 to the vector components (default: false).positions_unit: Sets the unit of the particle positions (default: code units); usefull in paraview to select regionschunk_size::Int = 50000: Size of data chunks for processing (reserved for future optimizations).compress: Iffalse(default), disable compression.max_particles: Maximum number of particles to export (caps output if exceeded), (default: 100000000)verbose: Iftrue(default), print detailed progress and diagnostic messages.
Mera.getclumps — MethodRead the clump-data:
- selected variables
- limited to a spatial range
- print the name of each data-file before reading it
- toggle verbose mode
- pass a struct with arguments (myargs)
getclumps( dataobject::InfoType;
vars::Array{Symbol,1}=[:all],
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
print_filenames::Bool=false,
verbose::Bool=true,
myargs::ArgumentsType=ArgumentsType() )Returns an object of type ClumpDataType, containing the clump-data table, the selected options and the simulation ScaleType and summary of the InfoType
return ClumpDataType()
# get an overview of the returned fields:
# e.g.:
julia> info = getinfo(100)
julia> clumps = getclumps(info)
julia> viewfields(clumps)
#or:
julia> fieldnames(clumps)Arguments
Required:
dataobject: needs to be of type: "InfoType", created by the function getinfo
Predefined/Optional Keywords:
vars: Currently, the length of the loaded variable list can be modified *(see examples below).xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..print_filenames: print on screen the current processed particle file of each CPUverbose: print timestamp, selected vars and ranges on screen; default: truemyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of xrange, yrange, zrange, center, range_unit, verbose
Defined Methods - function defined for different arguments
- getclumps(dataobject::InfoType; ...) # no given variables -> all variables loaded
- getclumps(dataobject::InfoType, vars::Array{Symbol,1}; ...) # one or several given variables -> array needed
Examples
# read simulation information
julia> info = getinfo(420)
# Example 1:
# read clump data of all variables, full-box
julia> clumps = getclumps(info)
# Example 2:
# read clump data of all variables
# data range 20x20x4 kpc; ranges are given in kpc relative to the box (here: 48 kpc) center at 24 kpc
julia> clumps = getclumps( info,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[24., 24., 24.],
range_unit=:kpc )
# Example 3:
# give the center of the box by simply passing: center = [:bc] or center = [:boxcenter]
# this is equivalent to center=[24.,24.,24.] in Example 2
# the following combination is also possible: e.g. center=[:bc, 12., 34.], etc.
julia> clumps = getclumps( info,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[33., bc:, 10.],
range_unit=:kpc )
# Example 4:
# Load less than the found 12 columns from the header of the clump files;
# Pass an array with the variables to the keyword argument *vars*.
# The order of the variables has to be consistent with the header in the clump files:
julia> lumps = getclumps(info, [ :index, :lev, :parent, :ncell,
:peak_x, :peak_y, :peak_z ])
# Example 5:
# Load more than the found 12 columns from the header of the clump files.
# E.g. the list can be extended with more names if there are more columns
# in the data than given by the header in the files.
# The order of the variables has to be consistent with the header in the clump files:
julia> clumps = getclumps(info, [ :index, :lev, :parent, :ncell,
:peak_x, :peak_y, :peak_z,
Symbol("rho-"), Symbol("rho+"),
:rho_av, :mass_cl, :relevance,
:vx, :vy, :vz ])
...Mera.getextent — MethodGet the extent of the dataset-domain:
function getextent( dataobject::DataSetType;
unit::Symbol=:standard,
center::Array{<:Any,1}=[0., 0., 0.],
center_unit::Symbol=:standard,
direction::Symbol=:z)
return (xmin, xmax), (ymin ,ymax ), (zmin ,zmax )Arguments
Required:
dataobject: needs to be of type: "DataSetType"
Predefined/Optional Keywords:
center: in unit given by argumentcenter_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..center_unit: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)direction: todounit: return the variables in given unit
Defined Methods - function defined for different arguments
- getextent( dataobject::DataSetType; # one given variable
- getextent( dataobject::DataSetType, unit::Symbol; ...) # one given variable with its unit
Mera.getgravity — MethodRead the leaf-cells of the gravity-data:
- select variables
- limit to a maximum level
- limit to a spatial range
- print the name of each data-file before reading it
- toggle verbose mode
- toggle progress bar
- pass a struct with arguments (myargs)
getgravity( dataobject::InfoType;
lmax::Real=dataobject.levelmax,
vars::Array{Symbol,1}=[:all],
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
print_filenames::Bool=false,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )Returns an object of type GravDataType, containing the gravity-data table, the selected options and the simulation ScaleType and summary of the InfoType
return GravDataType()
# get an overview of the returned fields:
# e.g.:
julia> info = getinfo(100)
julia> grav = getgravity(info)
julia> viewfields(grav)
#or:
julia> fieldnames(grav)Arguments
Required:
dataobject: needs to be of type: "InfoType", created by the function getinfo
Predefined/Optional Keywords:
lmax: the maximum level to be read from the datavar(s): the selected gravity variables in arbitrary order: :all (default), :cpu, :epot, :ax, :ay, :azxrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..print_filenames: print on screen the current processed gravity file of each CPUverbose: print timestamp, selected vars and ranges on screen; default: trueshow_progress: print progress bar on screenmyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of lmax, xrange, yrange, zrange, center, rangeunit, verbose, showprogress
Defined Methods - function defined for different arguments
- getgravity( dataobject::InfoType; ...) # no given variables -> all variables loaded
- getgravity( dataobject::InfoType, var::Symbol; ...) # one given variable -> no array needed
- getgravity( dataobject::InfoType, vars::Array{Symbol,1}; ...) # several given variables -> array needed
Examples
# read simulation information
julia> info = getinfo(420)
# Example 1:
# read gravity data of all variables, full-box, all levels
julia> grav = getgravity(info)
# Example 2:
# read gravity data of all variables up to level 8
# data range 20x20x4 kpc; ranges are given in kpc relative to the box (here: 48 kpc) center at 24 kpc
julia> grav = getgravity( info,
lmax=8,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[24., 24., 24.],
range_unit=:kpc )
# Example 3:
# give the center of the box by simply passing: center = [:bc] or center = [:boxcenter]
# this is equivalent to center=[24.,24.,24.] in Example 2
# the following combination is also possible: e.g. center=[:bc, 12., 34.], etc.
julia> grav = getgravity( info,
lmax=8,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[33., bc:, 10.],
range_unit=:kpc )
# Example 4:
# read gravity data of the variables epot and the x-acceleration, full-box, all levels
julia> grav = getgravity( info, [:epot, :ax] ) # use array for the variables
# Example 5:
# read gravity data of the single variable epot, full-box, all levels
julia> grav = getgravity( info, :epot ) # no array for a single variable needed
...Mera.gethydro — MethodRead the leaf-cells of the hydro-data:
- select variables
- limit to a maximum level
- limit to a spatial range
- set a minimum density or sound speed
- check for negative values in density and thermal pressure
- print the name of each data-file before reading it
- toggle verbose mode
- toggle progress bar
- pass a struct with arguments (myargs)
gethydro( dataobject::InfoType;
lmax::Real=dataobject.levelmax,
vars::Array{Symbol,1}=[:all],
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
smallr::Real=0.,
smallc::Real=0.,
check_negvalues::Bool=false,
print_filenames::Bool=false,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )Returns an object of type HydroDataType, containing the hydro-data table, the selected options and the simulation ScaleType and summary of the InfoType
return HydroDataType()
# get an overview of the returned fields:
# e.g.:
julia> info = getinfo(100)
julia> gas = gethydro(info)
julia> viewfields(gas)
#or:
julia> fieldnames(gas)Arguments
Required:
dataobject: needs to be of type: "InfoType", created by the function getinfo
Predefined/Optional Keywords:
lmax: the maximum level to be read from the datavar(s): the selected hydro variables in arbitrary order: :all (default), :cpu, :rho, :vx, :vy, :vz, :p, :var6, :var7...xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..smallr: set lower limit for density; zero means inactivesmallc: set lower limit for thermal pressure; zero means inactivecheck_negvalues: check loaded data of "rho" and "p" on negative values; false by defaultprint_filenames: print on screen the current processed hydro file of each CPUverbose: print timestamp, selected vars and ranges on screen; default: trueshow_progress: print progress bar on screenmyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of lmax, xrange, yrange, zrange, center, rangeunit, verbose, showprogress
Defined Methods - function defined for different arguments
- gethydro( dataobject::InfoType; ...) # no given variables -> all variables loaded
- gethydro( dataobject::InfoType, var::Symbol; ...) # one given variable -> no array needed
- gethydro( dataobject::InfoType, vars::Array{Symbol,1}; ...) # several given variables -> array needed
Examples
# read simulation information
julia> info = getinfo(420)
# Example 1:
# read hydro data of all variables, full-box, all levels
julia> gas = gethydro(info)
# Example 2:
# read hydro data of all variables up to level 8
# data range 20x20x4 kpc; ranges are given in kpc relative to the box (here: 48 kpc) center at 24 kpc
julia> gas = gethydro( info,
lmax=8,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[24., 24., 24.],
range_unit=:kpc )
# Example 3:
# give the center of the box by simply passing: center = [:bc] or center = [:boxcenter]
# this is equivalent to center=[24.,24.,24.] in Example 2
# the following combination is also possible: e.g. center=[:bc, 12., 34.], etc.
julia> gas = gethydro( info,
lmax=8,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[33., bc:, 10.],
range_unit=:kpc )
# Example 4:
# read hydro data of the variables density and the thermal pressure, full-box, all levels
julia> gas = gethydro( info, [:rho, :p] ) # use array for the variables
# Example 5:
# read hydro data of the single variable density, full-box, all levels
julia> gas = gethydro( info, :rho ) # no array for a single variable needed
...Mera.getinfo — MethodGet the simulation overview from RAMSES info, descriptor and output header files
getinfo(; output::Real=1, path::String="", namelist::String="", verbose::Bool=true)
return InfoTypeKeyword Arguments
output: timestep number (default=1)path: the path to the output folder relative to the current folder or absolute pathnamelist: give the path to a namelist file (by default the namelist.txt-file in the output-folder is read)verbose:: informations are printed on the screen by default
Examples
# read simulation information from output `1` in current folder
julia> info = getinfo()
# read simulation information from output `420` in given folder (relative path to the current working folder)
julia> info = getinfo(output=420, path="../MySimFolder/")
# or simply use
julia> info = getinfo(420, "../MySimFolder/")
# get an overview of the returned field-names
julia> propertynames(info)
# a more detailed overview
julia> viewfields(info)
...
julia> viewallfields(info)
...
julia> namelist(info)
...
julia> makefile(info)
...
julia> timerfile(info)
...
julia> patchfile(info)
...Mera.getmass — MethodGet mass-array from the dataset (cells/particles/clumps/...):
getmass(dataobject::HydroDataType)
getmass(dataobject::PartDataType)
getmass(dataobject::ClumpDataType)
return Array{Float64,1}Mera.getparticles — MethodRead the particle-data
- select variables
- limit to a maximum range
- print the name of each data-file before reading it
- toggle verbose mode
- toggle progress bar
- pass a struct with arguments (myargs)
getparticles( dataobject::InfoType;
lmax::Real=dataobject.levelmax,
vars::Array{Symbol,1}=[:all],
stars::Bool=true,
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
presorted::Bool=true,
print_filenames::Bool=false,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )Returns an object of type PartDataType, containing the particle-data table, the selected and the simulation ScaleType and summary of the InfoType
return PartDataType()
# get an overview of the returned fields:
# e.g.:
julia> info = getinfo(100)
julia> particles = getparticles(info)
julia> viewfields(particles)
#or:
julia> fieldnames(particles)Arguments
Required:
dataobject: needs to be of type: "InfoType", created by the function getinfo
Predefined/Optional Keywords:
lmax: not definedstars: not definedvar(s): the selected particle variables in arbitrary order: :all (default), :cpu, :mass, :vx, :vy, :vz, :birth :metals, ...xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..presorted: presort data according to the key vars (by default)print_filenames: print on screen the current processed particle file of each CPUverbose: print timestamp, selected vars and ranges on screen; default: trueshow_progress: print progress bar on screenmyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of lmax not!, xrange, yrange, zrange, center, rangeunit, verbose, showprogress
Defined Methods - function defined for different arguments
- getparticles( dataobject::InfoType; ...) # no given variables -> all variables loaded
- getparticles( dataobject::InfoType, var::Symbol; ...) # one given variable -> no array needed
- getparticles( dataobject::InfoType, vars::Array{Symbol,1}; ...) # several given variables -> array needed
Examples
# read simulation information
julia> info = getinfo(420)
# Example 1:
# read particle data of all variables, full-box, all levels
julia> particles = getparticles(info)
# Example 2:
# read particle data of all variables
# data range 20x20x4 kpc; ranges are given in kpc relative to the box (here: 48 kpc) center at 24 kpc
julia> particles = getparticles( info,
xrange=[-10., 10.],
yrange=[-10., 10.],
zrange=[-2., 2.],
center=[24., 24., 24.],
range_unit=:kpc )
# Example 3:
# give the center of the box by simply passing: center = [:bc] or center = [:boxcenter]
# this is equivalent to center=[24.,24.,24.] in Example 2
# the following combination is also possible: e.g. center=[:bc, 12., 34.], etc.
julia> particles = getparticles( info,
xrange=[-10.,10.],
yrange=[-10.,10.],
zrange=[-2.,2.],
center=[33., bc:, 10.],
range_unit=:kpc )
# Example 4:
# read particle data of the variables mass and the birth-time, full-box, all levels
julia> particles = getparticles( info, [:mass, :birth] ) # use array for the variables
# Example 5:
# read particle data of the single variable mass, full-box, all levels
julia> particles = getparticles( info, :mass ) # no array for a single variable needed
...Mera.getpositions — MethodGet the x,y,z positions from the dataset (cells/particles/clumps/...):
getpositions( dataobject::DataSetType, unit::Symbol;
direction::Symbol=:z,
center::Array{<:Any,1}=[0., 0., 0.],
center_unit::Symbol=:standard,
mask::MaskType=[false])
return x, y, zArguments
Required:
dataobject: needs to be of type: "DataSetType"
Predefined/Optional Keywords:
center: in unit given by argumentcenter_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..center_unit: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)direction: todounit: return the variables in given unitmask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Defined Methods - function defined for different arguments
- getpositions( dataobject::DataSetType; ...) # one given dataobject
- getpositions( dataobject::DataSetType, unit::Symbol; ...) # one given dataobject and position unit
Mera.gettime — MethodGet physical time in selected units
returns Float
gettime(output::Real; path::String="./", unit::Symbol=:standard)
gettime(dataobject::DataSetType; unit::Symbol=:standard)
gettime(dataobject::InfoType, unit::Symbol=:standard)
return timeArguments Function 1
Required:
output: give the output-number of the simulation
Predefined/Optional Keywords:
path: the path to the output folder relative to the current folder or absolute pathunit: return the variable in given unit
Arguments Function 2
Required:
dataobject: needs to be of type: "DataSetType"
Predefined/Optional Keywords:
unit: return the variable in given unit
Arguments Function 3
Required:
dataobject: needs to be of type: "InfoType"
Predefined/Optional Keywords:
unit: return the variable in given unit
Mera.getvar — MethodGet variables or derived quantities from the dataset:
- overview the list of predefined quantities with: getinfo()
- select variable(s) and their unit(s)
- give the spatial center (with units) of the data within the box (relevant e.g. for radius dependency)
- relate the coordinates to a direction (x,y,z)
- pass a modified database
- pass a mask to exclude elements (cells/particles/...) from the calculation
getvar( dataobject::DataSetType, var::Symbol;
filtered_db::IndexedTables.AbstractIndexedTable=IndexedTables.table([1]),
center::Array{<:Any,1}=[0.,0.,0.],
center_unit::Symbol=:standard,
direction::Symbol=:z,
unit::Symbol=:standard,
mask::MaskType=[false],
ref_time::Real=dataobject.info.time)
return Array{Float64,1}
Arguments
Required:
dataobject: needs to be of type: "DataSetType"var(s): select a variable from the database or a predefined quantity (see field: info, function getvar(), dataobject.data)
Predefined/Optional Keywords:
filtered_db: pass a filtered or manipulated database together with the corresponding DataSetType object (required argument)center: in units given by argumentcenter_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..center_unit: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)direction: todounit(s): return the variable in given unitsmask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)ref_time: reference zero-time for particle age calculation
Defined Methods - function defined for different arguments
- getvar( dataobject::DataSetType, var::Symbol; ...) # one given variable -> returns 1d array
- getvar( dataobject::DataSetType, var::Symbol, unit::Symbol; ...) # one given variable with its unit -> returns 1d array
- getvar( dataobject::DataSetType, vars::Array{Symbol,1}; ...) # several given variables -> array needed -> returns dictionary with 1d arrays
- getvar( dataobject::DataSetType, vars::Array{Symbol,1}, units::Array{Symbol,1}; ...) # several given variables and their corresponding units -> both arrays -> returns dictionary with 1d arrays
- getvar( dataobject::DataSetType, vars::Array{Symbol,1}, unit::Symbol; ...) # several given variables that have the same unit -> array for the variables and a single Symbol for the unit -> returns dictionary with 1d arrays
Examples
```julia
read simulation information
julia> info = getinfo(420) julia> gas = gethydro(info)
Example 1: get the mass for each cell of the hydro data (1dim array)
mass1 = getvar(gas, :mass) # in [code units] mass = getvar(gas, :mass) * gas.scale.Msol # scale the result from code units to solar masses mass = getvar(gas, :mass, unit=:Msol) # unit calculation, provided by a keyword argument mass = getvar(gas, :mass, :Msol) # unit calculation provided by an argument
Example 2: get the mass and |v| (several variables) for each cell of the hydro data
quantities = getvar(gas, [:mass, :v]) # in [code units] returns: Dict{Any,Any} with 2 entries: :mass => [8.9407e-7, 8.9407e-7, 8.9407e-7, 8.9407e-7, 8.9407e-7, 8.9407e-7, 8… :v => [0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 … 2.28274e-7, 2.…
quantities = getvar(gas, [:mass, :v], units=[:Msol, :kms]) # unit calculation, provided by a keyword argument quantities = getvar(gas, [:mass, :v], [:Msol, :kms]) # unit calculation provided by an argument
Example 3: get several variables in the same units by providing a single argument
quantities = getvar(gas, [:vx, :vy, :vz], :km_s) ...
Mera.getvelocities — MethodGet the vx,vy,vz velocities from the dataset (cells/particles/clumps/...):
function getvelocities( dataobject::DataSetType, unit::Symbol;
mask::MaskType=[false])
return vx, vy, vzArguments
Required:
dataobject: needs to be of type: "DataSetType"
Predefined/Optional Keywords:
unit: return the variables in given unitmask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Defined Methods - function defined for different arguments
- getvelocities( dataobject::DataSetType; ...) # one given dataobject
- getvelocities( dataobject::DataSetType, unit::Symbol; ...) # one given dataobject and velocity unit
Mera.humanize — MethodConvert a value to human-readable astrophysical units and round to ndigits
(pass the value in code units and the quantity specification (length, time) )
function humanize(value::Float64, scale::ScalesType001, ndigits::Int, quantity::String)
return value, value_unitMera.infodata — MethodGet the simulation overview from RAMSES, saved in JLD2 == function getinfo
infodata(output::Int;
path::String="./",
fname = "output_",
datatype::Any=:nothing,
verbose::Bool=true)
return InfoTypeKeyword Arguments
output: timestep numberpath: the path to the output JLD2 file relative to the current folder or absolute pathfname: "output"-> filename = "output***.jld2" by default, can be changed to "myname***.jld2"verbose:: informations are printed on the screen by default
Examples
# read simulation information from output `1` in current folder
julia> info = infodata(1) # filename="output_00001.jld2"
# read simulation information from output `420` in given folder (relative path to the current working folder)
julia> info = infodata(420, path="../MySimFolder/")
# or simply use
julia> info = infodata(420, "../MySimFolder/")
# get an overview of the returned field-names
julia> propertynames(info)
# a more detailed overview
julia> viewfields(info)
...
julia> viewallfields(info)
...
julia> namelist(info)
...
julia> makefile(info)
...
julia> timerfile(info)
...
julia> patchfile(info)
...Mera.loaddata — MethodRead stored simulation data into a dataobject:
- supported datatypes: HydroDataType, PartDataType, GravDataType, ClumpDataType
- select a certain data range (data is fully loaded; the selected subregion is returned)
- toggle verbose mode
function loaddata(output::Int; path::String="./",
fname = "output_",
datatype::Symbol,
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
verbose::Bool=true,
myargs::ArgumentsType=ArgumentsType() )
return dataobject
Arguments
Required:
output: output numberdatatype: :hydro, :particles, :gravity or :clumps
Predefined/Optional Keywords:
path: path to the file; default is local path.fname: default name of the files "output_" and the running number is added. Change the string to apply a user-defined name.xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)myargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of xrange, yrange, zrange, center, range_unit, verboseverbose: print timestamp and further information on screen; default: true
Defined Methods - function defined for different arguments
- loaddata(output::Int64; ...) # opens first datatype in the file
- loaddata(output::Int64, datatype::Symbol; ...)
- loaddata(output::Int64, path::String; ...)
- loaddata(output::Int64, path::String, datatype::Symbol; ...)
Mera.makefile — MethodGet a printout of the makefile:
makefile(object::InfoType)Mera.msum — MethodCalculate the total mass of any ContainMassDataSetType:
msum(dataobject::ContainMassDataSetType; unit::Symbol=:standard, mask::MaskType=[false])
return Float64Arguments
Required:
dataobject: needs to be of type: "ContainMassDataSetType"
Optional Keywords:
unit: the unit of the result (can be used w/o keyword): :standard (code units) :Msol, :Mearth, :Mjupiter, :g, :kg (of typye Symbol) ..etc. ; see for defined mass-scales viewfields(info.scale)mask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)
Mera.notifyme — MethodGet an email notification, e.g., when your calculations are finished.
Mandatory:
- the email client "mail" needs to be installed
- put a file with the name "email.txt" in your home folder that contains your email address in the first line
julia> notifyme()or:
julia> notifyme("Calculation 1 finished!")Mera.patchfile — MethodGet a printout of the patchfile:
patchfile(object::InfoType)Mera.prep_cylindrical_shellranges — Methodcylindrical shell ranges convert given ranges + radius and print overview on screen used for shellregions...
Mera.prep_spherical_shellranges — Methodspherical shell ranges convert given ranges + radius and print overview on screen used for shellregions...
Mera.prepranges — Methodcuboid ranges convert given ranges and print overview on screen used for gethydro, getparticles, getgravity, getclumps..., subregions...
Mera.prepranges — Methodcylinder ranges (height != 0.), sphere ranges (height==0.) convert given ranges + radius and print overview on screen used for subregions...
Mera.printtime — Functionprint a Mera timestamp on the screen if the global variable: verbose_mode == true
function printtime(text::String="", verbose::Bool=verbose_mode)Mera.projection — MethodProject variables or derived quantities from the hydro-dataset:
- projection to an arbitrary large grid: give pixelnumber for each dimension = res
- overview the list of predefined quantities with: projection()
- select variable(s) and their unit(s)
- limit to a maximum range
- select a coarser grid than the maximum resolution of the loaded data (maps with both resolutions are created)
- give the spatial center (with units) of the data within the box (relevant e.g. for radius dependency)
- relate the coordinates to a direction (x,y,z)
- select arbitrary weighting: mass (default), volume weighting, etc.
- pass a mask to exclude elements (cells) from the calculation
- toggle verbose mode
- toggle progress bar
- pass a struct with arguments (myargs)
projection( dataobject::HydroDataType, vars::Array{Symbol,1};
units::Array{Symbol,1}=[:standard],
lmax::Real=dataobject.lmax,
res::Union{Real, Missing}=missing,
pxsize::Array{<:Any,1}=[missing, missing],
mask::Union{Vector{Bool}, MaskType}=[false],
direction::Symbol=:z,
weighting::Array{<:Any,1}=[:mass, missing],
mode::Symbol=:standard,
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
data_center::Array{<:Any,1}=[missing, missing, missing],
data_center_unit::Symbol=:standard,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )
return HydroMapsType
Arguments
Required:
dataobject: needs to be of type: "HydroDataType"var(s): select a variable from the database or a predefined quantity (see field: info, function projection(), dataobject.data)
Predefined/Optional Keywords:
unit(s): return the variable in given unitspxsize`: creates maps with the given pixel size in physical/code units (dominates over: res, lmax) : pxsize=[physical size (Number), physical unit (Symbol)]rescreate maps with the given pixel number for each deminsion; if res not given by user -> lmax is selected; (pixel number is related to the full boxsize)lmax: create maps with 2^lmax pixels for each dimensionxrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..weighting: select between:massweighting (default) and any other pre-defined quantity, e.g.:volume. Pass an array with the weighting=[quantity (Symbol), physical unit (Symbol)]data_center: to calculate the data relative to the datacenter; in units given by argument `datacenterunit`; by default the argument datacenter = center ;data_center_unit: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)direction: select between: :x, :y, :zmask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)mode: :standard (default) handles projections other than surface density. mode=:standard (default) -> weighted average; mode=:sum sums-up the weighted quantities in projection direction.show_progress: print progress bar on screenmyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of lmax, xrange, yrange, zrange, center, rangeunit, verbose, showprogress
Defined Methods - function defined for different arguments
- projection( dataobject::HydroDataType, var::Symbol; ...) # one given variable
- projection( dataobject::HydroDataType, var::Symbol, unit::Symbol; ...) # one given variable with its unit
- projection( dataobject::HydroDataType, vars::Array{Symbol,1}; ...) # several given variables -> array needed
- projection( dataobject::HydroDataType, vars::Array{Symbol,1}, units::Array{Symbol,1}; ...) # several given variables and their corresponding units -> both arrays
- projection( dataobject::HydroDataType, vars::Array{Symbol,1}, unit::Symbol; ...) # several given variables that have the same unit -> array for the variables and a single Symbol for the unit
Examples
...
Mera.projection — MethodProject variables or derived quantities from the particle-dataset:
- projection to a grid related to a given level
- overview the list of predefined quantities with: projection()
- select variable(s) and their unit(s)
- limit to a maximum range
- give the spatial center (with units) of the data within the box (relevant e.g. for radius dependency)
- relate the coordinates to a direction (x,y,z)
- select between mass (default) and volume weighting
- pass a mask to exclude elements (cells/particles/...) from the calculation
- toggle verbose mode
- toggle progress bar
- pass a struct with arguments (myargs)
projection( dataobject::PartDataType, vars::Array{Symbol,1};
units::Array{Symbol,1}=[:standard],
lmax::Real=dataobject.lmax,
res::Union{Real, Missing}=missing,
pxsize::Array{<:Any,1}=[missing, missing],
mask=[false],
direction::Symbol=:z,
weighting::Symbol=:mass,
xrange::Array{<:Any,1}=[missing, missing],
yrange::Array{<:Any,1}=[missing, missing],
zrange::Array{<:Any,1}=[missing, missing],
center::Array{<:Any,1}=[0., 0., 0.],
range_unit::Symbol=:standard,
data_center::Array{<:Any,1}=[missing, missing, missing],
data_center_unit::Symbol=:standard,
ref_time::Real=dataobject.info.time,
verbose::Bool=true,
show_progress::Bool=true,
myargs::ArgumentsType=ArgumentsType() )
return HydroMapsType
Arguments
Required:
dataobject: needs to be of type: "PartDataType"var(s): select a variable from the database or a predefined quantity (see field: info, function projection(), dataobject.data)
Predefined/Optional Keywords:
unit(s): return the variable in given unitspxsize`: creates maps with the given pixel size in physical/code units (dominates over: res, lmax) : pxsize=[physical size (Number), physical unit (Symbol)]rescreate maps with the given pixel number for each deminsion; if res not given by user -> lmax is selected; (pixel number is related to the full boxsize)lmax: create maps with 2^lmax pixels for each dimensionxrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible lengthrange_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..weighting: select between:massweighting (default) and:volumeweightingdata_center: to calculate the data relative to the datacenter; in units given by argument `datacenterunit`; by default the argument datacenter = center ;data_center_unit: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)direction: select between: :x, :y, :zmask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the database (rows)ref_time: the age quantity relative to a given time (code_units); default relative to the loaded snapshot timeshow_progress: print progress bar on screenmyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of lmax, xrange, yrange, zrange, center, rangeunit, verbose, showprogress
Defined Methods - function defined for different arguments
- projection( dataobject::PartDataType, var::Symbol; ...) # one given variable
- projection( dataobject::PartDataType, var::Symbol, unit::Symbol; ...) # one given variable with its unit
- projection( dataobject::PartDataType, vars::Array{Symbol,1}; ...) # several given variables -> array needed
- projection( dataobject::PartDataType, vars::Array{Symbol,1}, units::Array{Symbol,1}; ...) # several given variables and their corresponding units -> both arrays
- projection( dataobject::PartDataType, vars::Array{Symbol,1}, unit::Symbol; ...) # several given variables that have the same unit -> array for the variables and a single Symbol for the unit
Examples
...
Mera.savedata — MethodSave loaded simulation data into a compressed/uncompressed JLD2 format:
- write new file; add datatype to existing file
- running number is taken from original RAMSES folders
- use different compression methods
- add a string to describe the simulation
- toggle verbose mode
function savedata( dataobject::DataSetType;
path::String="./",
fname = "output_",
fmode::Any=nothing,
dataformat::Symbol=:JLD2,
compress::Any=nothing,
comments::Any=nothing,
merafile_version::Float64=1.,
verbose::Bool=true)
return
Arguments
Required:
dataobject: needs to be of type: "HydroDataType", "PartDataType", "GravDataType", "ClumpDataType"fmode: nothing is written/appended by default to avoid overwriting files by accident. Need: fmode=:write (new file or overwriting existing file); fmode=:append further datatype. (overwriting of existing datatypes is not possible)
Predefined/Optional Keywords:
path: path to save the file; default is local path.fname: default name of the files "output_" and the running number is added. Change the string to apply a user-defined name.dataformat: currently, only JLD2 can be selected.compress: by default compression is activated. compress=false (deactivate).
If necessary, choose between different compression types: LZ4FrameCompressor() (default), Bzip2Compressor(), ZlibCompressor(). Load the required package to choose the compression type and to see their parameters: CodecZlib, CodecBzip2 or CodecLz4
comments: add a string that includes e.g. a description about your simulationmerafile_version: default: 1.; current only versionverbose: print timestamp and further information on screen; default: true
Defined Methods - function defined for different arguments
- savedata( dataobject::DataSetType; ...) # note: fmode needs to be given for action!
- savedata( dataobject::DataSetType, fmode::Symbol; ...)
- savedata( dataobject::DataSetType, path::String; ...)
- savedata( dataobject::DataSetType, path::String, fmode::Symbol; ...)
Mera.shellregion — FunctionCutout sub-regions of the data base of DataSetType
- select shape of a shell-region
- select size of a region (with or w/o intersecting cells)
- give the spatial center (with units) of the data relative to the full box
- relate the coordinates to a direction (x,y,z)
- inverse the selected region
- pass a struct with arguments (myargs)
shellregion(dataobject::DataSetType, shape::Symbol=:cylinder;
radius::Array{<:Real,1}=[0.,0.], # cylinder, sphere;
height::Real=0., # cylinder
direction::Symbol=:z, # cylinder
center::Array{<:Any,1}=[0., 0., 0.], # all
range_unit::Symbol=:standard, # all
cell::Bool=true, # hydro and gravity
inverse::Bool=false, # all
verbose::Bool=true, # all
myargs::ArgumentsType=ArgumentsType() ) # allArguments
Required:
dataobject: needs to be of type: "DataSetType"shape: select between region shapes: :cylinder/:disc, :sphere
Predefined/Optional Keywords:
For cylindrical shell-region, related to a given center:
radius: the inner and outer radius of the shell in units given by argumentrange_unitand relative to the givencenterheight: the hight above and below a plane [-height, height] in units given by argumentrange_unitand relative to the givencenterdirection: todo
For spherical shell-region, related to a given center:
radius: the inner and outer radius of the shell in units given by argumentrange_unitand relative to the givencenter
Keywords related to all region shapes
range_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..inverse: inverse the region selection = get the data outside of the regioncell: take intersecting cells of the region boarder into account (true) or only the cells-centers within the selected region (false)verbose: print timestamp, selected vars and ranges on screen; default: truemyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of radius, height, direction, center, range_unit, verbose
Mera.storageoverview — MethodPrint overview of the used storage per file type for a given timestep
function storageoverview(dataobject::InfoType, verbose::Bool=true)
function storageoverview(dataobject::InfoType; verbose::Bool=true)
return dictionary in bytesMera.subregion — FunctionCutout sub-regions of the data base of DataSetType
- select shape of a region
- select size of a region (with or w/o intersecting cells)
- give the spatial center (with units) of the data relative to the full box
- relate the coordinates to a direction (x,y,z)
- inverse the selected region
- pass a struct with arguments (myargs)
subregion(dataobject::DataSetType, shape::Symbol=:cuboid;
xrange::Array{<:Any,1}=[missing, missing], # cuboid
yrange::Array{<:Any,1}=[missing, missing], # cuboid
zrange::Array{<:Any,1}=[missing, missing], # cuboid
radius::Real=0., # cylinder, sphere
height::Real=0., # cylinder
direction::Symbol=:z, # cylinder
center::Array{<:Any,1}=[0.,0.,0.], # all
range_unit::Symbol=:standard, # all
cell::Bool=true, # hydro and gravity
inverse::Bool=false, # all
verbose::Bool=true, # all
myargs::ArgumentsType=ArgumentsType() ) # allArguments
Required:
dataobject: needs to be of type: "DataSetType"shape: select between region shapes: :cuboid, :cylinder/:disc, :sphere
Predefined/Optional Keywords:
For cuboid region, related to a given center:
xrange: the range between [xmin, xmax] in units given by argumentrange_unitand relative to the givencenter; zero length for xmin=xmax=0. is converted to maximum possible lengthyrange: the range between [ymin, ymax] in units given by argumentrange_unitand relative to the givencenter; zero length for ymin=ymax=0. is converted to maximum possible lengthzrange: the range between [zmin, zmax] in units given by argumentrange_unitand relative to the givencenter; zero length for zmin=zmax=0. is converted to maximum possible length
For cylindrical region, related to a given center:
radius: the radius between [0., radius] in units given by argumentrange_unitand relative to the givencenterheight: the hight above and below a plane [-height, height] in units given by argumentrange_unitand relative to the givencenterdirection: todo
For spherical region, related to a given center:
radius: the radius between [0., radius] in units given by argumentrange_unitand relative to the givencenter
Keywords related to all region shapes
range_unit: the units of the given ranges: :standard (code units), :Mpc, :kpc, :pc, :mpc, :ly, :au , :km, :cm (of typye Symbol) ..etc. ; see for defined length-scales viewfields(info.scale)center: in units given by argumentrange_unit; by default [0., 0., 0.]; the box-center can be selected by e.g. [:bc], [:boxcenter], [value, :bc, :bc], etc..inverse: inverse the region selection = get the data outside of the regioncell: take intersecting cells of the region boarder into account (true) or only the cells-centers within the selected region (false)verbose: print timestamp, selected vars and ranges on screen; default: truemyargs: pass a struct of ArgumentsType to pass several arguments at once and to overwrite default values of xrange, yrange, zrange, radius, height, direction, center, range_unit, verbose
Mera.timerfile — MethodGet a printout of the timerfile:
timerfile(object::InfoType)Mera.usedmemory — FunctionGet the memory that is used for an object in human-readable units
function usedmemory(object, verbose::Bool=true)
return value, unitMera.viewallfields — MethodGet a detailed overview of many fields from the MERA InfoType:
viewallfields(dataobject::InfoType)Mera.viewdata — MethodGet overview of stored datatypes:
- compression
- versions of the used/loaded compression
- MERA/MERA-file version
- compressed/uncompressed data size
- returns stored conversion statistics, when available (created by convertdata-function)
function viewdata(output::Int;
path::String="./",
fname = "output_",
showfull::Bool=false,
verbose::Bool=true)
return overview (dictionary)Arguments
Required:
output: output numberdatatype: :hydro, :particles, :gravity or :clumps
Predefined/Optional Keywords:
path: the path to the output JLD2 file relative to the current folder or absolute pathfname: "output"-> filename = "output***.jld2" by default, can be changed to "myname***.jld2"showfull: shows the full data tree of the datafileverbose:: informations are printed on the screen by default
Mera.viewfields — MethodGet an overview of the fields from MERA composite types:
viewfields(object)Mera.viewmodule — MethodGet a list of all exported Mera types and functions:
function viewmodule(modulename::Module)Mera.wstat — MethodCalculate statistical values w/o weighting of any Array:
wstat(array::Array{<:Real,1}; weight::Array{<:Real,1}=[1.], mask::MaskType=[false])
WStatType(mean, median, std, skewness, kurtosis, min, max)Arguments
Required:
array: Array needs to be of type: "<:Real"
Optional Keywords:
weight: Array needs to be of type: "<:Real" (can be used w/o keyword)mask: needs to be of type MaskType which is a supertype of Array{Bool,1} or BitArray{1} with the length of the Array