galaximview.snap module¶
This module contains:
a class Snap that loads a snapshot of a simulation and has methods to perform various computations on data from the snapshot.
a class IterOnSim for a simulation from which snapshots can be obtained.
- class galaximview.snap.IterOnSim(simu_directory, snapshot_basis='snapshot_', input_type='G1/G2', gadget2_header=True, ntypes=6, multiple_snaps=False, centre_on_com=False, store_pos_vel=True)[source]¶
Bases:
galaximview.units.SimCaracs
Class for simulation from which snapshots can be obtained.
- Parameters
simu_directory (str) – the directory containing the snapshots files.
snapshot_basis (str) – the basis of the snapshots names, to which formated numbers will be added, see
get_snap()
input_type ({'G1/G2', 'G2_str', 'G4_hdf5'}) – Snapshot data type, ‘G1/G2’ for usual Gadget-2 format.
gadget2_header – Only for Gadget1/Gadget2 snapshots types used in Gadget4. See
read_gadget()
ntypes – Only for Gadget1/Gadget2 snapshots types used in Gadget4. See
read_gadget()
centre_on_com (bool, default = False) – True to centre positions and velocities on the centre of mass
store_pos_vel (bool, default = True) – True to store positions and velocities before any rotation, useful for viewer
- check_snaps_present()[source]¶
Checks that snapshots are numbered as [min_nb, min_nb + 1, min_nb + 2, …, min_nb + N], if not, prints that they are not numbered this way. Returns the number of the first snapshot min_nb and the last one.
- Returns
Numbers of the first snapshot and of the last one
- Return type
int, int
- class galaximview.snap.Snap(snapshot_directory, snapshot_file_name, input_type='G1/G2', gadget2_header=True, ntypes=6, multiple_snaps=False, in_sim=True, centre_on_com=False, store_pos_vel=True)[source]¶
Bases:
galaximview.units.SimCaracs
This class loads a snapshot of a simulation and has methods to perform various computations on data from the snapshot.
- Parameters
snapshot_directory (str) – the directory containing the snapshot file
snapshot_file_name (str) – the name of the snapshot file
input_type ({'G1/G2', 'G2_str', 'G4_hdf5'}) – Snapshot data type, ‘G1/G2’ for usual Gadget-2 format.
in_sim (bool, default = True) – True if units are specified in snapshot directory in the case of G2 snapshots
gadget2_header – Only for Gadget1/Gadget2 snapshots types used in Gadget4. See
read_gadget()
ntypes – Only for Gadget1/Gadget2 snapshots types used in Gadget4. See
read_gadget()
centre_on_com (bool, default = False) – True to centre positions and velocities on the centre of mass
store_pos_vel (bool, default = True) – True to store positions and velocities before any rotation, useful for viewer
- data¶
Dictionary read from a functon of
readgadget
, according to ‘input_type’.- Type
dictionary
- listcomps¶
List of components that can be used in various methods calls.
- Type
list of str
- listindivcomps¶
List of components that can be used in various methods calls and that correspond to only one particle type.
- Type
list of str
- listquantities¶
List of quantities that can be used in various methods calls
- Type
list of str
- posini¶
(N,3) Positions as initially read, minus the centre of mass if centre_on_com is set to True.
- Type
ndarray
- velini¶
(N,3) Velocities as initially read, minus the velocity of the centre of mass if centre_on_com is set to True.
- Type
ndarray
- Raises
ValueError – If input type is unknown.
Examples
>>> from galaximview.snap import * >>> snap = Snap('path/to/simulation', 'snapshot_name') >>> hm2d, edy, edx = snap.surface_density_2d('total','xpos','ypos')
Notes
Methods often have keywords for units conversion. Possibilities are ‘cgs’, ‘astro’ and None:
‘cgs’ units are usual cgs units.
‘astro’ units are: lengths in kpc, velocities in km/s, time in Gyr, masses in sun mass, mass surface density in sun mass/pc^2, mass volume density in sun mass/pc^3. See
get_conversion_factor_from_string()
for computation of conversion factors.None for no conversion.
- binned_dispersion_of_quantity_in_2d(compstr, xstr, ystr, tobinstr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, weight_type='mass', quantity_conv_type='astro', cic_hist=False, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes and plot 2D weighted dispersion of array for component or particles from indexes or identifiers.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Member of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Member of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.tobinstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
quantity_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
cic_hist (bool, default = False) – If True, computes cloud in cloud histogram.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) 2D weighted dispersion histogram of array, (nby+1,) y axis bins edges and (nbx+1,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- binned_quantity_in_1d(compstr, xstr, ystr, inds=None, ids=None, funcofx=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, weight_type='mass', show_plot=False, funcofyforpl=<function identity_function>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Returns bins middle-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function of array. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xrange (optional) – See
HistOfX
nbins (optional) – See
HistOfX
edges (optional) – See
HistOfX
fill_range (optional) – See
HistOfX
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
show_plot (bool, default = ` False`) – True to show plot.
funcofyforpl (ufunc) – Function of weighted histogram to plot. Defaults to no transformation.
xlimit_pl – See
plot_1d()
ylimit_pl – See
plot_1d()
linestyle – See
plot_1d()
linewidth – See
plot_1d()
col – See
plot_1d()
leglab – See
plot_1d()
- Returns
(N,) Weighted histogram of y quantity at x and (N+1,) bins edges.
- Return type
ndarray, ndarray
- binned_quantity_in_2d(compstr, xstr, ystr, tobinstr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, weight_type='mass', quantity_conv_type='astro', remove_average=False, cic_hist=False, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes and plot 2D weighted histogram of array for component or particles from indexes or identifiers.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.tobinstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
quantity_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
remove_average (bool, default = False) – If True, removes mass-weighted average of histogram.
cic_hist (bool, default = False) – If True, computes cloud in cloud histogram.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Weighted 2D histogram of array, (nby+1,) y axis bins edges and (nbx+1,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- binned_quantity_in_2d_with_marginals(compstr, xstr, ystr, tobinstr, inds=None, ids=None, funcofx=<ufunc 'log10'>, funcofy=<ufunc 'log10'>, x_conv_type='astro', y_conv_type='astro', mass_conv_type=None, xlimit=(None, None), ylimit=(None, None), nbx=128, nby=128, weight_type='mass', funcforpl=<ufunc 'log10'>, cmap=None, fig=None, ax=None)[source]¶
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.tobinstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default = np.log10) – Function used to create bins on specified dimension.
funcofy (ufunc, default = np.log10) – Function used to create bins on specified dimension.
x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
mass_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
funcforpl – See
plot_hist2d_with_pdfs()
cmap – See
plot_hist2d_with_pdfs()
fig – See
plot_hist2d_with_pdfs()
ax – See
plot_hist2d_with_pdfs()
- Returns
(nby, nbx) Weighted 2D histogram, (nbx,) marginal x histogram, (nby,) marginal y histogram, (nby+1) y axis bins edges, (nbx+1) x axis bins edges and matplolib axes.
- Return type
ndarray, ndarray, ndarray, ndarray, ndarray, matplotlib axis, matplotlib axis, matplotlib axis
- change_viewing_angle(azim, elev)[source]¶
Projects positions and velocities on the right basis to obtain projection. This actually modifies self.data[‘pos’] and self.data[‘vel’]. Untouched positions and velocities are stored in self.posini and self.velini.
- Parameters
azim (float) – See ~basefuncs.get_basis_los
elev (float) – See ~basefuncs.get_basis_los
- entropy(inds)[source]¶
Returns the specific entropy, in internal units.
- Parameters
inds (ndarray) – (Ngas,) indexes to select
- Returns
(Ngas,) specific entropy, in internal units.
- Return type
ndarray
- get_ang_mom(inds)[source]¶
Returns angular momentum of the particles having the given indexes.
- Parameters
inds (ndarray) – (N,) indexes to select
- Returns
(3,) Angular momentum.
- Return type
ndarray
- get_array_from_string(arrstr, inds, conv_type=None)[source]¶
Associates array of data to the given string belonging to self.listquantities, for the given indexes. Not defined as an attribute dictionary because some arrays are not created by default (e.g. cylindrical radius).
- Parameters
arrstr (str) – Member of
listquantities
.inds (ndarray) – (N,) indexes to select
conv_type ({None, 'astro', 'cgs'}) – Type of units conversion.
- Returns
Array obtained from string identification.
- Return type
ndarray
- get_com(inds)[source]¶
Returns the centre of mass of the particles having the given indexes.
- Parameters
inds (ndarray) – (N,) indexes to select
- Returns
(3,) Position of centre of mass.
- Return type
ndarray
- get_dic_of_components_masses(conv_type='astro')[source]¶
Returns masses of components in code units and also in sun masses if self.insim=True.
- Parameters
conv_type ({'astro', 'cgs', None}) – Type of units conversion.
- Returns
Dictionary of total masses of components in self.listcomps.
- Return type
dict
- get_ids_from_indexes(inds)[source]¶
Returns the ids (Gadget-2 identifiers integers) of the given indexes.
- Parameters
inds (ndarray) – (N,) array of indexes.
- Returns
(N,) array of particle identifiers corresponding to indexes.
- Return type
ndarray
- get_indexes_from_ids(ids)[source]¶
Returns the indexes corresponding to the given ids (Gdget-2 integer identifiers). Does not keep the order of input ids ndarray.
- Parameters
ids (ndarray) – (N,) Particle identifiers.
- Returns
(N,) Indexes corresponding to identifiers. Does not keep the order of input ids ndarray.
- Return type
ndarray
- get_indexes_from_user_input(compstr, inds, ids)[source]¶
If inds and ids are None, returns indexes of component type compstr of self.listcomps. If ids is not None and inds is None, returns indexes corresponding to particle identifiers ids. If inds is not None and ids is None, returns input inds.
- Parameters
compstr (str) – Member of
listcomps
used to select indexes if inds and ids are both None.inds (ndarray or None) – (N,) Indexes returned as such if not None and if ids is not None.
ids (ndarray or None) – (N,) Identifiers for which indexes are selected if inds is None and ids is not None.
- Returns
(N,) Selected indexes.
- Return type
ndarray
- Raises
ValueError – if inds and ids are both not None.
- get_indexes_of_component(compstr)[source]¶
Selects indexes of component type of self.listcomps.
- Parameters
compstr (str) – Member of
listcomps
.- Returns
(N,) Indexes of the components selected by compstr.
- Return type
ndarray
- Raises
ValueError – If comsptr is not in self.listcomps.
- get_isnap()[source]¶
Returns the number of the snapshot, taken from numbers after last underscore of the file name.
- Returns
Number of the snapshot.
- Return type
int
- get_spin(inds)[source]¶
Returns spin (angular momentum of the particles having the given indexes, normalized to unity).
- Parameters
inds (ndarray) – (N,) indexes to select
- Returns
(3,) Angular momentum.
- Return type
ndarray
- get_temp_in_kelvin()[source]¶
Computes temperatrue in Kelvin, either from electron abundance, if present, or from some approximation.
- get_vel_of_com(inds)[source]¶
Returns the velocity of the center of mass of the particles having the given indexes.
- Parameters
inds (ndarray) – (N,) indexes to select
- Returns
(3,) Velocity of centre of mass.
- Return type
ndarray
- info()[source]¶
Prints the data dictionnary, components of the snapshot and components that can be used in methods.
- mass_in_radius(compstr, radstr, inds=None, ids=None, funcofx=<function identity_function>, rad_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, mass_conv_type='astro', show_plot=False, funcofyforpl=<function identity_function>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Mass inside radius (to be chosen between cylindrical or spherical). Returns bins middle-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.radstr (str) – radcyl’ for cylindrucal radius or ‘radspher’ for spherical radius.
inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function of radius used for x axis. Defaults to identity_function (no transfoormation).rad_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
mass_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xrange (optional) – See
HistOfX
nbins (optional) – See
HistOfX
edges (optional) – See
HistOfX
fill_range (optional) – See
HistOfX
show_plot (bool, default = ` False`) – True to show plot.
funcofyforpl (ufunc, default =
identity_function()
) – Function of mass in radius to plot. Defaults to no transformation.xlimit_pl (optional) – See
plot_1d()
ylimit_pl (optional) – See
plot_1d()
linestyle (optional) – See
plot_1d()
linewidth (optional) – See
plot_1d()
col (optional) – See
plot_1d()
leglab (optional) – See
plot_1d()
- Returns
(N,) Mass inside radius histogram and (N+1,) radius bins edges.
- Return type
ndarray, ndarray
- nhnucl(inds)[source]¶
Returns the number of H nuclei per unit of volume, in internal units.
- Parameters
inds (ndarray) – (Ngas,) indexes to select
- Returns
(Ngas,) number of H nuclei per unit of volume, in internal units.
- Return type
ndarray
- plot_3d_for_subset(compstr, frac, xstr='xpos', ystr='ypos', zstr='zpos', xlimit=(- 40, 40), ylimit=None, zlimit=None, x_conv_type='astro', y_conv_type='astro', z_conv_type='astro', inds=None, ids=None, seed=None, basis=None, zero_point=array([0, 0, 0]), ax=None, col=None, leglab=None)[source]¶
Plots a subset of particles of compstr (and returns the indexes of the subset). A seed can be provided for identical reproduction.
- Parameters
compstr (str) – Member of
listcomps
.frac (float) – Fraction of particles plotted.
xlimit ((float, float)) – Limits of x axis in units chosen by x_conv_type.
ylimit ((float, float) or None, default = None) – Limits of y axis in units chosen by y_conv_type. If None, set as same as xlimit.
zlimit ((float, float) or None, default = None) – Limits of z axis in units chosen by z_conv_type. If None, set as same as xlimit.
xstr (str, default = 'xpos') – Member of
listquantities
.ystr (str, default = 'ypos') – Member of
listquantities
.zstr (str, default = 'zpos') – Member of
listquantities
.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
z_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
seed (int, optional) – Seed for random generator, can be set to specific natural integer for reproduction.
ax (matplotlib axis, optional) – Specific matplotlib axis to plot on.
basis (ndarray or None, optional) – (3,3) basis to project positions on.
zero_point (ndarray, optional) – (3,) Removed from x, y, and z values for centering, after potential projection on basis. Defaults to np.array([0, 0, 0]).
col – See
plot_3d()
leglab – See
plot_3d()
- Returns
(N,) Indexes of the subset.
- Return type
ndarray
Examples
>>> from galaximview.snap import * >>> snap = Snap('path/to/simulation', 'snapshot_name') >>> indexes = snap.plot_3d_for_subset('all stars',0.1,xstr='xpos',ystr='ypos',zstr='zpos',xlimit=(-100,100))
- pressure(inds)[source]¶
Returns the pressure, in internal units.
- Parameters
inds (ndarray) – (Ngas,) indexes to select
- Returns
(Ngas,) pressure, in internal units.
- Return type
ndarray
- quiver(compstr, xstr, ystr, vxstr, vystr, inds=None, ids=None, x_conv_type='astro', y_conv_type='astro', xlimit=(- 40.0, 40.0), ylimit=None, nbx=128, nby=None, quantity_conv_type='astro', remove_average=False, color='k')[source]¶
Plots quiver plot (arrows), typically of 2D velocities with vxstr=’xvel’, vystr=’yvel’, xstr=’xpos’, ystr=’ypos’.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.vxstr (str) – Member of
listquantities
. Typically ‘vxpos’ if xstr is ‘xpos’.vystr (str) – Member of
listquantities
. Typically ‘vypos’ if ystr is ‘ypos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
quantity_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
remove_average (bool, default = False) – If True, removes mass-weighted average of histogram.
color (matplotlib color, default = 'k') – Matplotlib color of arrows. Defaults to black.
- reset_pos_vel()[source]¶
Can be used to reset positions and velocities of data dictionnary as read in snapshot, minus position and velocity of centre of masse if centre_on_com is True.
- Raises
ValueError – if store_pos_vel is False
- select_gas_indexes_with_temp_nh(ulimits, nhlimits, conv_type=None)[source]¶
Select gas indexes that have internal energy and number density in given limits.
- Parameters
ulimits ((float, float)) – limits of internal energy
nhlimits ((float, float)) – limits of number density
conv_type ({None, 'astro', 'cgs'}) – Type of units conversion.
- Returns
(N,) Gas particles indexes that have internal energy in ulimits and number density in nhlimits.
- Return type
ndarray
- select_gas_indexes_with_u_nh(ulimits, nhlimits, conv_type=None)[source]¶
Select gas indexes that have internal energy and number density in given limits.
- Parameters
ulimits ((float, float)) – limits of internal energy
nhlimits ((float, float)) – limits of number density
conv_type ({None, 'astro', 'cgs'}) – Type of units conversion.
- Returns
(N,) Gas particles indexes that have internal energy in ulimits and number density in nhlimits.
- Return type
ndarray
- select_new_stars_indexes_with_tbirth(tbirthlimits, conv_type='astro')[source]¶
Select indexes of new stars that have a time of birth between limits.
- Parameters
tbirthlimits ((float, float)) – limits of birth time
conv_type ({'astro', 'cgs', None}) – Type of units conversion.
- Returns
(N,) Indexes of new stars that have a time of birth in tbirthlimits.
- Return type
ndarray
- sort_data_by_id()[source]¶
Sorts all data dictionnary by ids (particle identifiers and adds a data[‘g2_type’] equal to the particle type.
- sph_binned_dispersion_of_quantity_in_2d(compstr, xstr, ystr, tobinstr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, quantity_conv_type='astro', nb_row_parts=4, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes SPH density.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.tobinstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
quantity_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
nb_row_parts (int) – Number of spatial rows parts used for speeding up of SPH treatment. Defaults to 4.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Surface density histogram, (nby+1,) y axis bins edges and (nbx,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- sph_binned_quantity_in_2d(compstr, xstr, ystr, tobinstr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, quantity_conv_type='astro', nb_row_parts=4, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes SPH density.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.tobinstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
quantity_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
nb_row_parts (int) – Number of spatial rows parts used for speeding up of SPH treatment. Defaults to 4.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Surface density histogram, (nby+1,) y axis bins edges and (nbx,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- sph_density_2d(compstr, xstr, ystr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, surfdens_conv_type='astro', nb_row_parts=4, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes SPH density.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
surfdens_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
nb_row_parts (int) – Number of spatial rows parts used for speeding up of SPH treatment. Defaults to 4.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Surface density histogram, (nby+1,) y axis bins edges and (nbx,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- stellar_formation(nb=128, tmin=0, tmax=None, show_plot=False)[source]¶
Returns formation time and histogram of amount of stellar particles.
- Parameters
nb (int, default = 128) – Number of time bins.
tmin (float, default = 0) – Minimum time.
tmax (float or None, default= None) – Maximum time. If None, set to snapshot time.
show_plot (bool, default = ` False`) – Shows plot if True.
- Returns
(nb,) time and SFR arrays, in Gyr and solar mass per yr, respectively.
- Return type
(ndarray, ndarray)
- surface_density_2d(compstr, xstr, ystr, inds=None, ids=None, funcofx=<function identity_function>, funcofy=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, volstrx='width', volstry='width', surfdens_conv_type='astro', cic_hist=False, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes and plots 2D surface density (mass weighted hist divided by surface of bins) of component or particles from indexes.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
surfdens_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
volstrx (str, optional) – Type of cell lengths used for 2D histogram. Set to ‘width’ (regular cells) if cic_hist=True. Defaults to ‘width’ (regular cells).
volstry (str, optional) – Type of cell lengths used for 2D histogram. Set to ‘width’ (regular cells) if cic_hist=True. Defaults to ‘width’ (regular cells).
cic_hist (bool, default = False) – If True, computes cloud in cloud histogram.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Surface density histogram, (nby+1,) y axis bins edges and (nbx,) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- surface_density_of_radius(compstr, inds=None, ids=None, funcofx=<function identity_function>, rad_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, surfdens_conv_type='astro', show_plot=False, funcofyforpl=<ufunc 'log10'>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Surface density as a function of cylindrical radius. Returns bins edges, middle-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function of radius used for x axis. Defaults to no transformation.rad_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
surfdens_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xrange (optional) – See
HistOfX
nbins (optional) – See
HistOfX
edges (optional) – See
HistOfX
fill_range (optional) – See
HistOfX
show_plot (bool, default = ` False`) – True to show plot.
funcofyforpl (ufunc, default = np.log10) – Function of surface density to plot.
xlimit_pl – See
plot_1d()
ylimit_pl – See
plot_1d()
linestyle – See
plot_1d()
linewidth – See
plot_1d()
col – See
plot_1d()
leglab – See
plot_1d()
- Returns
(N,) Surface density of cylindrical radius histogram and (N+1,) cylindrical radius bins edges.
- Return type
ndarray, ndarray
- velocity_dispersion_of_radius(compstr, xstr='radcyl', ystr='v_R', inds=None, ids=None, funcofx=<function identity_function>, x_conv_type='astro', y_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, show_plot=False, funcofyforpl=<ufunc 'log10'>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Default: radial velocity dispersion as a function of cylindrical radius. Returns bins middle-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function of array. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xrange (optional) – See
HistOfX
nbins (optional) – See
HistOfX
edges (optional) – See
HistOfX
fill_range (optional) – See
HistOfX
show_plot (bool, default = False) – True to show plot.
funcofyforpl (ufunc) – Function of surface density to plot. Defaults to np.log10.
xlimit_pl – See
plot_1d()
ylimit_pl – See
plot_1d()
linestyle – See
plot_1d()
linewidth – See
plot_1d()
col – See
plot_1d()
leglab – See
plot_1d()
- Returns
(N,) Velocity dispersion and (N+1,) radius bins edges.
- Return type
ndarray, ndarray
- volume_density_of_radius(compstr, inds=None, ids=None, funcofx=<function identity_function>, rad_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, voldens_conv_type='astro', show_plot=False, funcofyforpl=<ufunc 'log10'>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Volume density as a function of spherical radius. Returns bins middleR-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =) – Function of radius used for x axis. Defaults to identity_function (no transfoormation).
rad_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
voldens_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xrange (optional) – See
HistOfX
.nbins (optional) – See
HistOfX
.edges (optional) – See
HistOfX
.fill_range (optional) – See
HistOfX
.show_plot (bool, default = ` False`) – True to show plot.
funcofyforpl (ufunc, default = np.log10) – Function of volume density to plot.
xlimit_pl (optional) – See
plot_1d()
ylimit_pl (optional) – See
plot_1d()
linestyle (optional) – See
plot_1d()
linewidth (optional) – See
plot_1d()
col (optional) – See
plot_1d()
leglab (optional) – See
plot_1d()
- Returns
(N,) Volume density of spherical radius histogram and (N+1,) spherical radius bins edges.
- Return type
ndarray, ndarray
- weighted_hist_2d(compstr, xstr, ystr, inds=None, ids=None, funcofx=<ufunc 'log10'>, funcofy=<ufunc 'log10'>, x_conv_type='astro', y_conv_type='astro', mass_conv_type=None, xlimit=(-40.0, 40.0), ylimit=None, nbx=128, nby=None, weight_type='mass', cic_hist=False, show_plot=True, funcforpl=<ufunc 'log10'>, vmin_plot=None, vmax_plot=None, add_eps_in_plot=False, cmap=None, show_color_bar=True, fig=None, cb_loc=None)[source]¶
Computes and plot 2D mass-weighted hist of xstr and ystr value for component or particles from indexes. Outputs/plots the mass in each bin.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.funcofy (ufunc, default =
identity_function()
) – Function used to create bins on specified dimension. Defaults to no transformation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
mass_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
cic_hist (bool, default = False) – If True, computes cloud in cloud histogram.
show_plot (bool, default = True) – True to show plot.
funcforpl (optional) – See
plot_hist2d()
vmin_plot (optional) – See
plot_hist2d()
vmax_plot (optional) – See
plot_hist2d()
add_eps_in_plot (optional) – See
plot_hist2d()
cmap (optional) – See
plot_hist2d()
show_color_bar (optional) – See
plot_hist2d()
fig (optional) – See
plot_hist2d()
cb_loc (optional) – See
plot_hist2d()
- Returns
(nby, nbx) Weighted 2D histogram, (nby+1) y axis bins edges and (nbx+1) x axis bins edges.
- Return type
ndarray, ndarray, ndarray
- weighted_hist_2d_with_marginals(compstr, xstr, ystr, inds=None, ids=None, funcofx=<ufunc 'log10'>, funcofy=<ufunc 'log10'>, x_conv_type='astro', y_conv_type='astro', mass_conv_type=None, xlimit=(None, None), ylimit=(None, None), nbx=128, nby=128, weight_type='mass', funcforpl=<ufunc 'log10'>, cmap=None, fig=None, ax=None)[source]¶
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.ystr (str) – Members of
listquantities
. Typically ‘xpos’, ‘ypos’ or ‘zpos’.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default = np.log10) – Function used to create bins on specified dimension.
funcofy (ufunc, default = np.log10) – Function used to create bins on specified dimension.
x_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
y_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
mass_conv_type ({'astro', 'cgs', None}) – Type of units conversion for specified axis.
xlimit (optional) – See
get_hist2d()
ylimit (optional) – See
get_hist2d()
nbx (optional) – See
get_hist2d()
nby (optional) – See
get_hist2d()
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
funcforpl – See
plot_hist2d_with_pdfs()
cmap – See
plot_hist2d_with_pdfs()
fig – See
plot_hist2d_with_pdfs()
ax – See
plot_hist2d_with_pdfs()
- Returns
(nby, nbx) Weighted 2D histogram, (nbx,) marginal x histogram, (nby,) marginal y histogram, (nby+1) y axis bins edges, (nbx+1) x axis bins edges and matplolib axes.
- Return type
(ndarray, ndarray, ndarray, ndarray, ndarray, matplotlib axis, matplotlib axis, matplotlib axis)
- weighted_hist_of_array(compstr, xstr, inds=None, ids=None, funcofx=<function identity_function>, x_conv_type='astro', mass_conv_type='astro', xrange=(None, None), nbins=128, edges=None, fill_range=False, weight_type='mass', show_plot=False, funcofyforpl=<function identity_function>, xlimit_pl=None, ylimit_pl=None, linestyle='-', linewidth=1, col=None, leglab=None)[source]¶
Returns bins middle-points and histogram values.
- Parameters
compstr (str) – Member of
listcomps
.xstr (str) – Member of
listquantities
.inds (optional) – See
get_indexes_from_user_input()
ids (optional) – See
get_indexes_from_user_input()
funcofx (ufunc, default =
identity_function()
) – Function of array. Defaults to no transfoormation.x_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
mass_conv_type ({'astro', 'cgs', None}) – Type of units conversion.
xrange (optional) – See
HistOfX
nbins (optional) – See
HistOfX
edges (optional) – See
HistOfX
fill_range (optional) – See
HistOfX
weight_type ({'mass', 'number'}) – Type of weighting for 2D histogramm: ‘mass’ for mass-weighting, ‘number’ for number-weighting.
show_plot (bool, default = ` False`) – True to show plot.
funcofyforpl (ufunc, default =
identity_function()
) – Function of weighted histogram to plot. Defaults to no transformation.xlimit_pl – See
plot_1d()
ylimit_pl – See
plot_1d()
linestyle – See
plot_1d()
linewidth – See
plot_1d()
col – See
plot_1d()
leglab – See
plot_1d()
- Returns
(N,) Weighted histogram of quantity and (N+1,) bins edges.
- Return type
ndarray, ndarray
- galaximview.snap.get_hist2d(xcoords, ycoords, xrange, yrange, nbx, nby, cic_hist=False, funcofx=<function identity_function>, funcofy=<function identity_function>)[source]¶
Returns the histogram object (class Hist2DOfXY or CicHist2DOfXY). Sets nby=nbx if nby not given.
- Parameters
xcoords (ndarray) – (N,) Array to bin for the ‘x direction’ of the 2D histogram.
ycoords (ndarray) – (N,) Array to bin for the ‘y direction’ of the 2D histogram.
xrange ((float, float)) – Limits of x axis.
yrange ((float, float) or (None, None)) – Limits of y axis. If (None, None), set as equal to xlimit.
nbx (int) – Number of bins for x axis.
nby (int or None) – Number of bins for y axis. If is None, computed so that y bins have same length as x bins.
cic_hist (bool, default = False) – If True, computes cloud in cloud histogram.
funcofx (ufunc, default =
identity_function()
) – Function used for x axis. Defaults to identity_function (no transfoormation).funcofy (ufunc, default =
identity_function()
) – Function used for y axis. Defaults to identity_function (no transfoormation).
- Returns
Hist2DOfXY object.
- Return type
Hist2DOfXY
object.