API Reference¶
Main Routines (projected space)¶
These are the main routines for counting pairs and estimating correlation functions. They encapsulate all required steps to provide a true end-user experience : you only need to provide data+parameters to receive final results, and perhaps that amazing plot of your next paper
-
gundam.
pcf
(tab, tab1, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given two astropy tables corresponding to data and random samples, calculate the two-point projected auto-correlation function (pcf)
All input parameters that control binning, cosmology, estimator, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (pcf) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for DD/RR/DR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.rr, counts.wrp, etc. It also stores the complete log and all input parameters. See Output Dictionary (pcf) for a detailed description.
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data rans = Table.read('redgal_rand.fits') # read randoms par = gun.packpars(kind='pcf',outfn='redgal') # generate default parameters cnt = gun.pcf(gals, rans, par, write=True, plot=True) # get pcf and plot
-
gundam.
pccf
(tab, tab1, tab2, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given three astropy tables corresponding to data, random, and cross samples, this routine calculates the two-point projected cross-correlation function (pccf)
All input parameters that control binning, cosmology, estimator, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- tab2 : astropy table
- Table with cross sample particles (C)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (pccf) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for CD/CR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.cd, counts.cr, counts.wrp, etc. It also stores the complete log and all input parameters. See Output Dictionary (pccf) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data files rans = Table.read('redgal_rand.fits') qsos = Table.read('qso.fits') par = gun.packpars(kind='pccf',outfn='redgal_qso') # generate default parameters cnt = gun.pccf(gals, rans, qsos, par, write=True ) # get pcf
-
gundam.
rppi_A
(tab, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given an astropy data table, count pairs in the projected-radial space (\(r_p\) vs \(\pi\) space).
All input parameters that control binning, cosmology, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (rppiA) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=True
- True : generate a plot of counts (integrated radially) vs proj. radius
- False : do not generate any plots
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.bdd, etc. It also stores the complete log and all input parameters. See Output Dictionary (rppiA) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='rppiA',outfn='redgalpairs') # generate default parameters cnt = gun.rppiA(gals, par, write=True, plot=True) # get pair counts and plot
-
gundam.
rppi_C
(tab, tab1, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given two astropy data tables, cross-count pairs in the projected-radial space (\(r_p\) vs \(\pi\) space).
All input parameters that control binning, cosmology, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with particles (sample D)
- tab1 : astropy table
- Table with particles (sample R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (rppiC) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=True
- True : generate a plot of counts (integrated radially) vs proj. radius
- False : do not generate any plots
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dr, counts.bdr, etc. It also stores the complete log and all input parameters. See Output Dictionary (rppiC) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table qsos = Table.read('qso.fits') # read data gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='rppiC',outfn='qso_rg_pairs') # generate default parameters cnt = gun.rppiC(qsos, gals, par, write=True) # get pair counts
Main Routines (redshift space)¶
These are the main routines for counting pairs and estimating correlation functions. They encapsulate all required steps to provide a true end-user experience : you only need to provide data+parameters to receive final results, and even that amazing plot of your next paper
-
gundam.
rcf
(tab, tab1, par, nthreads=-1, write=True, para=False, plot=False, **kwargs)¶ Given two astropy tables corresponding to data and random samples, this routine calculates the two-point redshift space auto-correlation function (rcf)
All input parameters that control binning, cosmology, estimator, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (rcf) for a detailed description
- write : bool. Default=True
- True : write several output files for DD/RR/DR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- para : bool. Default=False
- True : run in parallel over all available ipyparallel engines
- False : run serially. No need for any ipyparallel cluster running
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.rr, counts.xis, etc. It also stores the complete log and all input parameters. See Output Dictionary (rcf) for a detailed description.
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data rans = Table.read('redgal_rand.fits') # read randoms par = gun.packpars(kind='rcf',outfn='redgal') # generate default parameters cnt = gun.rcf(gals, rans, par, write=True, plot=True) # get rcf and plot
-
gundam.
rccf
(tab, tab1, tab2, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given three astropy tables corresponding to data, random, and cross samples, this routine calculates the two-point redshift space cross-correlation function (rccf)
All input parameters that control binning, cosmology, estimator, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- tab2 : astropy table
- Table with cross sample particles (C)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (rccf) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for CD/CR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.cd, counts.cr, counts.xis, etc. It also stores the complete log and all input parameters. See Output Dictionary (rccf) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data files rans = Table.read('redgal_rand.fits') qsos = Table.read('qso.fits') par = gun.packpars(kind='pccf',outfn='redgal_qso') # generate default parameters cnt = gun.rccf(gals, rans, qsos, par, write=True ) # get rcf
-
gundam.
s_A
(tab, par, nthreads=-1, write=True, para=False, plot=False, **kwargs)¶ Given an astropy table, count pairs in redshift space
All input parameters that control binning, cosmology, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (sA) for a detailed description
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- para : bool. Default=False
- True : run in parallel over all available ipyparallel engines
- False : run serially. No need for any ipyparallel cluster running
- plot : bool. Default=False
- True : generate a plot of counts vs redshift separation
- False : do not generate any plots
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.bdd, etc. It also stores the complete log and all input parameters. See Output Dictionary (sA) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='sA',outfn='redgalpairs') # generate default parameters cnt = gun.s_A(gals, par, write=True, plot=True ) # get counts and plot
-
gundam.
s_C
(tab, tab1, par, nthreads=-1, write=True, para=False, plot=False, **kwargs)¶ Given two astropy tables, cross-count pairs in redshift space
All input parameters that control binning, cosmology, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with particles (sample D)
- tab1 : astropy table
- Table with particles (sample R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (sC) for a detailed description
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- para : bool. Default=False
- True : run in parallel over all available ipyparallel engines
- False : run serially. No need for any ipyparallel cluster running
- plot : bool. Default=False
- True : generate a plot of counts vs redshift separation
- False : do not generate any plots
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dr, counts.bdr, etc. It also stores the complete log and all input parameters. See Output Dictionary (sC) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table qsos = Table.read('qso.fits') # read data gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='sC',outfn='qso_rg_pairs') # generate default parameters cnt = gun.rppiC(qsos, gals, par, write=True) # get pair counts
Main Routines (angular space)¶
These are the main routines for counting pairs and estimating correlation functions. They encapsulate all required steps to provide a true end-user experience : you only need to provide data+parameters to receive final results, and even that amazing plot of your next paper
-
gundam.
acf
(tab, tab1, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given two astropy tables corresponding to data and random samples, this routine calculates the two-point angular space auto-correlation function (acf)
All input parameters that control binning, cosmology, estimator, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (acf) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for DD/RR/DR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.rr, counts.xis, etc. It also stores the complete log and all input parameters. See Output Dictionary (acf) for a detailed description.
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data rans = Table.read('redgal_rand.fits') # read randoms par = gun.packpars(kind='rcf',outfn='redgal') # generate default parameters cnt = gun.rcf(gals, rans, par, write=True, plot=True) # get rcf and plot
-
gundam.
accf
(tab, tab1, tab2, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given three astropy tables corresponding to data, random, and cross samples, this routine calculates the two-point angular space cross-correlation function (accf)
All input parameters that control binning, estimator, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles (D)
- tab1 : astropy table
- Table with random particles (R)
- tab2 : astropy table
- Table with cross sample particles (C)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (accf) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for CD/CR counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate the CF plot on screen (and saved to disk when
write=True
) - False : do not generate any plots
- True : generate the CF plot on screen (and saved to disk when
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.cd, counts.cr, counts.wth, etc. It also stores the complete log and all input parameters. See Output Dictionary (accf) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data files rans = Table.read('redgal_rand.fits') qsos = Table.read('qso.fits') par = gun.packpars(kind='accf',outfn='redgal_qso') # generate default parameters cnt = gun.accf(gals, rans, qsos, par, write=True ) # get accf
-
gundam.
th_A
(tab, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given an astropy data table, count pairs in angular space
All input parameters that control binning, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with data particles
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (thA) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate a plot of counts vs angular separation
- False : do not generate any plots
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dd, counts.bdd, etc. It also stores the complete log and all input parameters. See Output Dictionary (thA) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='thA',outfn='redgalpairs') # generate default parameters cnt = gun.th_A(gals, par, write=True, plot=True ) # get counts and plot
-
gundam.
th_C
(tab, tab1, par, nthreads=-1, write=True, plot=False, **kwargs)¶ Given two astropy data tables, cross-count pairs in angular space
All input parameters that control binning, column names, etc. are passed in a single dictionary
par
, which can be easily generated with default values usinggundam.packpars()
and customized laterParameters
- tab : astropy table
- Table with particles (sample D)
- tab1 : astropy table
- Table with particles (sample R)
- par : Munch dictionary
- Input parameters. See Input Parameters Dictionary (thC) for a detailed description
- nthreads : integer. Default=-1
- Number of threads to use. Default to -1, which means all available cpu cores
- write : bool. Default=True
- True : write several output files for counts, correlations, etc.
- False : do not write any output files (except for logs, which are always saved)
- plot : bool. Default=False
- True : generate a plot of counts vs angular separation
- False : do not generate any plots
- kwargs : dict
- Extra keywords passed to
gundam.plotcf()
Returns
- counts : Munch dictionary
- Munch dictionary, containing all counts and correlations accesible by field keys, e.g. counts.dr, counts.bdr, etc. It also stores the complete log and all input parameters. See Output Dictionary (thC) for a detailed description
Examples
import gundam as gun ; from astropy.table import Table qsos = Table.read('qso.fits') # read data gals = Table.read('redgal.fits') # read data par = gun.packpars(kind='thC',outfn='qso_rg_pairs') # generate default parameters cnt = gun.th_C(qsos, gals, par, write=True) # get pair counts
Auxiliary Routines¶
These are helper routines that make it easy to work with correlations, visualize
count data, optimize grid size, sort data, etc. For examples of how to use them
check their individual help, or even better, the source code of gundam.pcf()
-
gundam.
bestSKgrid2d
(par, npts, ras, dens=None)¶ Try to find the optimum size (mxh1,mxh2) of a 2D skip grid, for an arbitrary sample of (ra,dec) points.
This is far from trivial, so for the moment, this routine works as follows:
- Choose an “optimum” target cell density (around 22 ppcell in many tests)
- Find the best
mxh1
from an empirical best fit relation of npts vs mxh1 - Adjust
mxh2
to reach the target cell density
Parameters
- par : Munch dictionary
- Input parameters dictionary
- npts : integer or list of integers
- Nr. of objects in the sample(s). For cross-correlations it should
be a list of the form [sample_D_size, sample_R,size] and
the effective
npts
adopted is their sum - ras : float array or list of arrays
- Right ascention of objects in the samples(s). For cross-correlations, the two samples are concatenated
- dens : float. Default=None
- Target cell density. Note that dens=None implicates a default value of 22 if npts>50000 and 16 otherwise.
Returns
- mxh1, mxh2 : integer
- Nr. of DEC and RA cells of the SK table, respectively
- dens : float
- The effective cell density adopted
-
gundam.
bestSKgrid3d
(par, npts, ras, dens=None)¶ Try to find the optimum size (mxh1,mxh2,mxh3) of a 3D skip grid, for an arbitrary sample of (ra,dec,dcom) points.
This is far from trivial, so for the moment, this routine works as follows:
- Choose an “optimum” target cell density according to the type of correlation
- Choose the best
mxh3
asmxh3=int((dcmax-dcmin)/rvmax)
- Find the best
mxh1
from an empirical best fit relation of npts vs mxh1 - Adjust
mxh2
to reach the target cell density
Parameters
- par : Munch dictionary
- Input parameters dictionary
- npts : integer or list of integers
- Nr. of objects in the sample(s). For cross-correlations it should
be a list of the form [sample_D_size, sample_R,size] and
the effective
npts
adopted is their sum - ras : float array or list of arrays
- Right ascention of objects in the samples(s). For cross-correlations, the two samples are joined together
- dens : float. Default=None
- Target cell density. Note that dens=None implicates different default values. See function code.
Returns
- mxh1, mxh2, mxh3 : integer
- Nr. of DEC, RA, DCOM cells of the SK table, respectively
- dens : float
- The effective cell density adopted
-
gundam.
cfindex
(path='./')¶ List file_name + description of all count (.cnt) files in a given path. Useful to quickly check a dir with dozens of correlation function runs and therefore hundreds of files.
Parameters
- path : string
- Path to list descriptions of .cnt files inside
-
gundam.
cnttable
(cfi, fmt1=None, fmt2=None, write=False, browser=True)¶ Shows a nicely formatted tabular view of the count data stored in a counts output dictionary. The table is printed in stdout and optionally displayed in the default web browser.
Parameters
- cfi : string or Munch dictionary
- Filepath for the counts (.cnt) file, or the counts dictionary itself
- fmt1 : string
- Ouput format of numeric fields (bins, correlations and errors). Default=’.4f’
- fmt2 : string
- Ouput format of numeric fields (counts). Default=’.2f’
- write : bool
- If
write=True
, save table to disk. Filename will be asked for. Default=False - browser : bool
- If
browser=True
, display HTML table in browser. Default=True
Returns
- tab : astropy table
- Single table with all relevant counts as columns. Use
print(tab)
ortab.show_in_browser()
Examples
# show info for a projected correlation from file on disk cnttable('/proj/galred.cnt') # show info a from variable in the session's memory cnttable(galred)
-
gundam.
makebins
(nsep, sepmin, dsep, logsep)¶ Create arrays of bins in which Gundam will count pairs and estimate correlation functions, given the number of bins desired, the minimum bin value and the chosen bin width.
Note units are not needed, but should be interpreted according to the input parameters
Parameters
- nsep : integer
- Number of bins
- sepmin : float
- Minimum bin location
- dsep : float
- Bin width (dex if
logsep=1
) - logsep : bool
- If
True
, do log-space binning. IfFalse
, do linear-space binning
Returns
- sep : float array
- Bin locations used by Gundam (left-side + extra bin at right limit)
- sepout : float array
- Left, middle and right-side points of each bin. Useful to plot more easily
Examples
# Create 18 log bins of size=0.2 dex in redshift and projected space seps,sepsout = makebins(18,0.01,0.2,1) sepp,seppout = makebins(18,0.01,0.2,1) # Create 25 bins of size 2 Mpc in radial space, from 0 to 50 Mpc sepv = makebins(25,0.,2.,0)[0] # Create instead 1 bin of size 50 Mpc, e.g. to work out the pcf integrated from 0 to 50 Mpc sepv = makebins(1,0.,50.,0)[0]
-
gundam.
pixsort
(tab, colnms, par)¶ Arrange an astropy table with (ra,dec) or (ra,dec,z) data into a grid of SK pixels and sort the rows of the table according to the pixel index, ordered by a given methodology. This way data in a given pixel sits closer in memory, largely increasing the efficiency of the cache.
Several (experimental) orders such as Morton and Hilbert are available
Parameters
- tab: astropy table
- Data table with ra/dec or ra/dec/z points
- colnms: list of strings
- Colum names for coordinates and redshift, e.g. [‘ra’,’dec’,’z’], [‘RAJ2000’,’DEJ2000’,’redshift’]
- par: Munch dictionary
- Input parameters dictionary. Must contain
mxh1
,mxh2
,mxh3
(when appropiate), the survey boundariessbound
and the desired method for orderingpxorder
. For samples that cross the RA=0 limit, it can is useful to specify a custom RA boundary. See Custom RA boundaries
Returns
- sidx : array
- Sort index. Use
tab=tab[sidx]
to actually sort the data
-
gundam.
qprint
(self)¶ Prints a quick, nicely formatted version of Munch dictionaries, such as the counts ouput dictionary or the par input dictionary used by Gundam routines. Very useful to quickly explore what is stored inside.
This routine is added dynamically to the Munch class, so it can be accessed as
any_munch_obj.qprint()
Note counts output dictionaries can also be explored using
gundam.cnttable()
to display a more elaborated tabular view of (exclusively) the counts calculated.Examples
import gundam as gun # Explore the many arrays and counts of a typical run cred = gun.readcounts('redgalaxies.cnt') cred.qprint() # Check the parameters used to create the run cred.par.qprint()
-
gundam.
radec2xyz
(ra, dec, r=0.5)¶ Converts a (ra,dec) coordinates to rectangular (x,y,z) coordinates in a sphere of radius r.
For
r=0.5
, this allows to speed up subsecuent haversine distance calculations between two points, by simply computing \(dhav^2 = (x1-x2)^2 + (y1-y2)^2 + (z1-z2)^2\)Parameters
- ra,dec : float arrays
- Right ascention an declination coordinates
- r : float. Default=0.5
- Sphere radius
Returns
- (x,y,z) : tuple
- Rectangular coordinates
-
gundam.
tpcf
(npt, nrpt, dd, bdd, rr, dr, estimator)¶ Return the (auto)correlation function for a given estimator and arrays of data and random counts
If DR counts are not needed (e.g. the ‘NAT’ estimator), just set
dr=0
If boostrap errors are not needed or available, just set
bdd
to a zero-valued array with null 2nd dimension, e.g.bdd=np.zeros([len(dd),0])
Parameters
- npt,nrpt : integer
- Number of data and random particles
- dd : float array
- DD counts
- bdd : float array
- Bootstrap DD counts
- rr : float array
- RR counts in projected and radial separations
- dr : float array
- DR counts in projected and radial separations
- estimator : string
Statistical estimator of the correlation function. Default=`NAT`
- ‘NAT’ : Natural -> \(DD/RR - 1\)
- ‘HAM’ : Hamilton -> \(DD*RR/DR^{2} - 1\)
- ‘LS’ : Landy-Szalay -> \((DD - 2DR + RR) / RR\)
- ‘DP’ : Davis-Peebles -> \(DD/DR - 1\)
Returns
- xi : float array
- Correlation function
- xierr : float array
- Boostrap error estimate. Set to zero if
bdd
is nulled as explained above
Notes
See this paper for a nice review on estimators and their normalization factors. Here, the normalization factors are derived to : (1) keep estimator formulae clean, (2) avoid having operations such as (npt*(npt-1)) * dd, where counts are multiplied/divided by very big numbers when npt is large.
Examples
# Calculate the angular CF using the Landy-Szalay estimator acf, acferr = gun.tpcf(npt,nrpt,dd,bdd,rr,dr,estimator='LS')
-
gundam.
tpccf
(npt, nrpt, cd, bcd, cr, estimator)¶ Return the (cross)correlation function for a given estimator and count arrays for data (D), random (R) and cross (C) samples.
For the moment the only estimator implemented is the Davis-Peebles : \(\xi=CD/CR-1\)
If bootstrap errors are not needed or available, just set
bdd
to a zero-valued array, e.g.bdd=np.zeros([len(dd),0])
Parameters
- npt,nrpt : integer
- Number of particles in data (D) and random (R) samples
- cd : float array
- CD counts
- bcd : float array
- Bootstrap CD counts
- cr : float array
- CR counts
- estimator : string
- ‘DP’ : Davis-Peebles -> \(CD/CR-1\)
Notes
C and D are data samples while R is random sample corresponding to D
Returns
- fxi : float array
- Cross-correlation function
- fxierr : float array
- Boostrap error estimates
Examples
import gundam as gun c = gun.readcounts('qso_gal.cnt') (ccf,ccferr) = tpccf(c.npt, c.nrpt, c.cd, c.bcd, c.cr, estimator='DP')
-
gundam.
tpcf_wrp
(npt, nrpt, dd, bdd, rr, dr, dsepv, estimator)¶ Return the projected (auto)correlation function for a given estimator and arrays of data and random counts
If DR counts are not needed (e.g. the ‘NAT’ estimator), just set
dr=0
If boostrap errors are not needed or available, just set
bdd
to a zero-valued array with null 2nd dimension, e.g.bdd=np.zeros([len(dd),0])
Parameters
- npt,nrpt : integer
- Number of data and random particles
- dd : float array
- DD counts in projected and radial separations
- bdd : float array
- Bootstrap DD counts in projected and radial separations
- rr : float array
- RR counts in projected and radial separations
- dr : float array
- DR counts in projected and radial separations
- dsepv : float
- Bin size in radial direction
- estimator : string
Statistical estimator of the correlation function
- ‘NAT’ : Natural -> \(DD/RR - 1\)
- ‘HAM’ : Hamilton -> \(DD*RR/DR^{2} - 1\)
- ‘LS’ : Landy-Szalay -> \((DD - 2DR + RR) / RR\)
- ‘DP’ : Davis-Peebles -> \(DD/DR - 1\)
Returns
- wrp : float array
- Correlation function
- wrperr : float array
- Boostrap error estimate. Set to zero if
bdd
is nulled as explained above
Notes
See this paper for a nice review on estimators and their normalization factors. Here, the normalization factors are derived to : (1) keep estimator formulae clean, (2) avoid having operations such as (npt*(npt-1)) * dd, where counts are multiplied/divided by very big numbers when npt is large.
Examples (xxxx TODO)
(wrp,wrperr) = tpcf_wrp(npt,nrpt,ddpv,bddpv,rrpv,drpv,dsepv,estimator='HAM')
-
gundam.
tpccf_wrp
(npt, nrpt, cd, bcd, cr, dsepv, estimator)¶ Return the projected (cross)correlation function for a given estimator and count arrays
If boostrap errors are not needed or available, just set
bcd
to a zero-valued array with null 2nd dimension, e.g.bcd=np.zeros([len(ddXXX),0])
Parameters
- npt : integer
- Number of data particles
- nrpt : integer
- Number of random particles
- cd : float array
- CD counts in projected and radial separations
- bcd : float array
- Bootstrap CD counts in projected and radial separations
- cr : float array
- CR counts in projected and radial separations
- dsepv : float
- Radial bin size
- estimator : string
Statistical estimator of the correlation function
- ‘DP’ : Davis-Peebles -> \(CD/CR - 1\) (C,D data samples, R is random of D)
Returns
- wrp : float array
- Projected cross-correlation function
- wrperr : float array
- Boostrap error estimate
Examples
# remains to do XXXXX (wrp,wrperr) = tpccf_wrp(npt,nrpt,cd,bcd,cr,dsepv,estimator='DP')
Plotting Routines¶
-
gundam.
cntplot
(cnt, **kwargs)¶ Plot a correlation function from a counts output dictionary (either read from disk or passed directly). Both axes are set to log-space and axes labels are selected automatically according to the type of correlation (i.e. given by
par.kind
)This is a wrapper for
gundam.plotcf()
, so all of its parameters can be specified too.Parameters
- cnt : string or Munch dictionary
- Filepath for the counts (.cnt) file, or the counts dictionary itself
- kwargs : keyword list
- Any extra [key]=value pairs are passed to the underlying
gundam.plotcf()
routine
Examples
import gundam as gun # Read a pcf run and plot the correlation function cnt1 = gun.readcounts('/p01/redgalsP.cnt') cntplot(cnt1) # Plot the correlation function from a .cnt file cntplot('/p01/redgalsA.cnt', label='angcf of redgals', fill=True)
-
gundam.
cntplot2D
(cnt, estimator=None, slevel=5, write=False, figfile=None, xlabel=None, ylabel=None, cmap='jet', **kwargs)¶ Plot the 2D correlation function in the projected-radial space (\(r_p\) vs \(\pi\) space) with optional gaussian smoothing and contour levels
Parameters
- cnt : string or Munch dictionary
- Filepath for the counts (.cnt) file or the counts output dictionary
- estimator : string. Default=None
- Estimator for the correlation function. Any of (‘NAT’,’LS’,’HAM’,’DP’).
If
estimator=None
, then it is taken fromcnt.par.estimator
- slevel : float. Default=5
- Smoothing level (namely the size of the Gaussian smothing kernel)
- write : bool. Default=False
- Save the figure to disk (default format is pdf). See Notes to save in other graphic formats
- figfile : string. Default=None
- Specify an alternative file name for the figure. If
None
, then choosecnt.par.outfn
as default. Do not add extension. - xlabel, ylabel : string. Default=None
- X-axis and Y-axis labels. If supplied, they override the default labels (\(r_p \ [h^{-1} Mpc]\) and \(\pi \ [h^{-1} Mpc]\))
- cmap : string. Default=’jet’
- Colormap for the plot
- kwargs : keyword list
- Any extra [key]=value pairs are passed to
matplolib.pyplot.pcolor()
Use this to customize shading, edges, alpha, etc.
Notes
- The graphic format can be changed by passing the
figformat
key inkwargs
, e.g.figformat='pdf'
. Any format supported by matplotlib is valid.
Examples
# Check some nice Fingers of God and the Kaiser squashing cntplot2D('lum_red_gals.cnt', cmap='viridis')
-
gundam.
comparecf
(clist1, clist2=None, shift=0.0, fac=1.0, ploterrbar=True, fill=False, filtneg=False, label=None, plotratio=False, ratioxrange=None, color1=None, marker1=None, markers1=None, linestyle1=None, linewidth1=None, color2=None, marker2=None, markers2=None, linestyle2=None, linewidth2=None, f=None, ax1=None, ax2=None)¶ Plot multiple correlation functions in a single figure for easy comparison. Optionally show an additional lower panel displaying the ratio of each function respect to a single “control” correlation (many-to-one) or to multiple correlations (one-to-one).
Parameters
- clist1 : list of Munch dictionaries / list of strings
- Count dictionaries or filepaths of .cnt files of correlations
- clist2 : list of Munch dictionaries / list of strings. Default=None
- List of control correlations. When
plotratio=True
, the y-values of each correlation curve inclist1
are divided by those inclist2
(one-to-one) and plotted in a lower panel. Ifclist2
has a single element, the ratios are from allclist1
divided by the singleclist1
. See Notes for more details - shift : float. Default=0.0
- Fraction of bin size by which
x
values are shifted. Useful to slightly separate overlapping curves - fac : float. Default=1.0
- Multiplication factor for
y
andyerr
- ploterrbar : bool. Default=True
- If
ploterrbar=True
, plot error bars according toyerr
- fill : bool. Default=False
- If
fill=True
, plot a filled semi-transparent error region instead of the usual error bars - filtneg : bool. Default=False
- If
filtneg=True
, filter out points where (y-yerr)<0, i.e. those with large errors in a log plot - label: list
- Optional list of strings to label each correlation function. If ommited,
the values are taken from the
outfn
key stored in the count objects - plotratio : bool. Default=False
- If
plotratio=True
, plot also a lower panel with ratios ofclist1
respect toclist2
- ratioxrange : list of form [xmin,xmax]
- Only plot the ratio between xmin and xmax
- color1,marker1,markers1,linestyle1,linewidth1 : lists
- List of colors, marker symbols, marker sizes, line styles and line widths
for curves in
clist1
- color2,marker2,markers2,linestyle2,linewidth2 : lists
- List of colors, marker symbols, marker sizes, line styles and line widths
for control curves in
clist2
- f : figure instance
- Handle of existing Figure instance
- ax1,ax2 : axes instances
- Handles of correlation plot and ratio plot axes
Notes
The correlation curves in
clist2
are not plotted in the correlation function panel while curves present in both clists are not shown in the ratio panel (i.e. to avoid ratios of curves respect to themselves).Returns
- (f,ax1,ax2) or (f,ax1) : tuple of handles
- Handles of figure, correlation axis (ax1) and ratio axis (ax2), if present
Examples
# Compare two w(rp) correlations comparecf(['galred.cnt', 'galblue.cnt']) # Compare one acf on disk with another passed as a counts ouput dictionary comparecf(['/proj/galred.cnt', qso]) # Compare multiple samples and plot sample/control ratios f,ax1,ax2 = comparecf(['galred.cnt', 'galblue.cnt', 'galgreen.cnt'], clist2=['allgals.cnt'], fill=True, plotratio=True) # Add another curve to previous plot comparecf(['qso.cnt'], clist2=['control_qso.cnt'], color2=['k'], f=f, ax1=ax1, ax2=ax2, plotratio=True)
-
gundam.
fitpowerlaw
(x, y, yerr, iguess=[1.0, -1.0], fitrange=None, plot=False, markfitrange=False, **kwargs)¶ Fit a power-law of the form \(ax^{\gamma}\) to a correlation function over a given x-coordinate range. Optionally plot the fitted curve
Parameters
- x : float array
- x-coordinates such as cnt.rpm, cnt.thm, etc.
- y : float array
- x-coordinates such as cnt.wrp, cnt.wth, etc.
- yerr : float array
- Errors in y-coordinates such as cnt.wrperr, cnt.wtherr, etc.
- iguess : list of floats. Default=[1., -1.]
- Initial guesses for \(a\) and \(\gamma\)
- fitrange : float array of form [xmin,xmax]. Default=None
- Restrict fit to points inside the given interval
- plot : bool. Default=False
- Plot the fitted power-law curve
- markfitrange : bool. Default=False
- Overlay marks for points actually used in the fitting
- kwargs : keyword list
- Any extra [key]=value pairs are passed to
matplolib.pyplot.plot()
Use this to customize colors, linestyles, markers, etc.
Examples
import gundam as gun c1 = gun.readcounts('galaxies.cnt') cntplot(c1) gun.fitpowerlaw(c1.rpm, c1.wrp, c1.wrperr, plot=True)
-
gundam.
plotcf
(x, y, yerr, fac=1.0, write=False, figfile=None, par=None, angunit='deg', xlabel=None, ylabel=None, label=None, shift=0.0, ploterrbar=True, fill=False, filtneg=False, **kwargs)¶ Plot a correlation function from arrays of x, y and yerr values. Both axes are set to log-space and axes labels are selected automatically according to the type of correlation (i.e. given by
par.kind
)Parameters
- x,y,yerr : float arrays
- x, y coordinates and corresponding errors of the correlation function.
If
yerr=0
or all elements of yerr are <=0 orploterrbar=False
, no errorbar is plotted - fac : float. Default=1.0
- Multiplication factor for
y
andyerr
- write : bool. Default=False
- Save the figure to disk (default format is pdf). See Notes to save in other graphic formats
- figfile : string. Default=None
- Specify an alternative file name for the figure. If specified, overrides the
default which is to take it from
par.outfn
. Do not add extension. - par : dictionary of type Munch. Default=None
- Used to pass
outfn
name to name saved figures whenwrite=True
- angunit : string. Default=’deg’
- ‘arcsec’ : set ouput axis in arcsec (
x
values are unscaled) - ‘arcmin’ : set ouput axis in arcmin (
x
values are scaled asx
/60) - ‘deg’ : set ouput axis in degrees (
x
values are scaled asx
/3600)
- ‘arcsec’ : set ouput axis in arcsec (
- xlabel, ylabel : string. Default=None
- X-axis and Y-axis labels. If supplied, they override the default labels
deduced from
par.kind
- label : string. Default=None
- Label for the legend of the curve. If supplied, it will override the
default label which is the basename of
par.outfn
. Note you have to issue at least a plt.legend() to actually display the legend box - shift : float. Default=0.0
- Fraction of bin size by which
x
values are shifted. Useful to slightly separate overlapping curves - ploterrbar : bool. Default=True
- If
ploterrbar=True
, plot error bars according toyerr
- fill : bool. Default=False
- If
fill=True
, plot a filled semi-transparent error region instead of the usual error bars - filtneg : bool. Default=False
- If
filtneg=True
, filter out points where (y-yerr)<0, i.e. those with large errors in a log plot - kwargs : keyword list
- Any extra [key]=value pairs are passed to the underlying
matplotlib.pyplot.plot()
routine, except foralpha
which is passed topyplot.fill_between()
,capsize
which is passed topyplot.errorbar()
, andfigformat
which is passed topyplot.savefig()
. Use this to customize colors, linestyles, markers, etc.
Notes
- Sucessive calls cycle between 4 predefined styles (for color, marker,
linewidth, etc.) that can be overrriden by passing the corresponding
[key]=value pairs in
kwargs
- The output graphic format can be changed by passing the
figformat
key inkwargs
, e.g.figformat='pdf'
. Any format supported by matplotlib is valid.
Examples
import gundam as gun c1 = gun.readcounts('redgalPCF.cnt') c2 = gun.readcounts('redgalRCF.cnt') plt.figure(1) # Plot w(rp) gun.plotcf(c1.rpm,c1.wrp,c1.wrperr,par=c1.par) plt.figure(2) # Plot xi(s) gun.plotcf(c2.sm,c2.xis,c2.xiserr,par=c2.par,color='yellow')
Comprehensive Auxiliary Routines¶
These are helper routines that perform many common tasks (in the Python side)
during a correlation run, such as collecting counts, i/o functions, initialization,
logging, etc. For a nice example of how to use them, see source code of
gundam.pcf()
-
gundam.
addlog
(file, key, m)¶ Read a text file and dump it into a key of a given dictionary. Useful to add entire logs from disk into Munch objects
Parameters
- file : string
- Complete path of file, usually any log file
- key : string
- Name of the key to be created
- m : Munch dictionary
- The dictionary where the key
m.key
will be created
-
gundam.
allequal
(v)¶ Fast way to check if all elements of a 1D-array have the same value. Useful to detect when all weights are set to 1.0, and hence to call faster versions of the counting routines
Parameters
- v : array_like
- Array to be checked
Returns
- res : bool
- True if all elements have the same value
-
gundam.
addpvcols
(x, cfo, basecolname, **kwargs)¶ Auxiliary function used by
gundam.cnttable()
to append columns to a table populated with the fields of counts ouput dictionary that store pair counts, all with nice column names. Works with 1d counts or 2d counts arrays (e.g. those from a pcf run whennsepv>1
).Parameters
- x : astropy table
- Table to add data
- cfo : Much dictionary
- Counts dictionary with the count arrays, e.g. cfo.dd, cfo.rr, etc.
- basecolname : string
- The name of the field to add, e.g. dd, and also the prefix for the column name, which if needed will be appended with _001, _002, etc. for each radial bin
- kwargs :
- Any [key]=value pair to pass to the astropy Column constructor. Intended
to pass a format specification for the column, such as
format='.4f'
Returns
None, it modifies the input table
x
-
gundam.
bound2d
(decs)¶ Return the (maximal) survey boundaries in (RA,DEC) for multiple samples. Note in the RA direction, the limits are always set as ramin=0. and ramax=360. to make the pair counting valid for surveys that cross the origin of coordinates
Parameters
- decs : float array or list of arrays
- DEC coordinates of one or more samples [deg]
Returns
- bound : tuple
- The limits (ramin,ramax,decmin,decmax) that enclose all input samples
-
gundam.
bound3d
(decs, dcs)¶ Return the (maximal) survey boundaries in (RA,DEC,DCOM) for multiple samples. Note in the RA direction, the limits are always set as ramin=0. and ramax=360. to make the pair counting valid for surveys that cross the origin of coordinates
Parameters
- decs : float array or list of arrays
- DEC coordinates of one or more samples [deg]
- dcs : float array or list of arrays
- Comoving distances of one or more samples [Mpc/h]
Returns
- bound : tuple
- The limits (ramin,ramax,decmin,decmax,dcmin,dcmax) that enclose all input samples
-
gundam.
cross0guess
(ra)¶ Guess if a set of RA coordinates cross the RA=0 division (by finding one source 1deg to the left and another 1deg to the right of RA=0, at least)
Parameters
- ra : array
- Right ascention coordiantes
Returns
- res : bool
- True if the sample seems to cross the RA=0 boundary
-
gundam.
buildoutput
(par, npts=[], binslmr=[], dd=None, rr=None, dr=None, cd=None, cr=None, bootc=None, cf=None, cferr=None)¶ Given a set of pair count arrays, assemble the final counts output dictionary, adding also bins, input parameters and relevant sample sizes.
This is for the main Gundam functions that calculate a multiple pair counts, i.e.
gundam.pcf()
,gundam.pccf()
,gundam.rcf()
,gundam.rccf()
,gundam.acf()
,gundam.accf()
Parameters
- par : Munch dictionary
- Input parameters
- npts : list of int. Default=[]
- A list with the nr. of points in each sample, i.e. [npt, npt1]. For example npts=[400, 4000] in cross-count cases; npts=[400] for auto-count cases
- binslmr : list. Default=[]
- The left, mid and right-side locations of each bin, returned of example
by
gundam.makebins()
- dd,rr,dr : array. Default=None
- The “dd”, “rr” and “dr” count arrays, if available
- cd,cr : array. Default=None
- The “cd” and “cr” count arrays, if available
- bootc : array. Default=None
- The boostrap count array, if available
- cf,cferr : array. Default=None
- The correlation function and its error
Returns
- counts : Munch dictionary
- The ouput dictionary
-
gundam.
buildoutputC
(par, npts=[], binslmr=[], dd=None, dr=None, bootc=None, intpi=None, intpib=None)¶ Given a set of pair count arrays, assemble the final counts output dictionary, adding also bins, input parameters and relevant sample sizes.
This is for the main Gundam functions that calculate a single pair count, i.e.
gundam.rppi_A()
,gundam.rppi_C()
,gundam.s_A()
,gundam.s_C()
,gundam.th_A()
,gundam.th_C()
Parameters
- par : Munch dictionary
- Input parameters
- npts : list of int. Default=[]
- A list with the nr. of points in each sample, i.e. [npt, npt1]. For example npts=[400, 4000] in cross-count cases; npts=[400] for auto-count cases
- binslmr : list. Default=[]
- The left, mid and right-side locations of each bin, returned of example
by
gundam.makebins()
- dd,dr : array. Default=None
- The “dd” and “dr” count arrays, if available
- bootc : array. Default=None
- The boostrap count array, if available
- intpi, intpib : array. Default=None
- “dd” counts and boostrap “dd” counts integrated along all radial bins, if available
Returns
- counts : Munch dictionary
- The ouput dictionary
-
gundam.
check_kind
(par, kind)¶ Check that
par.kind = kind
. Useful to test, for example, if you are passing the right par object before really counting pairsParameters
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- kind : string
- The kind to check against, e.g. ‘pcf’, ‘accf’, ‘rppiA’, etc.
-
gundam.
closelog
(log, runspyder=True)¶ Close the log machinery used by the main counting routines by removing the handlers
Parameters
- log : logging object
- The log object
- runspyder : bool. Default=True
- If
runspyder=True
, remove the extra handler for stdout added when running under a Spyder console
-
gundam.
finalize
(log, logf, logff, Ltime, t0, counts)¶ Perform some finalization tasks common to all correlation runs, by logging loop/total times and adding the contents of log files to the counts output dictionary
Parameters
- log : logging object
- The log object
- logfile : string
- The complete path of the .log file
- Ltime, Ttime : float
- Loop time and total compute time of a correlation run
- counts : Munch dictionary
- Output dictionary containing all counts and correlations
-
gundam.
initialize
(kind, par, nthreads=None, write=None, plot=None)¶ Perform some initialization tasks common to all correlation function runs, such as reading par from disk file if needed, check the parameter kind, find out if running under Spyder, initialize the logs, etc.
Parameters
- kind : string
- The kind to check against, e.g. ‘pcf’, ‘accf’, ‘rppiA’, etc.
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- nthreads : integer.
- Number of threads to use. Passed just to store it under par
- write : bool. Default=True
- Flag to generate output count files. Passed just to store it under par
- plot : bool. Default=True
- Flag to generate plot. Passed just to store it under par
Returns
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- log : logging object
- The log object
- logfile, logfilefort : string
- The complete path of the .log file and the .fotran.log file
- runspyder : bool. Default=True
- Detect if running under a Spyder console
-
gundam.
logcallinfo
(log, par, npts=[])¶ Write to log some useful runtime parameters of the main counting routines
Parameters
- log : logging object
- The log object
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- npts : list of 1, 2 or 3 integers
- The nr of objects in the input data table, random table, and cross sample table
-
gundam.
logtimming
(log, cntid, t)¶ Write a message to a log instance, showing the compute time for the counts identified with a certain ID string
Parameters
- log : logging object
- Log object
- cntid : string
- String to ID a certain type of counts, e.g. ‘DD’, ‘RR’, DR’, etc.
- t : float
- The time elapsed [s]
-
gundam.
readcounts
(cfile, silent=False)¶ Read from disk the counts dictionary generated by the main counting routines.
Parameters
- cfile : string
- Filepath for the counts (.cnt) file
- silent : bool
- If False, print a status message indicating the file was read. Default=False
Returns
- counts : Munch dictionary
- The counts object
-
gundam.
readpars
(filename)¶ Load from a JSON (.par) file the input parameters dictionary used by many Gundam routines.
Parameters
- filename : string
- Filepath of .par file
-
gundam.
savecounts
(cnt, altname=None)¶ Save to disk the counts output dictionary returned by the main counting routines.
The default file name is
cnt.par.outfn
+ .cnt, which can be overriden asaltname
+ .cnt, if suppliedParameters
- cnt : Munch dictionary
- The counts object
- altname : string. Default=None
- If supplied, use an alternative file name instead of
cnt.par.outfn
-
gundam.
savepars
(par, altname=None)¶ Save the parameters dictionary par, such as the one generated by
gundam.packpars()
, in a JSON file. By default it is named aspar.outfn
+ .parParameters
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- altname : string. Default=None
- If supplied, use an alternative file name instead of
par.outfn
+ .par
Examples
import gundam as gun # Get default values for an angular CF run and save to disk par = gun.packpars(kind='acf', outfn='/proj/acfrun01') gun.savepars(par)
-
gundam.
setlogs
(par, runspyder=True)¶ Set up the log machinery used by the main counting routines by creating the logger object, adding the required handlers and cleaning previous logs if present
Parameters
- par : Munch dictionary
- Input parameters dictionary for Gundam routines
- runspyder : bool. Default=True
- If
runspyder=True
, add an extra handler for stdout when running under a Spyder console
Returns
- log : logging object
- The log object
- logfile : string
- The complete path of the .log file
-
gundam.
writeasc_cf
(lb, mb, rb, f, ferr, par, fmt='%17.5f', altname=None)¶ Write an ASCII file for w(rp) / xi(s) / w(th) ouput counts produced by the code
Parameters
- lb,mb,rb : float arrays
- Left, mid and rigth-side of bins
- f, ferr : float arrays
- Correlation function and its error
- par : Munch dictionary
- Used to pass
par.outfn
to name the output file - fmt : string. Default=’%17.5f’
- Numeric formatting string
- altname : string. Default=None
- If supplied, use an alternative file name instead of
par.outfn
Examples
import gundam as gun c1 = gun.readcounts('redgals.cnt') writeasc_cf(c1.rpl, c1.rpm, c1.rpr, c1.wrp, c1.wrperr, c1.par)
-
gundam.
writeasc_rppicounts
(lb, mb, rb, rppi, par, fmt='%17.5f', cntid=None, altname=None)¶ Write an ASCII file for a rp-pi count array produced by the code. This is a 2D array of counts in the projected (rp) and radial (pi) directions.
The columns in the output will be [lb mb rb tot_counts rppi] where the first 3 are the left, mid and right-side of bins, tot_counts are the counts integrated for all radial bins, and rppi has one column for each radial bin
Parameters
- lb,mb,rb : float
- Left, mid and rigth-side of bins
- rppi : float array
- 2-dimensional count array. Usually this is one of the fields cnt.dd, cnt.rr, etc. of a projected correlation run
- par : Munch dictionary
- Used to pass various data, including
par.outfn
to name the output file - fmt : string. Default=’%17.5f’
- Numeric formatting string
- cntid : string. Default=None
- ID string for column headers. Usually can be ‘dd’, ‘rr’, ‘dr’, etc.
Also appended as the extension of the ouput file (when
altname=None
) - altname : string. Default=None
- If supplied, use an alternative file name instead of
par.outfn
+ .cntid
Examples
import gundam as gun c1 = gun.readcounts('redgals.cnt') c1.par.outfn '/home/myuser/sdss/redgals' # Write the DD counts in rp-pi dimensions gun.writeasc_rppicounts(c1.rpl, c1.rpm, c1.rpr, c1.dd, c1.par, cntid='dd') # Inspect the output file with open('redgals.dd', 'r') as f: print(f.read(), end="") # lb mb rb dd dd_001 dd_002 dd_003 ... # 0.10000 0.12417 0.14835 11509.00 2082.00 1500.00 1168.00 ... # 0.14835 0.18421 0.22007 20273.00 3122.00 2378.00 1899.00 ... # 0.22007 0.27327 0.32647 36169.00 4940.00 3845.00 3283.00 ... # 0.32647 0.40539 0.48431 64866.00 8453.00 6302.00 5236.00 ... # ...
-
gundam.
writeasc_counts
(lb, mb, rb, c, par, fmt='%17.5f', cntid=None, altname=None)¶ Write an ASCII file for a (1-dimensional) counts array produced by the code.
The columns in the output will be [lb mb rb c] where the first 3 are the left, mid and right-side of bins and c are the counts
Parameters
- lb,mb,rb : float
- Left, mid and rigth-side of bins
- c : float array
- Counts array. Usually this is one of the fields cnt.dd, cnt.dr, cnt.rr, etc. of a correlation run
- par : Munch dictionary
- Used to pass
par.outfn
to name the output file - fmt : string. Default=’%17.5f’
- Numeric formatting string
- cntid : string. Default=None
- ID string for column header. Usually can be ‘dd’, ‘rr’, ‘dr’, etc.
Also appended as the extension of the ouput file (when
altname=None
) - altname : string. Default=None
- If supplied, use an alternative file name instead of
par.outfn
+ .cntid
Examples
import gundam as gun c1 = gun.readcounts('bluegals.cnt') # Write the DD counts in angular dimensions. Use an alternative file name gun.writeasc_counts(c1.thl, c1.thm, c1.thr, c1.dd, c1.par, cntid='dd', altname='akounts') # Inspect the output file with open('akounts', 'r') as f: print(f.read(), end="") # lb mb rb dd # 0.01000 0.01206 0.01413 3178.00 # 0.01413 0.01704 0.01995 6198.00 # 0.01995 0.02407 0.02818 12765.00 # 0.02818 0.03400 0.03981 24888.00 # 0.03981 0.04802 0.05623 49863.00 # 0.05623 0.06783 0.07943 98883.00 ...
Fortran Wrapper Routines¶
These are the wrappers that call Fortran pair counting routines. They are intended to: (1) choose the fastest counting routine depending whether weights, bootstrap errors, etc. are requested or not; and (2): concentrate all mayor Fortran calling in one place
Unless you are building your own custom script, heavily modifying, or extending the core algorithms, you normally should not need to use these functions
-
gundam.
pairs_auto
(par, wunit, logff, tab, x, y, z, sk, ll, dc=None)¶ Wrapper for calling Fortran counting routines (for auto-pairs)
This function isolates the call of external Fortran routines, choosing the correct one based on geometry, and choosing the fastest depending whether, weights, bootstrap errors, etc. are requested or not.
Parameters
- par : Munch dictionary
- Input parameters
- wunit : bool
- True : all weights are equal to 1
- False : at least one weight is not 1
- logff : string
- File for logging the output of Fortran routines
- tab : astropy table
- Table with data particles
- x,y,z : arrays
- Rectangular coordinates of data particles
- sk,ll : arrays
- Skip table (SK) and linked list table (see
gundam.skll2d()
orgundam.skll3d()
) - dc : array [optional]
- Array of comoving distances. Not needed for angular counts
Returns
- tt : list of ndarray
- Ouput counts as returned by Fortran counting routines
-
gundam.
pairs_cross
(par, wunit, logff, tab, x, y, z, tab1, x1, y1, z1, sk1, ll1, dc=None, dc1=None)¶ Wrapper for calling Fortran counting routines (for cross-pairs among two tables called D(data) and Random(R))
This function isolates the call of external Fortran routines, choosing the correct one based on geometry, and choosing the fastest depending whether, weights, bootstrap errors, etc. are requested or not.
Parameters
- par : Munch dictionary
- Input parameters
- wunit : bool
- True : all weights are equal to 1
- False : at least one weight is not 1
- logff : string
- File for logging the output of Fortran routines
- tab,tab1 : astropy tables
- Table with data particles (D) and random particles (R), respectively
- x,y,z,x1,y1,z1 : arrays
- Rectangular coordinates of particles in D table and R rable, respectively
- tab1 : astropy table
- Table with data particles (D)
- sk1,ll1 : arrays
- Skip table (SK) and linked list table (see
gundam.skll2d()
orgundam.skll3d()
) - dc,dc1 : arrays [optional]
- Array of comoving distances of particles in D and R tables. Not needed for angular counts
Returns
- tt : list of ndarray
- Ouput counts as returned by Fortran counting routines