Common Functionalities

Base Profiling Class

The cosmic_base_class module defines the CosmicBase class that constitutes the parent class of many profiling classes in the CosmicProfiles package.

class common.cosmic_base_class.CosmicBase(unicode SNAP, float L_BOX, int MIN_NUMBER_PTCS, unicode CENTER, unicode VIZ_DEST, unicode CAT_DEST, unicode SUFFIX)[source]

Bases: object

Parent class governing high-level cosmic shape calculations

Its public methods are _getMassesCentersBase(), _getShapeCatLocalBase(), _getShapeCatGlobalBase(), _getShapeCatVelLocalBase(), _getShapeCatVelGlobalBase(), _dumpShapeCatLocalBase(), _dumpShapeCatGlobalBase(), _dumpShapeCatVelLocalBase(), _dumpShapeCatVelGlobalBase(), _plotShapeProfsBase(), _plotLocalTHistBase(), _plotGlobalTHistBase(), _getDensProfsBestFitsBase(), _getConcentrationsBase(), _getDensProfsSphDirectBinningBase(), _getDensProfsEllDirectBinningBase(), _getDensProfsKernelBasedBase(), _getObjInfoBase()

Parameters
  • SNAP (string) – snapshot identifier, e.g. ‘024’

  • L_BOX (float) – simulation box side length in config.InUnitLength_in_cm

  • MIN_NUMBER_PTCS (int) – minimum number of particles for object to qualify for morphology calculation

  • CENTER (str) – shape quantities will be calculated with respect to CENTER = ‘mode’ (point of highest density) or ‘com’ (center of mass) of each halo

  • VIZ_DEST (string) – visualization folder

  • CAT_DEST (string) – catalogue destination

  • SUFFIX (string) – either ‘_dm_’ or ‘_gx_’ or ‘’ (latter for CosmicProfsDirect)

__init__()[source]
Parameters
  • SNAP (string) – snapshot identifier, e.g. ‘024’

  • L_BOX (float) – simulation box side length in config.InUnitLength_in_cm

  • MIN_NUMBER_PTCS (int) – minimum number of particles for object to qualify for morphology calculation

  • CENTER (str) – shape quantities will be calculated with respect to CENTER = ‘mode’ (point of highest density) or ‘com’ (center of mass) of each halo

  • VIZ_DEST (string) – visualization folder

  • CAT_DEST (string) – catalogue destination

  • SUFFIX (string) – either ‘_dm_’ or ‘_gx_’ or ‘’ (latter for CosmicProfsDirect)

estDensProfs(self, r_over_r200, obj_numbers, bool direct_binning=True, bool spherical=True, bool reduced=False, bool shell_based=False)[source]

Estimate density profiles

Parameters
  • r_over_r200 ((r_res,) floats) – normalized radii at which to-be-estimated density profiles are defined

  • obj_numbers (list of int) – list of object indices of interest

  • direct_binning (boolean) – whether or not direct binning approach or kernel-based approach should be used

  • spherical (boolean) – whether or not spherical shell-based or ellipsoidal shell-based should be used, ignored if direct_binning = False

  • reduced (boolean) – whether or not reduced shape tensor (1/r^2 factor) should be used, ignored if direct_binning = False

  • shell_based (boolean) – whether shell-based or ellipsoid-based algorithm should be run

Returns

density profiles in units of config.OutUnitMass_in_g/config.OutUnitLength_in_cm**3

Return type

(N2, r_res) floats

Cosmo Tools

The cosmo_tools module provides simple Python tools that are useful for the analysis of cosmological datasets, such as KS tests, chi-square tests etc.

common.cosmo_tools.KSTest(x, y)[source]

2-Sample KS Test of Independence for Continuous Distributions

Get KS test statistic and p-value, assuming x and y are drawn from the same distribution. Q we could ask: Is the PDF independent of mass?

Parameters
  • x ((# sample values) floats, can be both np.array and list) – 1st sample, all sample points, not binned, whether normalized or not is irrelevant

  • y ((# sample values) floats) – 2nd sample, all sample points, not binned, whether normalized or not is irrelevant

Returns

KS-stat, p-value

Return type

float, float

common.cosmo_tools.M_split(m, center, start_time, v=None, NB_BINS=2)[source]

Mass-splitting of quantities center, v

Mass units are determined by m

Parameters
  • m (list of floats) – total mass of objects

  • center (list of (3,) float arrays) – centers of objects

  • start_time (float) – time of start of shape analysis

  • v (list of (3,) float arrays) – major axis of objects (optional) or any other vectorial quantity

  • NB_BINS (float) – Number of mass bins

Returns

max_min_m (mass bin edges), m_groups (total mass in bins), center_groups (center average in each bin), v_groups (v average in each bin), idx_groups (indices of all objects in each bin)

Return type

(N,) floats, (N-1,) floats, (N-1,) floats, (N-1,) floats, list of lists (for each bin) of ints

common.cosmo_tools.calcMode(xyz, masses, rad)[source]

Find mode of point distribution xyz

Parameters
  • xyz ((N,3) floats) – coordinates of particles of type 1 or type 4

  • masses ((N,) floats) – masses of particles of type 1 or type 4

  • rad (float) – initial radius to consider away from COM of object

Returns

mode of point distribution

Return type

(3,) floats

common.cosmo_tools.chiSquareToICore(occs)[source]

Chi Square Test of Independence for Binned Data

Get chi square and p-value, assuming occs[0,:] and occs[1,:] are drawn from the same distribution. Q we could ask: Is the PDF independent of mass?

Parameters

occs ((2 = # mass groups, # mass or sep bins) floats) – spine of shaded curves, but rescaled to be cell counts (!or weighted cell counts!), not PDF

Returns

chi-square, p-value, degrees of freedom

Return type

float, float, int

common.cosmo_tools.chiSquareToIMultiDim(x, y)[source]

Calculate chi-square test for multidimensional categorical data

Both x and y are multi-dimensional (D_BINS+1, varying numbers of elements in each bin), binning into 1D-bins is done before handing over to chiSquareToICore() :param x: sample 1, multi-dimensional :type x: list of float lists :param y: sample 1, multi-dimensional :type y: list of float lists :return: chi-square, p-value, degrees of freedom :rtype: float, float, int

common.cosmo_tools.chiSquareToIShapeCurve(d_cdm, param_interest_cdm, d_fdm, param_interest_fdm, D_LOGSTART, D_LOGEND, D_BINS)[source]

Chi-square test for shape curves

Hands over to chiSquareToIMultiDim()

Parameters
  • d_cdm (float array) – array of ellipsoidal radii, sample 1

  • param_interest_cdm (float array) – shape quantitity of interest, sample 1

  • d_fdm (float array) – array of ellipsoidal radii, sample 2

  • param_interest_fdm (float array) – shape quantitity of interest, sample 2

  • D_LOGSTART (int) – logarithm of minimum ellipsoidal radius of interest, in units of R200 of parent halo

  • D_LOGEND (int) – logarithm of maximum ellipsoidal radius of interest, in units of R200 of parent halo

  • D_BINS (int) – number of ellipsoidal radii of interest minus 1 (i.e. number of bins)

Returns

chi-square, p-value, degrees of freedom

Return type

float, float, int

common.cosmo_tools.getBlocks(lst, n)[source]

Yield successive n-sized blocks from lst (list or np.array).

Parameters
  • lst (list or np.array) – array-like object to draw successive blocks from

  • n (int) – size of blocks, last block might be smaller

Returns

blocks of size n

Return type

generator function

common.cosmo_tools.getMeanOrMedianAndError(y, ERROR_METHOD, N_REAL=10)[source]

Return mean (if ERROR_METHOD == “bootstrap” or “SEM”) or median (if ERROR_METHOD == “median_quantile”) and the +- 1 sigma error attached

Parameters
  • y (1D float array) – data

  • ERROR_METHOD (string) – error method, either ‘bootstrap’, ‘SEM’ or ‘median_quantile’

  • N_REAL (int) – number of bootstraps, only relevant if ERROR_METHOD == ‘bootstrap’

Returns

mean_median, err_low, err_high

Return type

float, float, float

common.cosmo_tools.respectPBC(xyz, L_BOX)[source]

Return positions xyz that respect the box periodicity

If point distro xyz has particles separated in any Cartesian direction by more than L_BOX/2, reflect those particles along L_BOX/2

Parameters
  • xyz ((N,3) floats) – coordinates of particles of type 1 or type 4

  • L_BOX (float, units: Mpc/h) – simulation box side length

Returns

updated coordinates of particles of type 1 or type 4

Return type

(N,3) floats

Python Routines

The python_routines module helps in calculating the center of mass of a point cloud, drawing uniformly from an ellipsoid, drawing uniformly from an ellipsoidal shell, drawing Fibonacci sphere samples etc.

common.python_routines.calcCoM(xyz, masses)[source]

Calculate center of mass of point distribution

Parameters
  • xyz ((N,3) floats) – coordinates of particles of type 1 or type 4

  • masses ((N,3) floats) – masses of the particles

Returns

com, center of mass

Return type

(3,) floats

common.python_routines.calcMode(xyz, masses, rad)[source]

Find mode (point of highest local density) of point distribution xyz

Parameters
  • xyz ((N,3) floats) – coordinates of particles of type 1 or type 4

  • masses ((N,) floats) – masses of the particles

  • rad (float) – initial radius to consider from CoM of object

Returns

mode of distro

Return type

(3,) floats

common.python_routines.checkDensFitMethod(method)[source]

Check validity of density profile fitting method

Parameters

method (dictionary, method[‘profile’] is either einasto, alpha_beta_gamma, hernquist, nfw, minimum requirement) – describes density profile model assumed for fitting, if parameter should be kept fixed during fitting then it needs to be provided, e.g. method[‘alpha’] = 0.18

common.python_routines.checkKatzConfig(katz_config)[source]

Check (for types etc) and return configuration parameters for Katz algorithm

Parameters

katz_config (dictionary) – dictionary with parameters to the Katz algorithm, with fields ‘ROverR200’, ‘IT_TOL’, ‘IT_WALL’, ‘IT_MIN’, ‘REDUCED’, ‘SHELL_BASED’

Return ROverR200, IT_TOL, IT_WALL, IT_MIN, REDUCED, SHELL_BASED

configuration parameters

Return type

(r_res,) doubles, double, int, int, boolean, boolean

common.python_routines.drawUniformFromEllipsoid(dims, a, b, c, Nptc, ell_nb)[source]

Draw points uniformly from an ellipsoid volume, centered at the origin

This function is primarily used to generate synthetic halos. The approach is taken from Rubinstein & Kroese 2007:

Part 1: Generating uniform random vectors inside the 3-ball

Part 2: Lower Cholesky decomposition

Part 3: Draw uniformly from ellipsoid

Parameters
  • dims (int) – number of dimensions of ‘ellipsoid’, either 2, in which case c remains ‘None’, or 3

  • a (float array, units are Mpc/h) – major axis array

  • b (float array, units are Mpc/h) – intermediate axis array

  • c (float array, units are Mpc/h) – minor axis array

  • Nptc (int array) – number of particles in each ellipsoid Ell(a[idx], b[idx], c[idx]).

  • ell_nb (int) – ellipsoid number of interest

Returns

points drawn uniformly from an ellipsoid volume, ellipsoid number

Return type

(N_reals,3) floats, 1 int

common.python_routines.drawUniformFromShell(dims, a, b, c, Nptc, shell_nb)[source]

Draw points uniformly from an ellipsoidal shell volume, centered at the origin

This function is primarily used to generate synthetic halos. The approach is taken from Rubinstein & Kroese 2007:

Part 1: Generating uniform random vectors inside the 3-ball

Part 2: Lower Cholesky decomposition

Part 3: Draw uniformly from Ball(a[shell_nb]) Part 4: Move ptcs from ball Ball(a[shell_nb]) into ellipsoidal shell Shell(a[shell_nb-1],b[shell_nb-1],c[shell_nb-1],a[shell_nb],b[shell_nb],c[shell_nb])

Parameters
  • dims (int) – number of dimensions of ‘ellipsoid’, either 2, in which case c remains ‘None’, or 3

  • a (float array, units are Mpc/h) – major axis array

  • b (float array, units are Mpc/h) – intermediate axis array

  • c (float array, units are Mpc/h) – minor axis array

  • Nptc (int array) – number of particles in each shell Shell(a[shell_nb-1],b[shell_nb-1],c[shell_nb-1],a[shell_nb],b[shell_nb],c[shell_nb]).

  • shell_nb (int) – shell number of interest

Returns

points drawn uniformly from an ellipsoidal shell volume, shell number

Return type

(N_reals,3) floats, 1 int

common.python_routines.eTo10(st)[source]

Replace e+{xy} by “10^{xy}” etc..

Note that “st” should be handed over with {:.2E}. Works for exponents 00 to 19, plus and minus, e and E.

Example

>>> area = 2
>>> print("E.g. internally we write '²' for 'square'. The area of your rectangle is {} cm²".format(area))
>>> test = "2.88E-19"
>>> print("The string", test, "becomes", eTo10(test))
common.python_routines.fibonacci_ellipsoid(a, b, c, samples=1)[source]

Creating “evenly” distributed points on an ellipsoid surface

Here, this distribution of points on the ellipsoid surface is called a Fibonacci ellipsoid

Parameters
  • a (float) – major axis of ellipsoid surface

  • b (float) – intermediate axis of ellipsoid surface

  • c (float) – minor axis of ellipsoid surface

  • samples (int) – number of points to be put onto the unit sphere

Returns

N points on the unit sphere

Return type

list of N (3,) floats

common.python_routines.fibonacci_sphere(samples=1)[source]

Creating “evenly” distributed points on a unit sphere

This distribution of points on the unit sphere is called a Fibonacci sphere

Parameters

samples (int) – number of points to be put onto the unit sphere

Returns

N points on the unit sphere

Return type

list of N (3,) floats

common.python_routines.getCatWithinFracR200(cat_in, obj_size_in, xyz, masses, L_BOX, CENTER, r200, frac_r200)[source]

Cleanse index catalogue cat_in of particles beyond R200 r200

Parameters
  • cat_in ((N3) integers) – contains indices of particles belonging to an object

  • obj_size_in ((N1,) integers) – indicates how many particles are in each object

  • xyz ((N2,3) floats) – coordinates of particles of type 1 or type 4

  • masses ((N2,3) floats) – masses of the particles

  • L_BOX (float) – periodicity of the box

  • CENTER (str) – shape quantities will be calculated with respect to CENTER = ‘mode’ (point of highest density) or ‘com’ (center of mass) of each halo

  • r200 ((N1,) floats) – R_200 radii of the parent halos

  • frac_r200 (float) – depth of objects to plot ellipticity, in units of R200

Return cat_out, obj_size_out

updated cat_in and obj_size_out, with particles beyond R200 r200 removed

Return type

list of length N1

common.python_routines.getDelta(z, OMEGA_M, OMEGA_L)[source]

Calculates the overdensity required for spherical virialization relative to the mean background of the universe

The result is taken from Bryan, Norman, 1998, ApJ (Vir. Theorem + Spher. Collapse). Note that this quantity is z-dependent.

Parameters
  • z (float) – redshift of interest

  • OMEGA_M (float) – fractional matter density of the universe

  • OMEGA_L (float) – fractional dark energy density of the universe

Returns

spherical virialization overdensity

Return type

float

common.python_routines.getMassDMParticle(N, L)[source]

Returns the mass of each DM particle, assuming it is constant, given N and L

Parameters
  • N (int) – Simulation resolution parameter, usually power of 2

  • L (float) – Simulation box size in cMpc/h

Returns

Mass of each DM particle in solar masses*h^-1

Return type

float

common.python_routines.getRhoCrit()[source]

Returns critical comoving density of the universe

Return value is in units of [solar masses/(cMpc)^3*h^2]

Returns

critical comoving density

Return type

float

common.python_routines.getSubSetIdxCat(idx_cat, obj_size, obj_numbers)[source]

Get the indices from idx_cat that correspond to object numbers obj_numbers

Parameters
  • idx_cat ((N3,) integers) – contains indices of particles belonging to an object

  • obj_size ((N1,) integers) – indicates how many particles are in each object

  • obj_numbers (list of int) – list of object indices of interest

Return subset_idx_cat

indices from idx_cat that correspond to requested object numbers

Return type

(N4,) integers

common.python_routines.inShell(X, a, b, c, shell_nb)[source]

Determine whether point X is in Shell(a[shell_nb], b[shell_nb], c[shell_nb])

Parameters
  • X ((3,) float array) – point of interest

  • a (float array, units are Mpc/h) – major axis array

  • b (float array, units are Mpc/h) – intermediate axis array

  • c (float array, units are Mpc/h) – minor axis array

  • shell_nb (int) – ellipsoid number of interest

Returns

True if X is in Shell(a[shell_nb], b[shell_nb], c[shell_nb]), False otherwise

Return type

boolean

common.python_routines.isValidSelection(obj_numbers, nb_objects)[source]

Trivial function to check whether selection of objects is valid

Parameters
  • obj_numbers (list of int) – list of object indices of interest

  • nb_objects (integer) – number of objects in inventory

Return valid

True if selection of objects is valid one, False otherwise

Return type

boolean

Raises

ValueError if selection of objects is invalid

common.python_routines.print_status(rank, start_time, message, **kwargs)[source]

Routine to print script status to command line, with elapsed time

Parameters
  • rank (int) – rank of process

  • start_time (float) – time of start of script (or specific function)

  • message (string) – message to be printed

  • color (string) – for optional color argument, we allow for ‘green’, ‘blue’ and ‘red’

common.python_routines.recentreObject(xyz, L_BOX)[source]

Recentre object if fallen outside [L_BOX]^3 due to e.g. respectPBCNoRef()

Parameters
  • xyz ((N,3) floats) – coordinates of particles of type 1 or type 4

  • L_BOX (float) – periodicity of box (0.0 if non-periodic)

Returns

updated coordinates of particles

Return type

(N^3x3) floats

common.python_routines.respectPBCNoRef(xyz, L_BOX=None)[source]

Return modified positions xyz_out of an object that respect the box periodicity

If point distro xyz has particles separated in any Cartesian direction by more than L_BOX/2, translate those particles accordingly.

Parameters
  • xyz ((N^3x3) floats) – coordinates of particles

  • L_BOX (float) – periodicity of box (0.0 if non-periodic)

Returns

updated coordinates of particles

Return type

(N^3x3) floats

common.python_routines.set_axes_equal(ax: Axes)[source]

Set 3D plot axes to equal scale.

Make axes of 3D plot have equal scale so that spheres appear as spheres and cubes as cubes. Required since ax.axis('equal') and ax.set_aspect('equal') don’t work on 3D.

Caching Routines

The caching module implements the caching of function outputs to avoid recalculating/reloading repeatedly. A modified decorator has been added that in addition to the standard maxsize argument takes the optional argument use_memory_up_to. If set (i.e. use_memory_up_to is not False), the cache will be considered full if there are fewer than use_memory_up_to free bytes of memory available (according to psutil.virtual_memory().available). Note that the maxsize argument will have no effect in that case. The motivation behind a memory-aware LRU caching function is that caching too many values causes thrashing, which should be avoided if possible.

The user has the liberty to update use_memory_up_to and maxsize by calling updateCachingMaxGBs() or updateCachingMaxSize(), respectively.

common.config.LMVLabel()[source]

Return strings that can be used for plots, i.e. in out units

common.config.LengthInToInternal(length_arr)[source]

Rescale length_arr from in length units to internal length units (Mpc/h)

Parameters

length_arr (floats) – length array in in units

common.config.LengthInternalToOut(length_arr)[source]

Rescale length_arr from internal length units (Mpc/h) to out length units

Parameters

length_arr (floats) – length array in internal units

common.config.MassInToInternal(mass_arr)[source]

Rescale mass_arr from in mass units to internal mass units (E+10M_sun/h)

Parameters

mass_arr (floats) – mass array in in units

common.config.MassInternalToOut(mass_arr)[source]

Rescale mass_arr from internal mass units (E+10M_sun/h) to out mass units

Parameters

mass_arr (floats) – mass array in internal units

common.config.VelInToInternal(vel_arr)[source]

Rescale vel_arr from in velocity units to internal velocity units (km/s)

Parameters

vel_arr (floats) – velocity array in in units

common.config.VelInternalToOut(vel_arr)[source]

Rescale vel_arr from internal velocity units (km/s) to out velocity units

Parameters

vel_arr (floats) – velocity array in internal units

common.config.getAdmissibleLMV()[source]
common.config.getLMVInternal()[source]
common.config.initialize(gbs, maxsize, in_length_in_cm, in_mass_in_g, in_velocity_in_cm_per_s, little_h, out_length_in_cm, out_mass_in_g, out_velocity_in_cm_per_s)[source]

Initialize all global variables, should only be called once.

Parameters
  • gbs (int, or float) – new memory limit in GBs

  • maxsize (int) – new maximum cache size (note that multiple functions are decorated, each will have new maximum cache size)

  • in_length_in_cm (float) – unit for length in cm as used in provided data

  • in_mass_in_g (float) – unit for mass in g as used in provided data

  • in_velocity_in_cm_per_s (float) – unit for velocity in cm/s as used in provided data

  • little_h (float) – scaled Hubble constant

  • out_length_in_cm (float) – unit for length in cm to be used by CosmicProfiles

  • out_mass_in_g (float) – unit for mass in g to be used by CosmicProfiles

  • out_velocity_in_cm_per_s (float) – unit for velocity in cm/s to be used by CosmicProfiles

common.config.updateCachingMaxGBs(gbs)[source]

Update use_memory_up_to argument of LRU-caching decorator

Parameters

gbs (int, or float) – new memory limit in GBs

common.config.updateCachingMaxSize(maxsize)[source]

Update maxsize argument of LRU-caching decorator

Parameters

maxsize (int) – new maximum cache size (note that multiple functions are decorated, each will have new maximum cache size)

common.config.updateInUnitSystem(length_in_cm, mass_in_g, velocity_in_cm_per_s, little_h)[source]

Inform CosmicProfiles about system of units employed in snapshot data

Parameters
  • length_in_cm (float) – unit for length in cm as used in provided data

  • mass_in_g (float) – unit for mass in g as used in provided data

  • velocity_in_cm_per_s (float) – unit for velocity in cm/s as used in provided data

  • little_h (float) – scaled Hubble constant

common.config.updateOutUnitSystem(length_in_cm, mass_in_g, velocity_in_cm_per_s, little_h)[source]

Update unit system to be used by CosmicProfiles in its outputs

Parameters
  • length_in_cm (float) – unit for length in cm to be used by CosmicProfiles

  • mass_in_g (float) – unit for mass in g to be used by CosmicProfiles

  • velocity_in_cm_per_s (float) – unit for velocity in cm/s to be used by CosmicProfiles

  • little_h (float) – scaled Hubble constant

The default caching parameters are:

CACHE_MAXSIZE = 128
GBs = 2

To disable caching, use:

from cosmic_profiles import updateCachingMaxSize, updateCachingMaxGBs
updateCachingMaxGBs(False)
updateCachingMaxSize(None)