Auxiliary Functions
GADGETPlotting.center_of_mass
— Methodcenter_of_mass(
position_data::Matrix{<:Real},
mass_data::Vector{<:Real},
)::Union{NTuple{3, Float64}, Nothing}
Compute the center of mass as
\[R_c = \frac{1}{M} \sum_n m_n \, r_n \, ,\]
where $M = \sum_n m_n$ and $m_n$ and $r_n$ are the mass and distance from the origin of the $n$-th particle.
If the length of $R$ is less than $10^-3$ the length of the larger position vector in position_data
, nothing is returned.
Arguments
position_data::Matrix{<:Real}
: Positions of the particles.mass_data::Vector{<:Real}
: Masses of the particles.
Returns
- The center of mass in the unis of
position_data
.
GADGETPlotting.comparison
— Methodcomparison(x, y; atol::Float64 = 1e-5, rtol::Float64 = 1e-5)::Bool
Determine if two elements are equal, as per the isequal function.
Arguments
x
: First element to be compared.y
: Second element to be compared.atol::Float64 = 1e-5
: Absolute tolerance (for compatibility).rtol::Float64 = 1e-5
: Relative tolerance (for compatibility).
Returns
- Returns
isequal(x, y)
.
GADGETPlotting.comparison
— Methodcomparison(x, y; atol::Float64 = 1e-5, rtol::Float64 = 1e-5)::Bool
Determine if two elements are equal, as per the isequal function.
Arguments
x
: First element to be compared.y
: Second element to be compared.atol::Float64 = 1e-5
: Absolute tolerance (for compatibility).rtol::Float64 = 1e-5
: Relative tolerance (for compatibility).
Returns
- Returns
isequal(x, y)
.
GADGETPlotting.comparison
— Methodcomparison(x, y; atol::Float64 = 1e-5, rtol::Float64 = 1e-5)::Bool
Determine if two elements are equal, as per the isequal function.
Arguments
x
: First element to be compared.y
: Second element to be compared.atol::Float64 = 1e-5
: Absolute tolerance (for compatibility).rtol::Float64 = 1e-5
: Relative tolerance (for compatibility).
Returns
- Returns
isequal(x, y)
.
GADGETPlotting.compute_cmdf
— Methodcompute_cmdf(
mass_data::Vector{Float64},
metallicity_data::Vector{Float64},
max_Z::Float64,
bins::Int64;
<keyword arguments>
)::NTuple{2, Vector{Float64}}
Compute the cumulative metallicity distribution function up to a metallicity of max_Z
.
The CMDF is calculated separating in bins
windows the stellar metallicity, within the range [0, max_Z
]. For the $n$-th window the CMDF is
\[\sum_{i = 1}^n \dfrac{m_i}{M_T} \quad \mathrm{vs.} \quad \bar{Z}_n \, ,\]
or, for x_norm = true
,
\[\sum_{i = 1}^n \dfrac{m_i}{M_T} \quad \mathrm{vs.} \quad \dfrac{\bar{Z}_n}{\mathrm{max}(\bar{Z}_n)} \, ,\]
where $M_T$ is the total stellar mass, $m_i$ the stellar mass of the $i$-th window and $\bar{Z}_n$ the mean stellar metallicity of the $n$-th window.
mass_data
and metallicity_data
must be in the same mass units.
Arguments
mass_data::Vector{Float64}
: Masses of the particles.metallicity_data::Vector{Float64}
: Metallicities of the particles.max_Z::Float64
: Maximum metallicity up to which the CMDF will be calculated.bins::Int64
: Number of subdivisions of [0,max_Z
] to be used for the CMDF.x_norm::Bool = false
: If the x axis will be normalized to its maximum value.
Returns
- A Tuple of two Arrays. The first with the metallicities and the second with the accumulated masses, of each window.
GADGETPlotting.deep_comparison
— Methoddeep_comparison(
x::Dict,
y::Dict;
atol::Float64 = 1e-5,
rtol::Float64 = 1e-5,
)::Bool
Determine if two dictionaries are approximately equal.
Numeric elements are compared with the comparison
function, everything else with the isequal function.
Arguments
x::Dict
: First dictionary to be compared.y::Dict
: Second dictionary to be compared.atol::Float64 = 1e-5
: Absolute tolerance for numeric elements.rtol::Float64 = 1e-5
: Relative tolerance for numeric elements.
Returns
- Return
true
if every pair of elements within the dictionaries pass the equality tests.
GADGETPlotting.deep_comparison
— Methoddeep_comparison(
x::Union{AbstractArray, Tuple},
y::Union{AbstractArray, Tuple};
atol::Float64 = 1e-5,
rtol::Float64 = 1e-5,
)::Bool
Determines if two arrays or tuples are approximately equal.
Numeric elements are compared with the comparison
function, everything else with the isequal function.
Arguments
x::Union{AbstractArray, Tuple}
: First array to be compared.y::Union{AbstractArray, Tuple}
: Second array to be compared.atol::Float64 = 1e-5
: Absolute tolerance for numeric elements.rtol::Float64 = 1e-5
: Relative tolerance for numeric elements.
Returns
- Return
true
if every pair of elements pass the equality tests.
GADGETPlotting.density_profile
— Methoddensity_profile(
mass_data::Vector{Float64},
distance_data::Vector{Float64},
max_radius::Float64,
bins::Int64,
)::NTuple{2, Vector{Float64}}
Compute a density profile up to a radius max_radius
.
It divides a sphere of radius max_radius
, centered at (0, 0, 0), in bins
spherical shells of equal width max_radius
/ bins
. This results in a volume for the $n$-th shell of
\[V_n = \frac{4}{3}\,\pi\,\mathrm{width}^3\,(3\,n^2 - 3\,n + 1) \, .\]
So, the assigned density for that shell is $\rho_n = M_n / V_n$ where $M_n$ is the total mass within the shell.
max_radius
and distance_data
must be in the same length units.
Arguments
mass_data::Vector{Float64}
: Masses of the particles.distance_data::Vector{Float64}
: Radial distances of the particles.max_radius::Float64
: Maximum distance up to which the profile will be calculated.bins::Int64
: Number of subdivisions of [0,max_radius
] to be used for the profile.
Returns
- A Tuple with two arrays. The first with the radial distances and the second with the densities, of each shell.
GADGETPlotting.energy_integrand
— Methodenergy_integrand(header::GadgetIO.SnapshotHeader, a::Float64)::Float64
Give the integrand of the scale factor to physical time function
\[\frac{1}{H\,\sqrt{\epsilon}} \, ,\]
where
\[\epsilon = \Omega_\lambda + \frac{1 - \Omega_\lambda - \Omega_0}{a^2} + \frac{\Omega_0}{a^3} \, , \]
\[H = H_0 \, a \, .\]
Arguments
header::GadgetIO.SnapshotHeader
: Header of the relevant snapshot file.a::Float64
: Dimensionless scale factor.
Returns
- The integrand in Gyr evaluated in
a
.
GADGETPlotting.format_error
— Methodformat_error(mean::Float64, error::Float64)::String
Format the mean and error values.
It follows the traditional rules for error presentation. The error has only one significant digit, unless such digit is a one, in which case, two significant digits are used. The mean will have a number of digits such as to match the last significant position of the error.
Arguments
mean::Float64
: Mean value.error::Float64
: Error value. It must be positive.
Returns
- A Tuple with the formatted mean and error values.
Examples
julia> format_error(69.42069, 0.038796)
(69.42, 0.04)
julia> format_error(69.42069, 0.018796)
(69.421, 0.019)
julia> format_error(69.42069, 0.0)
(69.42069, 0.0)
julia> format_error(69.42069, 73.4)
(0.0, 70.0)
GADGETPlotting.kennicutt_schmidt_law
— Methodkennicutt_schmidt_law(
gas_mass_data::Vector{Float64},
gas_distance_data::Vector{Float64},
temperature_data::Vector{Float64},
star_mass_data::Vector{Float64},
star_distance_data::Vector{Float64},
age_data::Vector{Float64},
temp_filter::Float64,
age_filter::Float64,
max_r::Float64;
<keyword arguments>
)::Union{Nothing, Dict{String, Any}}
Compute the mass area density and the SFR area density for the Kennicutt-Schmidt law.
The area densities are calculated by projecting the positions of the stars and of the particles of gas to the x-y plane. Then the space is subdivided in bins
concentric rings from 0 to max_r
, each of equal width max_r
/ bins
. This results in an area for the $n$-th ring of
\[A_n = π \, \mathrm{width}^2 \, (2 \, n - 1) \, .\]
So, the assigned SFR for that ring is
\[\Sigma_\mathrm{SFR}^n = \frac{M_*^n}{\mathrm{age\_filter}\,A_n} \, ,\]
where $M_*^n$ is the total mass of stars younger than age_filter
within the ring.
Equivalently, the mass area density of the gas is given by
\[\Sigma_\rho^n = \frac{M_\rho^n}{A_n} \, ,\]
where $M_\rho^n$ is the total mass of gas colder than temp_filter
within the ring.
temp_filter
and temperature_data
must be in the same temperature units, and age_filter
and age_data
must be in the same time units.
Arguments
gas_mass_data::Vector{Float64}
: Masses of the gas particles.gas_distance_data::Vector{Float64}
: 2D radial distances of the gas particles.temperature_data::Vector{Float64}
: Temperatures of the gas particles.star_mass_data::Vector{Float64}
: Masses of the stars.star_distance_data::Vector{Float64}
: 2D radial distances of the stars.age_data::Vector{Float64}
: Ages of the stars.temp_filter::Float64
: Maximum temperature allowed for the gas particles.age_filter::Float64
: Maximum age allowed for the stars.max_r::Float64
: Maximum distance up to which the parameters will be calculated.bins::Int64 = 50
: Number of subdivisions of [0,max_r
] to be used. It has to be at least 5.
Returns
- A dictionary with three entries.
- Key "RHO" => Logarithm of the area mass densities.
- Key "SFR" => Logarithm of the SFR area densities.
- Key "LM" => Linear model given by GLM.jl.
- Or
nothing
if there are less than 5 data points in the end result.
GADGETPlotting.make_video
— Methodmake_video(
source_path::String,
source_format::String,
output_path::String,
output_filename::String,
frame_rate::Int64,
)::Nothing
Make an MP4 video from a series of images.
The H.264 codec with no compression is used and the source images can be in any format available in ImageIO.jl.
Arguments
source_path::String
: Path to the directory containing the images.source_format::String
: File format of the source images, e.g. ".png", ".svg", etc.output_path::String
: Path to the directory where the resulting video will be saved.output_filename::String
: Name of the video to be generated without extension.frame_rate::Int64
: Frame rate of the video to be generated.
GADGETPlotting.mass_profile
— Methodmass_profile(
mass_data::Vector{Float64},
distance_data::Vector{Float64},
max_radius::Float64,
bins::Int64,
)::NTuple{2, Vector{Float64}}
Compute an cumulative mass profile up to a radius max_radius
.
It divides a sphere of radius max_radius
, centered at (0, 0, 0), in bins
spherical shells of equal width max_radius / bins
. This results in a cumulative mass for the $n$-th shell of
\[M_n = \sum_{i = 1}^n m_i \, ,\]
where $m_i$ is the total mass within the $i$-th shell.
max_radius
and distance_data
must be in the same length units.
Arguments
mass_data::Vector{Float64}
: Masses of the particles.distance_data::Vector{Float64}
: Radial distances of the particles.max_radius::Float64
: Maximum distance up to which the profile will be calculated.bins::Int64
: Number of subdivisions of [0,max_radius
] to be used for the profile.
Returns
- A Tuple of two arrays. The first with the radial distances and the second with the accumulated masses, of each shell.
GADGETPlotting.max_length
— Methodmax_length(data::Matrix{<:Real})::Float64
Maximum norm of the positions in data
.
data
must be a matrix with three rows (x, y, and z coordinates respectively) and where each column is a position.
Arguments
data::Matrix{<:Real}
: Positions of the particles.
Returns
- The maximum norm of the position vectors in
data
.
GADGETPlotting.metallicity_profile
— Methodmetallicity_profile(
mass_data::Vector{Float64},
distance_data::Vector{Float64},
z_data::Vector{Float64},
max_radius::Float64,
bins::Int64,
)::NTuple{2, Vector{Float64}}
Compute a metallicity profile up to a radius max_radius
(normalized to solar metallicity).
It divides a sphere of radius max_radius
, centered at (0, 0, 0), in bins
spherical shells of equal width max_radius / bins
. This results in a relative metallicity for the $n$-th shell of
\[\rho_n = \dfrac{z_n}{M_n\,Z_\odot} \, ,\]
where $M_n$ is the total mass and $z_n$ the total content of metals, within the shell.
max_radius
and distance_data
must be in the same length units, and z_data
and mass_data
must be in the same mass units.
Arguments
mass_data::Vector{Float64}
: Masses of the particles.distance_data::Vector{Float64}
: Radial distances of the particles.z_data::Vector{Float64}
: Metal content of the particles in mass units.max_radius::Float64
: Maximum distance up to which the profile will be calculated.bins::Int64
: Number of subdivisions of [0,max_radius
] to be used for the profile.
Returns
- A Tuple of two arrays. The first with the radial distances and the second with the metallicities, of each shell.
GADGETPlotting.num_integrate
— Functionnum_integrate(
func::Function,
inf_lim::Float64,
sup_lim::Float64,
steps::Int64 = 200,
)::Float64
Give the numerical integral of func
between inf_val
and sup_val
.
The result is given by
\[\int_\mathrm{inf\_lim}^\mathrm{sup\_lim} f(x) \mathrm{dx} \approx \sum_{i = 1}^\mathrm{steps} f(\mathrm{inf\_lim} + \mathrm{width}\,i ) \, ,\]
where $\mathrm{width} = (\mathrm{sup\_lim} - \mathrm{inf\_lim}) / \mathrm{steps}$.
Arguments
func::Function
: 1D function to be integrated.inf_lim::Float64
: Lower limit of the integral.sup_lim::Float64
: Upper limit if the integral.steps::Int64
: Number of subdivisions to be used for the discretization of the [inf_lim
,sup_lim
] region.
Returns
- The value of the integral.
Examples
julia> num_integrate(sin, 0, 3π)
1.9996298761360816
julia> num_integrate(x -> x^3 + 6 * x^2 + 9 * x + 2, 0, 4.69)
438.9004836958452
julia> num_integrate(x -> exp(x^x), 0, 1.0)
2.1975912134624904
GADGETPlotting.pass_all
— Methodpass_all(snap_file::String, type::String)::Vector{Int64}
Default filter function for the read_blocks_filtered function.
It does not filter out any particles, allowing the data acquisition functions to gather all the data.
Arguments
snap_file::String
: Snapshot file path.type::String
: Particle type.- "gas" -> Gas particle.
- "dark_matter" -> Dark matter particle.
- "stars" -> Star particle.
Returns
- A Vector with the indices of the allowed particles.
GADGETPlotting.quantities_2D
— Methodquantities_2D(
gas_mass_data::Vector{Float64},
gas_distance_data::Vector{Float64},
temperature_data::Vector{Float64},
star_mass_data::Vector{Float64},
star_distance_data::Vector{Float64},
age_data::Vector{Float64},
metal_mass_data::Matrix{Float64},
fmol_data::Vector{Float64},
temp_filter::Float64,
age_filter::Float64,
max_r::Float64;
bins::Int64 = 50,
)::Dict{String, Vector}
Compute the surface density of several quantities.
The area densities are calculated by projecting the positions of the stars and the particles of gas to the x-y plane. Then the space is subdivided in bins
concentric rings from 0 to max_r
, each of equal width max_r
/ bins
. This results in an area for the $n$-th ring of
\[A_n = π \, \mathrm{width}^2 \, (2 \, n - 1) \, .\]
So, the assigned SFR for that ring is
\[\Sigma_\mathrm{SFR}^n = \frac{M_*^n}{\mathrm{age\_filter}\,A_n} \, ,\]
where $M_*^n$ is the total mass of stars younger than age_filter
within the ring.
Equivalently, the mass area density of the gas and stars is given by
\[\Sigma_\rho^n = \frac{M_\rho^n}{A_n} \, ,\]
\[\Sigma_*^n = \frac{M_*^n}{A_n} \, ,\]
where $M_\rho^n$ is the total mass of gas colder than temp_filter
within the ring.
temp_filter
and temperature_data
must be in the same temperature units, and age_filter
and age_data
must be in the same time units.
The rest of the parameters are define as follows
\[\mathrm{SSFR}^n = \frac{\Sigma_\mathrm{SFR}^n}{\Sigma_{*\mathrm{(total)}}^n} \, ,\]
\[\mathrm{SFE}^n = \frac{\Sigma_\mathrm{SFR}^n}{\Sigma_{\rho\mathrm{(total)}}^n} \, ,\]
\[\mathrm{P}^n = \Sigma_{\rho\mathrm{(total)}}^n\left(\Sigma_{\rho\mathrm{(total)}}^n + \Sigma_{*\mathrm{(total)}}^n\right) \, ,\]
\[\mathrm{Ψ/H_2}^n = \frac{\Sigma_\mathrm{SFR}^n \, A_n}{M_{H_2}^n} \, ,\]
where $\Sigma_{\rho\mathrm{(total)}}^n$ is the mass density of all the gas, $\Sigma_{*\mathrm{(total)}}^n$ is the stellar density of all the stars and $M_{H_2}^n$ is the total mass of molecular Hydrogen.
Arguments
gas_mass_data::Vector{Float64}
: Masses of the gas particles.gas_distance_data::Vector{Float64}
: 2D radial distances of the gas particles.temperature_data::Vector{Float64}
: Temperatures of the gas particles.star_mass_data::Vector{Float64}
: Masses of the stars.star_distance_data::Vector{Float64}
: 2D radial distances of the stars.age_data::Vector{Float64}
: Ages of the stars.metal_mass_data::Matrix{Float64}
: Masses of the individual elements (H, O, He, C, etc.) within the gas particles.fmol_data::Vector{Float64}
: Fraction of molecular Hydrogen of the gas particles.temp_filter::Float64
: Maximum temperature allowed for the gas particles.age_filter::Unitful.Quantity
: Maximum age allowed for the stars.max_r::Float64
: Maximum distance up to which the parameters will be calculated.bins::Int64 = 50
: Number of subdivisions of [0,max_r
] to be used. It has to be at least 5.
Returns
- A dictionary with nine entries.
- Key "GAS" => Surface mass density of gas
- Key "COLD_GAS" => Surface mass density of cold gas
- Key "STARS" => Surface mass density of stars
- Key "OH" => 12 + log10(oxygenmass / hydrogenmass)
- Key "SFR" => Star formation rate surface density
- Key "SSFR" => Specific star formation rate surface density
- Key "SFE" => Star formation efficiency surface density
- Key "P" => Proportional to the pressure
- Key "Psi_FMOL" => Star formation rate per unit of molecular gas
GADGETPlotting.relative
— Functionrelative(
p::Plots.Plot,
rx::Float64,
ry::Float64,
rz::Union{Float64, Nothing} = nothing;
<keyword arguments>
)::Union{NTuple{2, Float64}, NTuple{3, Float64}}
Give the absolute coordinates within a plot, from the relative ones.
If any of the axes are in a logarithmic scale, you should set log
appropriately to get correct results.
Arguments
p::Plots.Plot
: Plot for which the absolute coordinates will be calculated.rx::Float64
: relative x coordinate,rx
∈ [0, 1].ry::Float64
: relative y coordinate,ry
∈ [0, 1].rz::Union{Float64,Nothing} = nothing
: relative z coordinate,rz
∈ [0, 1].log::Union{NTuple{2, Bool}, NTuple{3, Bool}} = (false, false, false)
: If the x, y or z axes are in a logarithmic scale.
Returns
- A Tuple with the absolute coordinates: (x, y) or (x, y, z).
Examples
julia> GADGETPlotting.relative(plot(rand(100)), 0.5, 0.5)
(50.5, 0.5047114800322484)
julia> GADGETPlotting.relative(surface(rand(100, 100)), 0.5, 0.5, 0.5)
(50.5, 50.5, 0.5000284432744244)
GADGETPlotting.set_vertical_flags
— Methodset_vertical_flags(
flags::Union{Tuple{Vector{<:Real}, Vector{<:AbstractString}}, Nothing},
plot::Plots.Plot;
<keyword arguments>
)::Plots.Plot
Draw vertical lines at specified positions.
If you only want labels for the vertical lines, the original plot should have label = ""
.
Arguments
flags::Union{Tuple{Vector{<:Real}, Vector{<:AbstractString}}, Nothing} = nothing
: The first vector in the Tuple has the positions of the vetical lines. The second has the corresponding labels.plot::Plots.Plot
: Plot to which the vertical lines will be added.
Returns
- New plot with the vertical lines added.
GADGETPlotting.smooth_window
— Methodsmooth_window(
x_data::Vector{<:Real},
y_data::Vector{<:Real},
bins::Int64,
)::NTuple{2, Vector{Float64}}
Separate the values in x_data
in bins
contiguous windows, and replaces every x and y value within the window with the corresponding mean to smooth out the data.
Arguments
x_data::Vector{<:Real}
: x-axis data.y_data::Vector{<:Real}
: y-axis data.bins::Int64
: Number of windows to be used in the smoothing.log::Bool = false
: If the x-axis data will be separated using logarithmic bins.
Returns
- A Tuple with two arrays containing the smoothed out x and y data.