API: suave-specific routines

Corrfunc.theory

Wrapper for all clustering statistic calculations on galaxies in a simulation volume.

Corrfunc.theory.DDsmu(autocorr, nthreads, binfile, mu_max, nmu_bins, X1, Y1, Z1, weights1=None, periodic=True, X2=None, Y2=None, Z2=None, weights2=None, verbose=False, boxsize=0.0, output_savg=False, fast_divide_and_NR_steps=0, xbin_refine_factor=2, ybin_refine_factor=2, zbin_refine_factor=1, max_cells_per_dim=100, copy_particles=True, enable_min_sep_opt=True, c_api_timer=False, isa='fastest', weight_type=None, proj_type=None, ncomponents=None, projfn=None)[source]

Calculate the 2-D component values (e.g. pair-counts for the tophat basis) corresponding to the redshift-space correlation function, \(\xi(s, \mu)\) Pairs which are separated by less than the s bins (specified in binfile) in 3-D, and less than s*mu_max in the Z-dimension are counted.

If weights are provided, the mean pair weight is stored in the "weightavg" field of the results array. The raw pair counts in the "npairs" field are not weighted. The weighting scheme depends on weight_type.

To use the projection capability with suave, set the proj_type parameter for the desired basis functions, and set ncomponents and projfn accordingly.

Note

This module only returns pair counts and not the actual correlation function \(\xi(s, \mu)\). See the utilities Corrfunc.utils.convert_3d_counts_to_cf for computing \(\xi(s, \mu)\) from the pair counts.

New in version 2.1.0.

Parameters
  • autocorr (boolean, required) – Boolean flag for auto/cross-correlation. If autocorr is set to 1, then the second set of particle positions are not required.

  • nthreads (integer) – The number of OpenMP threads to use. Has no effect if OpenMP was not enabled during library compilation.

  • binfile (string or an list/array of floats) –

    For string input: filename specifying the s bins for DDsmu_mocks. The file should contain white-space separated values of (smin, smax) specifying each s bin wanted. The bins need to be contiguous and sorted in increasing order (smallest bins come first).

    For array-like input: A sequence of s values that provides the bin-edges. For example, np.logspace(np.log10(0.1), np.log10(10.0), 15) is a valid input specifying 14 (logarithmic) bins between 0.1 and 10.0. This array does not need to be sorted.

  • mu_max (double. Must be in range (0.0, 1.0]) –

    A double-precision value for the maximum cosine of the angular separation from the line of sight (LOS). Here, LOS is taken to be along the Z direction.

    Note: Only pairs with \(0 <= \cos(\theta_{LOS}) < \mu_{max}\) are counted (no equality).

  • nmu_bins (int) – The number of linear mu bins, with the bins ranging from from (0, \(\mu_{max}\))

  • X1/Y1/Z1 (array-like, real (float/double)) – The array of X/Y/Z positions for the first set of points. Calculations are done in the precision of the supplied arrays.

  • weights1 (array_like, real (float/double), optional) – A scalar, or an array of weights of shape (n_weights, n_positions) or (n_positions,). weight_type specifies how these weights are used; results are returned in the weightavg field. If only one of weights1 and weights2 is specified, the other will be set to uniform weights.

  • periodic (boolean) – Boolean flag to indicate periodic boundary conditions.

  • X2/Y2/Z2 (array-like, real (float/double)) – Array of XYZ positions for the second set of points. Must be the same precision as the X1/Y1/Z1 arrays. Only required when autocorr==0.

  • weights2 (array-like, real (float/double), optional) – Same as weights1, but for the second set of positions

  • verbose (boolean (default false)) – Boolean flag to control output of informational messages

  • boxsize (double) – The side-length of the cube in the cosmological simulation. Present to facilitate exact calculations for periodic wrapping. If boxsize is not supplied, then the wrapping is done based on the maximum difference within each dimension of the X/Y/Z arrays.

  • output_savg (boolean (default false)) – Boolean flag to output the average s for each bin. Code will run slower if you set this flag. Also, note, if you are calculating in single-precision, s will suffer from numerical loss of precision and can not be trusted. If you need accurate s values, then pass in double precision arrays for the particle positions.

  • fast_divide_and_NR_steps (integer (default 0)) – Replaces the division in AVX implementation with an approximate reciprocal, followed by fast_divide_and_NR_steps of Newton-Raphson. Can improve runtime by ~15-20% on older computers. Value of 0 uses the standard division operation.

  • (xyz)bin_refine_factor (integer (default (2,2,1) typical values in [1-3])) – Controls the refinement on the cell sizes. Can have up to a 20% impact on runtime.

  • max_cells_per_dim (integer (default 100, typical values in [50-300])) – Controls the maximum number of cells per dimension. Total number of cells can be up to (max_cells_per_dim)^3. Only increase if rmax is too small relative to the boxsize (and increasing helps the runtime).

  • copy_particles (boolean (default True)) –

    Boolean flag to make a copy of the particle positions If set to False, the particles will be re-ordered in-place

    New in version 2.3.0.

  • enable_min_sep_opt (boolean (default true)) –

    Boolean flag to allow optimizations based on min. separation between pairs of cells. Here to allow for comparison studies.

    New in version 2.3.0.

  • c_api_timer (boolean (default false)) – Boolean flag to measure actual time spent in the C libraries. Here to allow for benchmarking and scaling studies.

  • isa (string (default fastest)) –

    Controls the runtime dispatch for the instruction set to use. Options are: [fastest, avx512f, avx, sse42, fallback]

    Setting isa to fastest will pick the fastest available instruction set on the current computer. However, if you set isa to, say, avx and avx is not available on the computer, then the code will revert to using fallback (even though sse42 might be available). Unless you are benchmarking the different instruction sets, you should always leave isa to the default value. And if you are benchmarking, then the string supplied here gets translated into an enum for the instruction set defined in utils/defs.h.

  • weight_type (str, optional) – The type of pair weighting to apply. Options: “pair_product”, “pair_product_gradient”, None; Default: None.

  • proj_type (string (default None)) –

    Projection method to use; currently supported methods are [‘tophat’, ‘piecewise’, ‘generalr’, ‘gaussian_kernel’, ‘gradient’]. If using ‘gradient’, must set weight_type="pair_product_gradient".

    New in version suave.

  • ncomponents (int (default None)) –

    Number of basis functions; necessary if projection method proj_type is defined.

    New in version suave.

  • projfn (string (default None)) –

    File path of projection file; necessary for proj_type='generalr'.

    New in version suave.

Returns

  • results (A python list) – A python list containing nmu_bins of [smin, smax, savg, mu_max, npairs, weightavg] for each spatial bin specified in the binfile. There will be a total of nmu_bins ranging from [0, mu_max) per spatial bin. If output_savg is not set, then savg will be set to 0.0 for all bins; similarly for weight_avg. npairs contains the number of pairs in that bin.

  • api_time (float, optional) – Only returned if c_api_timer is set. api_time measures only the time spent within the C library and ignores all python overhead.

  • v_proj (array-like, double, optional) – Only returned if proj_type is not None. The projection vector, an array of length ncomponents.

  • t_proj (array-like, double, optional) – Only returned if proj_type is not None. The projection tensor, unrolled in the form of an array with length ncomponents``*``ncomponents.

Example

>>> from __future__ import print_function
>>> import numpy as np
>>> from os.path import dirname, abspath, join as pjoin
>>> import Corrfunc
>>> from Corrfunc.theory.DDsmu import DDsmu
>>> binfile = pjoin(dirname(abspath(Corrfunc.__file__)),
...                 "../theory/tests/", "bins")
>>> N = 10000
>>> boxsize = 420.0
>>> nthreads = 4
>>> autocorr = 1
>>> mu_max = 1.0
>>> seed = 42
>>> nmu_bins = 10
>>> np.random.seed(seed)
>>> X = np.random.uniform(0, boxsize, N)
>>> Y = np.random.uniform(0, boxsize, N)
>>> Z = np.random.uniform(0, boxsize, N)
>>> weights = np.ones_like(X)
>>> results = DDsmu(autocorr, nthreads, binfile, mu_max, nmu_bins,
...                  X, Y, Z, weights1=weights, weight_type='pair_product', output_savg=True)
>>> for r in results[100:]: print("{0:10.6f} {1:10.6f} {2:10.6f} {3:10.1f}"
...                               " {4:10d} {5:10.6f}".format(r['smin'], r['smax'],
...                               r['savg'], r['mu_max'], r['npairs'], r['weightavg']))
...                         
 5.788530   8.249250   7.148213        0.1        230   1.000000
 5.788530   8.249250   7.157218        0.2        236   1.000000
 5.788530   8.249250   7.165338        0.3        208   1.000000
 5.788530   8.249250   7.079905        0.4        252   1.000000
 5.788530   8.249250   7.251661        0.5        184   1.000000
 5.788530   8.249250   7.118536        0.6        222   1.000000
 5.788530   8.249250   7.083466        0.7        238   1.000000
 5.788530   8.249250   7.198184        0.8        170   1.000000
 5.788530   8.249250   7.127409        0.9        208   1.000000
 5.788530   8.249250   6.973090        1.0        206   1.000000
 8.249250  11.756000  10.149183        0.1        592   1.000000
 8.249250  11.756000  10.213009        0.2        634   1.000000
 8.249250  11.756000  10.192220        0.3        532   1.000000
 8.249250  11.756000  10.246931        0.4        544   1.000000
 8.249250  11.756000  10.102675        0.5        530   1.000000
 8.249250  11.756000  10.276180        0.6        644   1.000000
 8.249250  11.756000  10.251264        0.7        666   1.000000
 8.249250  11.756000  10.138399        0.8        680   1.000000
 8.249250  11.756000  10.191916        0.9        566   1.000000
 8.249250  11.756000  10.243229        1.0        608   1.000000
11.756000  16.753600  14.552776        0.1       1734   1.000000
11.756000  16.753600  14.579991        0.2       1806   1.000000
11.756000  16.753600  14.599611        0.3       1802   1.000000
11.756000  16.753600  14.471100        0.4       1820   1.000000
11.756000  16.753600  14.480192        0.5       1740   1.000000
11.756000  16.753600  14.493679        0.6       1746   1.000000
11.756000  16.753600  14.547713        0.7       1722   1.000000
11.756000  16.753600  14.465390        0.8       1750   1.000000
11.756000  16.753600  14.547465        0.9       1798   1.000000
11.756000  16.753600  14.440975        1.0       1828   1.000000
16.753600  23.875500  20.720406        0.1       5094   1.000000
16.753600  23.875500  20.735403        0.2       5004   1.000000
16.753600  23.875500  20.721069        0.3       5172   1.000000
16.753600  23.875500  20.723648        0.4       5014   1.000000
16.753600  23.875500  20.650621        0.5       5094   1.000000
16.753600  23.875500  20.688135        0.6       5076   1.000000
16.753600  23.875500  20.735691        0.7       4910   1.000000
16.753600  23.875500  20.714097        0.8       4864   1.000000
16.753600  23.875500  20.751836        0.9       4954   1.000000
16.753600  23.875500  20.721183        1.0       5070   1.000000

Corrfunc.mocks

Wrapper for all clustering statistic calculations on galaxies in a mock catalog.

Corrfunc.mocks.DDsmu_mocks(autocorr, cosmology, nthreads, mu_max, nmu_bins, binfile, RA1, DEC1, CZ1, weights1=None, RA2=None, DEC2=None, CZ2=None, weights2=None, is_comoving_dist=False, verbose=False, output_savg=False, fast_divide_and_NR_steps=0, xbin_refine_factor=2, ybin_refine_factor=2, zbin_refine_factor=1, max_cells_per_dim=100, copy_particles=True, enable_min_sep_opt=True, c_api_timer=False, isa='fastest', weight_type=None, proj_type=None, ncomponents=None, projfn=None)[source]

Calculate the 2-D pair-counts corresponding to the projected correlation function, \(\xi(s, \mu)\). The pairs are counted in bins of radial separation and cosine of angle to the line-of-sight (LOS). The input positions are expected to be on-sky co-ordinates. This module is suitable for calculating correlation functions for mock catalogs.

If weights are provided, the resulting pair counts are weighted. The weighting scheme depends on weight_type.

Returns a numpy structured array containing the pair counts for the specified bins.

To use the projection capability with suave, set the proj_type parameter for the desired basis functions, and set ncomponents and projfn accordingly.

Note

This module only returns pair counts and not the actual correlation function \(\xi(s, \mu)\). See the utilities Corrfunc.utils.convert_3d_counts_to_cf for computing \(\xi(s, \mu)\) from the pair counts.

New in version 2.1.0.

Parameters
  • autocorr (boolean, required) – Boolean flag for auto/cross-correlation. If autocorr is set to 1, then the second set of particle positions are not required.

  • cosmology (integer, required) –

    Integer choice for setting cosmology. Valid values are 1->LasDamas cosmology and 2->Planck cosmology. If you need arbitrary cosmology, easiest way is to convert the CZ values into co-moving distance, based on your preferred cosmology. Set is_comoving_dist=True, to indicate that the co-moving distance conversion has already been done.

    Choices:
    1. LasDamas cosmology. \(\Omega_m=0.25\), \(\Omega_\Lambda=0.75\)

    2. Planck cosmology. \(\Omega_m=0.302\), \(\Omega_\Lambda=0.698\)

    To setup a new cosmology, add an entry to the function, init_cosmology in ROOT/utils/cosmology_params.c and re-install the entire package.

  • nthreads (integer) – The number of OpenMP threads to use. Has no effect if OpenMP was not enabled during library compilation.

  • mu_max (double. Must be in range [0.0, 1.0]) –

    A double-precision value for the maximum cosine of the angular separation from the line of sight (LOS). Here, mu is defined as the angle between s and l. If \(v_1\) and \(v_2\) represent the vectors to each point constituting the pair, then \(s := v_1 - v_2\) and \(l := 1/2 (v_1 + v_2)\).

    Note: Only pairs with \(0 <= \cos(\theta_{LOS}) < \mu_{max}\) are counted (no equality).

  • nmu_bins (int) – The number of linear mu bins, with the bins ranging from from (0, \(\mu_{max}\))

  • binfile (string or an list/array of floats) –

    For string input: filename specifying the s bins for DDsmu_mocks. The file should contain white-space separated values of (smin, smax) specifying each s bin wanted. The bins need to be contiguous and sorted in increasing order (smallest bins come first).

    For array-like input: A sequence of s values that provides the bin-edges. For example, np.logspace(np.log10(0.1), np.log10(10.0), 15) is a valid input specifying 14 (logarithmic) bins between 0.1 and 10.0. This array does not need to be sorted.

  • RA1 (array-like, real (float/double)) –

    The array of Right Ascensions for the first set of points. RA’s are expected to be in [0.0, 360.0], but the code will try to fix cases where the RA’s are in [-180, 180.0]. For peace of mind, always supply RA’s in [0.0, 360.0].

    Calculations are done in the precision of the supplied arrays.

  • DEC1 (array-like, real (float/double)) –

    Array of Declinations for the first set of points. DEC’s are expected to be in the [-90.0, 90.0], but the code will try to fix cases where the DEC’s are in [0.0, 180.0]. Again, for peace of mind, always supply DEC’s in [-90.0, 90.0].

    Must be of same precision type as RA1.

  • CZ1 (array-like, real (float/double)) –

    Array of (Speed Of Light * Redshift) values for the first set of points. Code will try to detect cases where redshifts have been passed and multiply the entire array with the speed of light.

    If is_comoving_dist is set, then CZ1 is interpreted as the co-moving distance, rather than cz.

  • weights1 (array_like, real (float/double), optional) – A scalar, or an array of weights of shape (n_weights, n_positions) or (n_positions,). weight_type specifies how these weights are used; results are returned in the weightavg field. If only one of weights1 or weights2 is specified, the other will be set to uniform weights.

  • RA2 (array-like, real (float/double)) –

    The array of Right Ascensions for the second set of points. RA’s are expected to be in [0.0, 360.0], but the code will try to fix cases where the RA’s are in [-180, 180.0]. For peace of mind, always supply RA’s in [0.0, 360.0].

    Must be of same precision type as RA1/DEC1/CZ1.

  • DEC2 (array-like, real (float/double)) –

    Array of Declinations for the second set of points. DEC’s are expected to be in the [-90.0, 90.0], but the code will try to fix cases where the DEC’s are in [0.0, 180.0]. Again, for peace of mind, always supply DEC’s in [-90.0, 90.0].

    Must be of same precision type as RA1/DEC1/CZ1.

  • CZ2 (array-like, real (float/double)) –

    Array of (Speed Of Light * Redshift) values for the second set of points. Code will try to detect cases where redshifts have been passed and multiply the entire array with the speed of light.

    If is_comoving_dist is set, then CZ2 is interpreted as the co-moving distance, rather than cz.

    Must be of same precision type as RA1/DEC1/CZ1.

  • weights2 (array-like, real (float/double), optional) – Same as weights1, but for the second set of positions

  • is_comoving_dist (boolean (default false)) – Boolean flag to indicate that cz values have already been converted into co-moving distances. This flag allows arbitrary cosmologies to be used in Corrfunc.

  • verbose (boolean (default false)) – Boolean flag to control output of informational messages

  • output_savg (boolean (default false)) – Boolean flag to output the average s for each bin. Code will run slower if you set this flag. Also, note, if you are calculating in single-precision, savg will suffer from numerical loss of precision and can not be trusted. If you need accurate savg values, then pass in double precision arrays for the particle positions.

  • fast_divide_and_NR_steps (integer (default 0)) – Replaces the division in AVX implementation with an approximate reciprocal, followed by fast_divide_and_NR_steps of Newton-Raphson. Can improve runtime by ~15-20% on older computers. Value of 0 uses the standard division operation.

  • (xyz)bin_refine_factor (integer, default is (2,2,1); typically within [1-3]) – Controls the refinement on the cell sizes. Can have up to a 20% impact on runtime.

  • max_cells_per_dim (integer, default is 100, typical values in [50-300]) – Controls the maximum number of cells per dimension. Total number of cells can be up to (max_cells_per_dim)^3. Only increase if rpmax is too small relative to the boxsize (and increasing helps the runtime).

  • copy_particles (boolean (default True)) –

    Boolean flag to make a copy of the particle positions If set to False, the particles will be re-ordered in-place

    New in version 2.3.0.

  • enable_min_sep_opt (boolean (default true)) –

    Boolean flag to allow optimizations based on min. separation between pairs of cells. Here to allow for comparison studies.

    New in version 2.3.0.

  • c_api_timer (boolean (default false)) – Boolean flag to measure actual time spent in the C libraries. Here to allow for benchmarking and scaling studies.

  • isa (string, case-insensitive (default fastest)) –

    Controls the runtime dispatch for the instruction set to use. Options are: [fastest, avx512f, avx, sse42, fallback]

    Setting isa to fastest will pick the fastest available instruction set on the current computer. However, if you set isa to, say, avx and avx is not available on the computer, then the code will revert to using fallback (even though sse42 might be available).

    Unless you are benchmarking the different instruction sets, you should always leave isa to the default value. And if you are benchmarking, then the string supplied here gets translated into an enum for the instruction set defined in utils/defs.h.

  • weight_type (string, optional (default None)) – The type of weighting to apply. One of [“pair_product”, “pair_product_gradient”, None].

  • proj_type (string (default None)) –

    Projection method to use; currently supported methods are [‘tophat’, ‘piecewise’, ‘generalr’, ‘gaussian_kernel’, ‘gradient’]. If using ‘gradient’, must set weight_type="pair_product_gradient".

    New in version suave.

  • ncomponents (int (default None)) –

    Number of basis functions; necessary if projection method proj_type is defined.

    New in version suave.

  • projfn (string (default None)) –

    Path to projection file; necessary for proj_type='generalr'.

    New in version suave.

Returns

  • results (Numpy structured array) – A numpy structured array containing [smin, smax, savg, mumax, npairs, weightavg]. There are a total of nmu_bins in mu for each separation bin specified in the binfile, with mumax being the upper limit of the mu bin. If output_savg is not set, then savg will be set to 0.0 for all bins; similarly for weightavg. npairs contains the number of pairs in that bin and can be used to compute the actual \(\xi(s, \mu)\) by combining with (DR, RR) counts.

  • api_time (float, optional) – Only returned if c_api_timer is set. api_time measures only the time spent within the C library and ignores all python overhead.

  • v_proj (array-like, double, optional) – Only returned if proj_type is not None. The projection vector, an array of length ncomponents.

  • t_proj (array-like, double, optional) – Only returned if proj_type is not None. The projection tensor, unrolled in the form of an array with length ncomponents``*``ncomponents.

Corrfunc.utils

A set of utility routines

Corrfunc.utils.compute_amps(ncomponents, nd1, nd2, nr1, nr2, dd, dr, rd, rr, trr)[source]

Compute the amplitude vector for the continuous correlation function from the component vectors (a.k.a. pair counts in a given basis).

Parameters
  • ncomponents (int) – Number of basis functions

  • nd1 (int) – Number of particles in data catalog 1

  • nd2 (int) – Number of particles in data catalog 2

  • nr1 (int) – Number of particles in random1

  • nr2 (int) – Number of particles in random2

  • dd (array-like, double) – Projection vector for data-data cross-correlation

  • dr (array-like, double) – Projection vector for data-random cross-correlation

  • rd (array-like, double) – Projection vector for random-data cross-correlation

  • rr (array-like, double) – Projection vector for random-random cross-correlation

  • trr (array-like, double) – Projection tensor for random-random cross-correlation

Returns

amps – Vector of amplitudes, with length ncomponents

Return type

array-like, double

Corrfunc.utils.evaluate_xi(amps, rvals, proj_type, rbins=None, projfn=None, weights1=None, weights2=None, weight_type=None)[source]

Evaluate the correlation function for the given amplitudes and separation values.

Parameters
  • amps (array-like, double) – Vector of amplitudes, e.g. that returned by compute_amps

  • rvals (array-like, double) – Array of radial separation values at which to evaluate the correlation function

  • proj_type (string) – Projection method to use; currently supported methods are [‘tophat’, ‘piecewise’, ‘generalr’, ‘gaussian_kernel’, ‘gradient’]

  • rbins (array-like, double, default=None) – Bin edges; necessary for tophat or piecewise bases

  • projfn (string, default=None) – Path to projection file; necessary for proj_type=’generalr’

  • weights1 (array-like, double, default=None) – Weights/metadata for custom basis functions, for first set of imaginary galaxies

  • weights2 (array-like, double, default=None) – Weights/metadata for custom basis functions, for second set of imaginary galaxies

  • weight_type (string, default=None) – Name of weight function. If using ‘gradient’, must set weight_type="pair_product_gradient".

Returns

xi – Vector of xi values, same shape as rvals

Return type

array-like, double

Corrfunc.utils.trr_analytic(rmin, rmax, nd, volume, ncomponents, proj_type, rbins=None, projfn=None)[source]

Compute the T_RR tensor analytically for a periodic box, for the given set of basis functions.

Parameters
  • rmin (double) – Minimum r-value for integration

  • rmax (double) – Maximum r-value for integration

  • nd (int) – Number of points in data catalog

  • volume (double) – Volume of data cube, in the same units as the r-values (cubed)

  • ncomponents (int) – Number of basis functions

  • proj_type (string) – Projection method to use; currently supported methods are [‘tophat’, ‘piecewise’, ‘generalr’, ‘gaussian_kernel’, ‘gradient’]

  • rbins (array-like, double, default=None) – Edges of r-bins for tophat or piecewise projections

  • projfn (string, default=None) – Path to projection file if necessary

Returns

trr – Tensor of T_RR values, with size (ncomponents, ncomponents)

Return type

array-like, double

Corrfunc.bases

Corrfunc.bases.bao_bases(rmin, rmax, projfn, cosmo_base=None, ncont=2000, redshift=0.0, alpha_guess=1.0, dalpha=0.001, bias=1.0, k0=0.1, k1=10.0, k2=0.1, k3=0.001)[source]

Compute the 5-component BAO basis functions based on a cosmological model and linearized around the scale dilation parameter alpha.

Parameters
  • rmin (double) – Minimum r-value for basis functions

  • rmax (double) – Maximum r-value for basis functions

  • projfn (string, default=None) – Path to projection file if necessary

  • cosmo_base (nbodykit cosmology object, default=nbodykit.cosmology.Planck15) – Cosmology object for the BAO model.

  • ncont (int, default=2000) – Number of continuous r-values at which to write the basis function file

  • redshift (double, default=0.0) – Redshift at which to compute power spectrum

  • alpha_guess (double, default=1.0) – The alpha (scale dilation parameter) at which to compute the model (alpha=1.0 is no scale shift)

  • dalpha (double, default=0.001) – The change in alpha (scale dilation parameter) used to calculate the numerical partial derivative

  • bias (double, default=1.0) – The bias parameter by which to scale the model amplitude (bias=1.0 indicates no bias)

  • k0 (double, default=0.1) – The initial magnitude of the derivative term

  • k1 (double, default=1.0) – The initial magnitude of the s^2 nuisance parameter term

  • k2 (double, default=0.1) – The initial magnitude of the s nuisance parameter term

  • k3 (double, default=0.001) – The initial magnitude of the constant nuisance parameter term

Returns

bases – 2-d array of basis function values; first column is r-values

Return type

array-like, double

Corrfunc.bases.spline_bases(rmin, rmax, projfn, ncomponents, ncont=2000, order=3)[source]

Compute a set of spline basis functions for the given order.

Parameters
  • rmin (double) – Minimum r-value for basis functions

  • rmax (double) – Maximum r-value for basis functions

  • projfn (string, default=None) – Path to projection file if necessary

  • ncomponents (int) – Number of components (basis functions)

  • ncont (int, default=2000) – Number of continuous r-values at which to write the basis function file

  • order (int, default=3) – Order of spline to use; default is cubic spline

Returns

bases – 2-d array of basis function values; first column is r-values

Return type

array-like, double