pycbc.tmpltbank package¶
Submodules¶
pycbc.tmpltbank.bank_conversions module¶
This module is supplied to make a convenience function for converting into specific values from PyCBC template banks.
- pycbc.tmpltbank.bank_conversions.get_bank_property(parameter, bank, template_ids, duration_approximant='SEOBNRv4')[source]¶
Get a specific value from a hdf file object in standard PyCBC template bank format
- Parameters
parameter (str) – the parameter to convert to, must be in conversion_options
bank (h5py File object or dictionary of arrays) – Template bank containing the parameters for use in conversions must contain mass1, mass2, spin1z, spin2z as a minimum
template_ids (numpy array) – Array of template IDs for reading a set of templates from the bank
Parameters (Optional) –
------------------- –
duration_approximant (string) – The approximant used to calculate the duration of the template if not already given
- Returns
values – Array of whatever the requested parameter is calculated for the specified templates in the bank
- Return type
numpy array, same size as template_ids
pycbc.tmpltbank.bank_output_utils module¶
- pycbc.tmpltbank.bank_output_utils.calculate_ethinca_metric_comps(metricParams, ethincaParams, mass1, mass2, spin1z=0.0, spin2z=0.0, full_ethinca=True)[source]¶
Calculate the Gamma components needed to use the ethinca metric. At present this outputs the standard TaylorF2 metric over the end time and chirp times au_0 and au_3. A desirable upgrade might be to use the chi coordinates [defined WHERE?] for metric distance instead of au_0 and au_3. The lower frequency cutoff is currently hard-coded to be the same as the bank layout options fLow and f0 (which must be the same as each other).
- Parameters
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
ethincaParams (ethincaParameters instance) – Structure holding options relevant to the ethinca metric computation.
mass1 (float) – Mass of the heavier body in the considered template.
mass2 (float) – Mass of the lighter body in the considered template.
spin1z (float (optional, default=0)) – Spin of the heavier body in the considered template.
spin2z (float (optional, default=0)) – Spin of the lighter body in the considered template.
full_ethinca (boolean (optional, default=True)) – If True calculate the ethinca components in all 3 directions (mass1, mass2 and time). If False calculate only the time component (which is stored in Gamma0).
- Returns
fMax_theor (float) – Value of the upper frequency cutoff given by the template parameters and the cutoff formula requested.
gammaVals (numpy_array) – Array holding 6 independent metric components in (end_time, tau_0, tau_3) coordinates to be stored in the Gamma0-5 slots of a SnglInspiral object.
- pycbc.tmpltbank.bank_output_utils.convert_to_sngl_inspiral_table(params, proc_id)[source]¶
Convert a list of m1,m2,spin1z,spin2z values into a basic sngl_inspiral table with mass and spin parameters populated and event IDs assigned
- Parameters
params (iterable) – Each entry in the params iterable should be a sequence of [mass1, mass2, spin1z, spin2z] in that order
proc_id (int) – Process ID to add to each row of the sngl_inspiral table
- Returns
Bank of templates in SnglInspiralTable format
- Return type
SnglInspiralTable
- pycbc.tmpltbank.bank_output_utils.output_sngl_inspiral_table(outputFile, tempBank, metricParams, ethincaParams, programName='', optDict=None, outdoc=None, **kwargs)[source]¶
Function that converts the information produced by the various pyCBC bank generation codes into a valid LIGOLW xml file containing a sngl_inspiral table and outputs to file.
- Parameters
outputFile (string) – Name of the file that the bank will be written to
tempBank (iterable) – Each entry in the tempBank iterable should be a sequence of [mass1,mass2,spin1z,spin2z] in that order.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
ethincaParams ({ethincaParameters instance, None}) – Structure holding options relevant to the ethinca metric computation including the upper frequency cutoff to be used for filtering. NOTE: The computation is currently only valid for non-spinning systems and uses the TaylorF2 approximant.
(key-word-argument) (programName) – Name of the executable that has been run
argument) (outdoc (key-word) – Dictionary of the command line arguments passed to the program
argument) – If given add template bank to this representation of a xml document and write to disk. If not given create a new document.
kwargs (key-word arguments) – All other key word arguments will be passed directly to ligolw_process.register_to_xmldoc
pycbc.tmpltbank.brute_force_methods module¶
- pycbc.tmpltbank.brute_force_methods.find_xi_extrema_brute(xis, bestMasses, bestXis, direction_num, req_match, massRangeParams, metricParams, fUpper, find_minimum=False, scaleFactor=0.8, numIterations=3000)[source]¶
This function is used to find the largest or smallest value of the xi space in a specified dimension at a specified point in the higher dimensions. It does this by iteratively throwing points at the space to find extrema.
- Parameters
xis (list or array) – Position in the xi space at which to assess the depth. This can be only a subset of the higher dimensions than that being sampled.
bestMasses (list) – Contains [totalMass, eta, spin1z, spin2z]. Is a physical position mapped to xi coordinates in bestXis that is close to the xis point. This is aimed to give the code a starting point.
bestXis (list) – Contains the position of bestMasses in the xi coordinate system.
direction_num (int) – The dimension that you want to assess the depth of (0 = 1, 1 = 2 …)
req_match (float) – When considering points to assess the depth with, only consider points with a mismatch that is smaller than this with xis.
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper that was used when obtaining the xi_i coordinates. This lets us know how to rotate potential physical points into the correct xi_i space. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
find_minimum (boolean, optional (default = False)) – If True, find the minimum value of the xi direction. If False find the maximum value.
scaleFactor (float, optional (default = 0.8)) – The value of the scale factor that is used when calling pycbc.tmpltbank.get_mass_distribution.
numIterations (int, optional (default = 3000)) – The number of times to make calls to get_mass_distribution when assessing the maximum/minimum of this parameter space. Making this smaller makes the code faster, but at the cost of accuracy.
- Returns
xi_extent – The extremal value of the specified dimension at the specified point in parameter space.
- Return type
- pycbc.tmpltbank.brute_force_methods.get_mass_distribution(bestMasses, scaleFactor, massRangeParams, metricParams, fUpper, numJumpPoints=100, chirpMassJumpFac=0.0001, etaJumpFac=0.01, spin1zJumpFac=0.01, spin2zJumpFac=0.01)[source]¶
Given a set of masses, this function will create a set of points nearby in the mass space and map these to the xi space.
- Parameters
bestMasses (list) – Contains [ChirpMass, eta, spin1z, spin2z]. Points will be placed around tjos
scaleFactor (float) – This parameter describes the radius away from bestMasses that points will be placed in.
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper that was used when obtaining the xi_i coordinates. This lets us know how to rotate potential physical points into the correct xi_i space. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
numJumpPoints (int, optional (default = 100)) – The number of points that will be generated every iteration
chirpMassJumpFac (float, optional (default=0.0001)) – The jump points will be chosen with fractional variation in chirpMass up to this multiplied by scaleFactor.
etaJumpFac (float, optional (default=0.01)) – The jump points will be chosen with fractional variation in eta up to this multiplied by scaleFactor.
spin1zJumpFac (float, optional (default=0.01)) – The jump points will be chosen with absolute variation in spin1z up to this multiplied by scaleFactor.
spin2zJumpFac (float, optional (default=0.01)) – The jump points will be chosen with absolute variation in spin2z up to this multiplied by scaleFactor.
- Returns
Totmass (numpy.array) – Total mass of the resulting points
Eta (numpy.array) – Symmetric mass ratio of the resulting points
Spin1z (numpy.array) – Spin of the heavier body of the resulting points
Spin2z (numpy.array) – Spin of the smaller body of the resulting points
Diff (numpy.array) – Mass1 - Mass2 of the resulting points
Mass1 (numpy.array) – Mass1 (mass of heavier body) of the resulting points
Mass2 (numpy.array) – Mass2 (mass of smaller body) of the resulting points
new_xis (list of numpy.array) – Position of points in the xi coordinates
- pycbc.tmpltbank.brute_force_methods.get_physical_covaried_masses(xis, bestMasses, bestXis, req_match, massRangeParams, metricParams, fUpper, giveUpThresh=5000)[source]¶
This function takes the position of a point in the xi parameter space and iteratively finds a close point in the physical coordinate space (masses and spins).
- Parameters
xis (list or array) – Desired position of the point in the xi space. If only N values are provided and the xi space’s dimension is larger then it is assumed that any value in the remaining xi coordinates is acceptable.
bestMasses (list) – Contains [totalMass, eta, spin1z, spin2z]. Is a physical position mapped to xi coordinates in bestXis that is close to the desired point. This is aimed to give the code a starting point.
bestXis (list) – Contains the position of bestMasses in the xi coordinate system.
req_match (float) – Desired maximum mismatch between xis and the obtained point. If a point is found with mismatch < req_match immediately stop and return that point. A point with this mismatch will not always be found.
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper that was used when obtaining the xi_i coordinates. This lets us know how to rotate potential physical points into the correct xi_i space. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
giveUpThresh (int, optional (default = 5000)) – The program will try this many iterations. If no close matching point has been found after this it will give up.
- Returns
mass1 (float) – The heavier mass of the obtained point.
mass2 (float) – The smaller mass of the obtained point
spin1z (float) – The heavier bodies spin of the obtained point.
spin2z (float) – The smaller bodies spin of the obtained point.
count (int) – How many iterations it took to find the point. For debugging.
mismatch (float) – The mismatch between the obtained point and the input xis.
new_xis (list) – The position of the point in the xi space
- pycbc.tmpltbank.brute_force_methods.stack_xi_direction_brute(xis, bestMasses, bestXis, direction_num, req_match, massRangeParams, metricParams, fUpper, scaleFactor=0.8, numIterations=3000)[source]¶
This function is used to assess the depth of the xi_space in a specified dimension at a specified point in the higher dimensions. It does this by iteratively throwing points at the space to find maxima and minima.
- Parameters
xis (list or array) – Position in the xi space at which to assess the depth. This can be only a subset of the higher dimensions than that being sampled.
bestMasses (list) – Contains [totalMass, eta, spin1z, spin2z]. Is a physical position mapped to xi coordinates in bestXis that is close to the xis point. This is aimed to give the code a starting point.
bestXis (list) – Contains the position of bestMasses in the xi coordinate system.
direction_num (int) – The dimension that you want to assess the depth of (0 = 1, 1 = 2 …)
req_match (float) – When considering points to assess the depth with, only consider points with a mismatch that is smaller than this with xis.
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper that was used when obtaining the xi_i coordinates. This lets us know how to rotate potential physical points into the correct xi_i space. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
scaleFactor (float, optional (default = 0.8)) – The value of the scale factor that is used when calling pycbc.tmpltbank.get_mass_distribution.
numIterations (int, optional (default = 3000)) – The number of times to make calls to get_mass_distribution when assessing the maximum/minimum of this parameter space. Making this smaller makes the code faster, but at the cost of accuracy.
- Returns
xi_min (float) – The minimal value of the specified dimension at the specified point in parameter space.
xi_max (float) –
- The maximal value of the specified dimension at the specified point in
parameter space.
pycbc.tmpltbank.calc_moments module¶
- pycbc.tmpltbank.calc_moments.calculate_metric(Js, logJs, loglogJs, logloglogJs, loglogloglogJs, mapping)[source]¶
This function will take the various integrals calculated by get_moments and convert this into a metric for the appropriate parameter space.
- Parameters
Js (Dictionary) – The list of (log^0 x) * x**(-i/3) integrals computed by get_moments() The index is Js[i]
logJs (Dictionary) – The list of (log^1 x) * x**(-i/3) integrals computed by get_moments() The index is logJs[i]
loglogJs (Dictionary) – The list of (log^2 x) * x**(-i/3) integrals computed by get_moments() The index is loglogJs[i]
logloglogJs (Dictionary) – The list of (log^3 x) * x**(-i/3) integrals computed by get_moments() The index is logloglogJs[i]
loglogloglogJs (Dictionary) – The list of (log^4 x) * x**(-i/3) integrals computed by get_moments() The index is loglogloglogJs[i]
mapping (dictionary) – Used to identify which Lambda components are active in this parameter space and map these to entries in the metric matrix.
- Returns
metric – The resulting metric.
- Return type
numpy.matrix
- pycbc.tmpltbank.calc_moments.calculate_metric_comp(gs, unmax_metric, i, j, Js, logJs, loglogJs, logloglogJs, loglogloglogJs, mapping)[source]¶
Used to compute part of the metric. Only call this from within calculate_metric(). Please see the documentation for that function.
- pycbc.tmpltbank.calc_moments.calculate_moment(psd_f, psd_amp, fmin, fmax, f0, funct, norm=None, vary_fmax=False, vary_density=None)[source]¶
Function for calculating one of the integrals used to construct a template bank placement metric. The integral calculated will be
int funct(x) * (psd_x)**(-7./3.) * delta_x / PSD(x)
where x = f / f0. The lower frequency cutoff is given by fmin, see the parameters below for details on how the upper frequency cutoff is chosen
- Parameters
psd_f (numpy.array) – numpy array holding the set of evenly spaced frequencies used in the PSD
psd_amp (numpy.array) – numpy array holding the PSD values corresponding to the psd_f frequencies
fmin (float) – The lower frequency cutoff used in the calculation of the integrals used to obtain the metric.
fmax (float) – The upper frequency cutoff used in the calculation of the integrals used to obtain the metric. This can be varied (see the vary_fmax option below).
f0 (float) – This is an arbitrary scaling factor introduced to avoid the potential for numerical overflow when calculating this. Generally the default value (70) is safe here. IMPORTANT, if you want to calculate the ethinca metric components later this MUST be set equal to f_low.
funct (Lambda function) – The function to use when computing the integral as described above.
norm (Dictionary of floats) – If given then moment[f_cutoff] will be divided by norm[f_cutoff]
vary_fmax (boolean, optional (default False)) – If set to False the metric and rotations are calculated once, for the full range of frequency [f_low,f_upper). If set to True the metric and rotations are calculated multiple times, for frequency ranges [f_low,f_low + i*vary_density), where i starts at 1 and runs up until f_low + (i+1)*vary_density > f_upper. Thus values greater than f_upper are not computed. The calculation for the full range [f_low,f_upper) is also done.
vary_density (float, optional) – If vary_fmax is True, this will be used in computing the frequency ranges as described for vary_fmax.
- Returns
moment – moment[f_cutoff] will store the value of the moment at the frequency cutoff given by f_cutoff.
- Return type
Dictionary of floats
- pycbc.tmpltbank.calc_moments.determine_eigen_directions(metricParams, preserveMoments=False, vary_fmax=False, vary_density=None)[source]¶
This function will calculate the coordinate transfomations that are needed to rotate from a coordinate system described by the various Lambda components in the frequency expansion, to a coordinate system where the metric is Cartesian.
- Parameters
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric.
preserveMoments (boolean, optional (default False)) – Currently only used for debugging. If this is given then if the moments structure is already set within metricParams then they will not be recalculated.
vary_fmax (boolean, optional (default False)) – If set to False the metric and rotations are calculated once, for the full range of frequency [f_low,f_upper). If set to True the metric and rotations are calculated multiple times, for frequency ranges [f_low,f_low + i*vary_density), where i starts at 1 and runs up until f_low + (i+1)*vary_density > f_upper. Thus values greater than f_upper are not computed. The calculation for the full range [f_low,f_upper) is also done.
vary_density (float, optional) – If vary_fmax is True, this will be used in computing the frequency ranges as described for vary_fmax.
- Returns
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric. THIS FUNCTION ONLY RETURNS THE CLASS The following will be added to this structure
metricParams.evals (Dictionary of numpy.array) – Each entry in the dictionary corresponds to the different frequency ranges described in vary_fmax. If vary_fmax = False, the only entry will be f_upper, this corresponds to integrals in [f_low,f_upper). This entry is always present. Each other entry will use floats as keys to the dictionary. These floats give the upper frequency cutoff when it is varying. Each numpy.array contains the eigenvalues which, with the eigenvectors in evecs, are needed to rotate the coordinate system to one in which the metric is the identity matrix.
metricParams.evecs (Dictionary of numpy.matrix) – Each entry in the dictionary is as described under evals. Each numpy.matrix contains the eigenvectors which, with the eigenvalues in evals, are needed to rotate the coordinate system to one in which the metric is the identity matrix.
metricParams.metric (Dictionary of numpy.matrix) – Each entry in the dictionary is as described under evals. Each numpy.matrix contains the metric of the parameter space in the Lambda_i coordinate system.
metricParams.moments (Moments structure) – See the structure documentation for a description of this. This contains the result of all the integrals used in computing the metrics above. It can be used for the ethinca components calculation, or other similar calculations.
- pycbc.tmpltbank.calc_moments.get_moments(metricParams, vary_fmax=False, vary_density=None)[source]¶
This function will calculate the various integrals (moments) that are needed to compute the metric used in template bank placement and coincidence.
- Parameters
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric.
vary_fmax (boolean, optional (default False)) – If set to False the metric and rotations are calculated once, for the full range of frequency [f_low,f_upper). If set to True the metric and rotations are calculated multiple times, for frequency ranges [f_low,f_low + i*vary_density), where i starts at 1 and runs up until f_low + (i+1)*vary_density > f_upper. Thus values greater than f_upper are not computed. The calculation for the full range [f_low,f_upper) is also done.
vary_density (float, optional) – If vary_fmax is True, this will be used in computing the frequency ranges as described for vary_fmax.
- Returns
None (None) – THIS FUNCTION RETURNS NOTHING The following will be added to the metricParams structure
metricParams.moments (Moments structure) – This contains the result of all the integrals used in computing the metrics above. It can be used for the ethinca components calculation, or other similar calculations. This is composed of two compound dictionaries. The first entry indicates which moment is being calculated and the second entry indicates the upper frequency cutoff that was used.
In all cases x = f/f0.
For the first entries the options are:
moments[‘J%d’ %(i)][f_cutoff] This stores the integral of x**((-i)/3.) * delta X / PSD(x)
moments[‘log%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.))) x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**2 x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**3 x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**4 x**((-i)/3.) * delta X / PSD(x)
The second entry stores the frequency cutoff used when computing the integral. See description of the vary_fmax option above.
All of these values are nomralized by a factor of
x**((-7)/3.) * delta X / PSD(x)
The normalization factor can be obtained in
moments[‘I7’][f_cutoff]
- pycbc.tmpltbank.calc_moments.interpolate_psd(psd_f, psd_amp, deltaF)[source]¶
Function to interpolate a PSD to a different value of deltaF. Uses linear interpolation.
- Parameters
- Returns
new_psd_f (numpy.array) – Array of the frequencies contained within the interpolated PSD
new_psd_amp (numpy.array) – Array of the interpolated PSD values at the frequencies in new_psd_f.
pycbc.tmpltbank.coord_utils module¶
- pycbc.tmpltbank.coord_utils.calc_point_dist(vsA, entryA)[source]¶
This function is used to determine the distance between two points.
- Parameters
- Returns
val – The metric distance between the two points.
- Return type
- pycbc.tmpltbank.coord_utils.calc_point_dist_vary(mus1, fUpper1, mus2, fUpper2, fMap, norm_map, MMdistA)[source]¶
Function to determine if two points, with differing upper frequency cutoffs have a mismatch < MMdistA for both upper frequency cutoffs.
- Parameters
mus1 (List of numpy arrays) – mus1[i] will give the array of point 1’s position in the chi_j coordinate system. The i element corresponds to varying values of the upper frequency cutoff. fMap is used to map between i and actual frequencies
fUpper1 (float) – The upper frequency cutoff of point 1.
mus2 (List of numpy arrays) – mus2[i] will give the array of point 2’s position in the chi_j coordinate system. The i element corresponds to varying values of the upper frequency cutoff. fMap is used to map between i and actual frequencies
fUpper2 (float) – The upper frequency cutoff of point 2.
fMap (dictionary) – fMap[fUpper] will give the index needed to get the chi_j coordinates in the two sets of mus
norm_map (dictionary) – norm_map[fUpper] will give the relative frequency domain template amplitude (sigma) at the given value of fUpper.
MMdistA – The minimal mismatch allowed between the points
- Returns
True if the points have a mismatch < MMdistA False if the points have a mismatch > MMdistA
- Return type
Boolean
- pycbc.tmpltbank.coord_utils.estimate_mass_range(numPoints, massRangeParams, metricParams, fUpper, covary=True)[source]¶
This function will generate a large set of points with random masses and spins (using pycbc.tmpltbank.get_random_mass) and translate these points into the xi_i coordinate system for the given upper frequency cutoff.
- Parameters
numPoints (int) – Number of systems to simulate
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper to use when getting the mu coordinates from the lambda coordinates. This must be a key in metricParams.evals and metricParams.evecs (ie. we must know how to do the transformation for the given value of fUpper). It also must be a key in metricParams.evecsCV if covary=True.
covary (boolean, optional (default = True)) – If this is given then evecsCV will be used to rotate from the Cartesian coordinate system into the principal coordinate direction (xi_i). If not given then points in the original Cartesian coordinates are returned.
- Returns
xis – A list of the positions of each point in the xi_i coordinate system.
- Return type
numpy.array
- pycbc.tmpltbank.coord_utils.find_closest_calculated_frequencies(input_freqs, metric_freqs)[source]¶
Given a value (or array) of input frequencies find the closest values in the list of frequencies calculated in the metric.
- Parameters
input_freqs (numpy.array or float) – The frequency(ies) that you want to find the closest value in metric_freqs
metric_freqs (numpy.array) – The list of frequencies calculated by the metric
- Returns
output_freqs – The list of closest values to input_freqs for which the metric was computed
- Return type
numpy.array or float
- pycbc.tmpltbank.coord_utils.find_max_and_min_frequencies(name, mass_range_params, freqs)[source]¶
ADD DOCS
- pycbc.tmpltbank.coord_utils.get_conv_params(mass1, mass2, spin1z, spin2z, metricParams, fUpper, lambda1=None, lambda2=None, quadparam1=None, quadparam2=None)[source]¶
Function to convert between masses and spins and locations in the mu parameter space. Mu = Cartesian metric, but not principal components.
- Parameters
mass1 (float) – Mass of heavier body.
mass2 (float) – Mass of lighter body.
spin1z (float) – Spin of body 1.
spin2z (float) – Spin of body 2.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper to use when getting the mu coordinates from the lambda coordinates. This must be a key in metricParams.evals and metricParams.evecs (ie. we must know how to do the transformation for the given value of fUpper)
- Returns
mus – Position of the system(s) in the mu coordinate system
- Return type
list of floats or numpy.arrays
- pycbc.tmpltbank.coord_utils.get_cov_params(mass1, mass2, spin1z, spin2z, metricParams, fUpper, lambda1=None, lambda2=None, quadparam1=None, quadparam2=None)[source]¶
Function to convert between masses and spins and locations in the xi parameter space. Xi = Cartesian metric and rotated to principal components.
- Parameters
mass1 (float) – Mass of heavier body.
mass2 (float) – Mass of lighter body.
spin1z (float) – Spin of body 1.
spin2z (float) – Spin of body 2.
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper to use when getting the mu coordinates from the lambda coordinates. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
- Returns
xis – Position of the system(s) in the xi coordinate system
- Return type
list of floats or numpy.arrays
- pycbc.tmpltbank.coord_utils.get_covaried_params(mus, evecsCV)[source]¶
Function to rotate from position(s) in the mu_i coordinate system into the position(s) in the xi_i coordinate system
- Parameters
mus (list of floats or numpy.arrays) – Position of the system(s) in the mu coordinate system
evecsCV (numpy.matrix) – This matrix is used to perform the rotation to the xi_i coordinate system.
- Returns
xis – Position of the system(s) in the xi coordinate system
- Return type
list of floats or numpy.arrays
- pycbc.tmpltbank.coord_utils.get_mu_params(lambdas, metricParams, fUpper)[source]¶
Function to rotate from the lambda coefficients into position in the mu coordinate system. Mu = Cartesian metric, but not principal components.
- Parameters
lambdas (list of floats or numpy.arrays) – Position of the system(s) in the lambda coefficients
metricParams (metricParameters instance) – Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
fUpper (float) – The value of fUpper to use when getting the mu coordinates from the lambda coordinates. This must be a key in metricParams.evals and metricParams.evecs (ie. we must know how to do the transformation for the given value of fUpper)
- Returns
mus – Position of the system(s) in the mu coordinate system
- Return type
list of floats or numpy.arrays
- pycbc.tmpltbank.coord_utils.get_point_distance(point1, point2, metricParams, fUpper)[source]¶
Function to calculate the mismatch between two points, supplied in terms of the masses and spins, using the xi_i parameter space metric to approximate the mismatch of the two points. Can also take one of the points as an array of points and return an array of mismatches (but only one can be an array!)
- point1List of floats or numpy.arrays
point1[0] contains the mass(es) of the heaviest body(ies). point1[1] contains the mass(es) of the smallest body(ies). point1[2] contains the spin(es) of the heaviest body(ies). point1[3] contains the spin(es) of the smallest body(ies).
- point2List of floats
point2[0] contains the mass of the heaviest body. point2[1] contains the mass of the smallest body. point2[2] contains the spin of the heaviest body. point2[3] contains the spin of the smallest body.
- metricParamsmetricParameters instance
Structure holding all the options for construction of the metric and the eigenvalues, eigenvectors and covariance matrix needed to manipulate the space.
- fUpperfloat
The value of fUpper to use when getting the mu coordinates from the lambda coordinates. This must be a key in metricParams.evals, metricParams.evecs and metricParams.evecsCV (ie. we must know how to do the transformation for the given value of fUpper)
- Returns
dist (float or numpy.array) – Distance between the point2 and all points in point1
xis1 (List of floats or numpy.arrays) – Position of the input point1(s) in the xi_i parameter space
xis2 (List of floats) – Position of the input point2 in the xi_i parameter space
- pycbc.tmpltbank.coord_utils.get_random_mass(numPoints, massRangeParams)[source]¶
This function will generate a large set of points within the chosen mass and spin space, and with the desired minimum remnant disk mass (this applies to NS-BH systems only). It will also return the corresponding PN spin coefficients for ease of use later (though these may be removed at some future point).
- Parameters
numPoints (int) – Number of systems to simulate
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
- Returns
mass1 (float) – Mass of heavier body.
mass2 (float) – Mass of lighter body.
spin1z (float) – Spin of body 1.
spin2z (float) – Spin of body 2.
- pycbc.tmpltbank.coord_utils.get_random_mass_point_particles(numPoints, massRangeParams)[source]¶
This function will generate a large set of points within the chosen mass and spin space. It will also return the corresponding PN spin coefficients for ease of use later (though these may be removed at some future point).
- Parameters
numPoints (int) – Number of systems to simulate
massRangeParams (massRangeParameters instance) – Instance holding all the details of mass ranges and spin ranges.
- Returns
mass1 (float) – Mass of heavier body.
mass2 (float) – Mass of lighter body.
spin1z (float) – Spin of body 1.
spin2z (float) – Spin of body 2.
- pycbc.tmpltbank.coord_utils.outspiral_loop(N)[source]¶
Return a list of points that will loop outwards in a 2D lattice in terms of distance from a central point. So if N=2 this will be [0,0], [0,1], [0,-1],[1,0],[-1,0],[1,1] …. This is useful when you want to loop over a number of bins, but want to start in the center and work outwards.
- pycbc.tmpltbank.coord_utils.return_nearest_cutoff(name, mass_dict, freqs)[source]¶
Given an array of total mass values and an (ascending) list of frequencies, this will calculate the specified cutoff formula for each mtotal and return the nearest frequency to each cutoff from the input list. Currently only supports cutoffs that are functions of the total mass and no other parameters (SchwarzISCO, LightRing, ERD)
- Parameters
name (string) – Name of the cutoff formula to be approximated
mass_dict (Dictionary where the keys are used to call the functions) – returned by tmpltbank.named_frequency_cutoffs. The values can be numpy arrays or single values.
freqs (list of floats) – A list of frequencies (must be sorted ascending)
- Returns
The frequencies closest to the cutoff for each value of totmass.
- Return type
numpy.array
- pycbc.tmpltbank.coord_utils.rotate_vector(evecs, old_vector, rescale_factor, index)[source]¶
Function to find the position of the system(s) in one of the xi_i or mu_i directions.
- Parameters
evecs (numpy.matrix) – Matrix of the eigenvectors of the metric in lambda_i coordinates. Used to rotate to a Cartesian coordinate system.
old_vector (list of floats or numpy.arrays) – The position of the system(s) in the original coordinates
rescale_factor (float) – Scaling factor to apply to resulting position(s)
index (int) – The index of the final coordinate system that is being computed. Ie. if we are going from mu_i -> xi_j, this will give j.
- Returns
positions – Position of the point(s) in the resulting coordinate.
- Return type
float or numpy.array
- pycbc.tmpltbank.coord_utils.test_point_dist(point_1_chis, point_2_chis, distance_threshold)[source]¶
This function tests if the difference between two points in the chi parameter space is less than a distance threshold. Returns True if it is and False if it is not.
- Parameters
point_1_chis (numpy.array) – An array of point 1’s position in the chi_i coordinate system
point_2_chis (numpy.array) – An array of point 2’s position in the chi_i coordinate system
distance_threshold (float) – The distance threshold to use.
pycbc.tmpltbank.em_progenitors module¶
- pycbc.tmpltbank.em_progenitors.ISCO_eq(r, chi)[source]¶
Polynomial that enables the calculation of the Kerr inntermost stable circular orbit (ISCO) radius via its roots.
- pycbc.tmpltbank.em_progenitors.ISSO_eq_at_pole(r, chi)[source]¶
Polynomial that enables the calculation of the Kerr polar (inclination = +/- pi/2) innermost stable spherical orbit (ISSO) radius via its roots. Physical solutions are between 6 and 1+sqrt[3]+sqrt[3+2sqrt[3]].
- pycbc.tmpltbank.em_progenitors.PG_ISSO_eq(r, chi, incl)[source]¶
Polynomial that enables the calculation of a generic innermost stable spherical orbit (ISSO) radius via its roots. Physical solutions are between the equatorial ISSO (aka the ISCO) radius and the polar ISSO radius. See Stone, Loeb, Berger, PRD 87, 084053 (2013).
- Parameters
- Returns
r**8*Z+chi**2*(1-cos_incl**2)*(chi**2*(1-cos_incl**2)*Y-2*r**4*X)
whereX=chi**2*(chi**2*(3*chi**2+4*r*(2*r-3))+r**2*(15*r*(r-4)+28))-6*r**4*(r**2-4)
Y=chi**4*(chi**4+r**2*(7*r*(3*r-4)+36))+6*r*(r-2)*(chi**6+2*r**3*(chi**2*(3*r+2)+3*r**2*(r-2)))
Z=ISCO_eq(r,chi)
- Return type
- pycbc.tmpltbank.em_progenitors.PG_ISSO_solver(chi, incl)[source]¶
Function that determines the radius of the innermost stable spherical orbit (ISSO) for a Kerr BH and a generic inclination angle between the BH spin and the orbital angular momentum. This function finds the appropriat root of PG_ISSO_eq.
- pycbc.tmpltbank.em_progenitors.load_ns_sequence(eos_name)[source]¶
Load the data of an NS non-rotating equilibrium sequence generated using the equation of state (EOS) chosen by the user. [Only the 2H 2-piecewise polytropic EOS is currently supported. This yields NSs with large radiss (15-16km).]
- Parameters
eos_name (string) – NS equation of state label (‘2H’ is the only supported choice at the moment)
- Returns
ns_sequence (3D-array) –
- contains the sequence data in the form NS gravitational
mass (in solar masses), NS baryonic mass (in solar masses), NS compactness (dimensionless)
max_ns_g_mass (float) – the maximum NS gravitational mass (in solar masses) in the sequence (this is the mass of the most massive stable NS)
- pycbc.tmpltbank.em_progenitors.ns_g_mass_to_ns_b_mass(ns_g_mass, ns_sequence)[source]¶
Determines the baryonic mass of an NS given its gravitational mass and an NS equilibrium sequence.
- Parameters
ns_g_mass (float) – NS gravitational mass (in solar masses)
ns_sequence (3D-array) –
- contains the sequence data in the form NS gravitational
mass (in solar masses), NS baryonic mass (in solar masses), NS compactness (dimensionless)
- Returns
The NS baryonic mass (in solar massesr**3*(r**2*(r-6)+chi**2*(3*r+4))+ chi**4*(3*r*(r-2)+chi**2))
- Return type
- pycbc.tmpltbank.em_progenitors.ns_g_mass_to_ns_compactness(ns_g_mass, ns_sequence)[source]¶
Determines the compactness of an NS given its gravitational mass and an NS equilibrium sequence.
- Parameters
ns_g_mass (float) – NS gravitational mass (in solar masses)
ns_sequence (3D-array) –
- contains the sequence data in the form NS gravitational
mass (in solar masses), NS baryonic mass (in solar masses), NS compactness (dimensionless)
- Returns
The NS compactness (dimensionless)
- Return type
- pycbc.tmpltbank.em_progenitors.remnant_mass(eta, ns_g_mass, ns_sequence, chi, incl)[source]¶
Function that determines the remnant disk mass of an NS-BH system using the fit to numerical-relativity results discussed in Foucart, Hinderer, Nissanke PRD 98, 081501(R) (2018).
- Parameters
eta (float) – the symmetric mass ratio of the binary
ns_g_mass (float) – NS gravitational mass (in solar masses)
ns_sequence (3D-array) –
- contains the sequence data in the form NS gravitational
mass (in solar masses), NS baryonic mass (in solar masses), NS compactness (dimensionless)
chi (float) – the BH dimensionless spin parameter
incl (float) – the inclination angle between the BH spin and the orbital angular momentum in radians
- Returns
remnant_mass – The remnant mass in solar masses
- Return type
pycbc.tmpltbank.lambda_mapping module¶
- pycbc.tmpltbank.lambda_mapping.ethinca_order_from_string(order)[source]¶
Returns the integer giving twice the post-Newtonian order used by the ethinca calculation. Currently valid only for TaylorF2 metric
- Parameters
order (string) –
- Return type
- pycbc.tmpltbank.lambda_mapping.generate_inverse_mapping(order)[source]¶
Genereate a lambda entry -> PN order map.
This function will generate the opposite of generate mapping. So where generate_mapping gives dict[key] = item this will give dict[item] = key. Valid PN orders are:
zeroPN: Will only include the dominant term (proportional to chirp mass)
onePN: Will only the leading orbit term and first correction at 1PN
onePointFivePN: Will include orbit and spin terms to 1.5PN.
twoPN: Will include orbit and spin terms to 2PN.
twoPointFivePN: Will include orbit and spin terms to 2.5PN.
threePN: Will include orbit terms to 3PN and spin terms to 2.5PN.
threePointFivePN: Include orbit terms to 3.5PN and spin terms to 2.5PN
- Parameters
order (string) – A string containing a PN order. Valid values are given above.
- Returns
mapping – An inverse mapping between the active Lambda terms and index in the metric
- Return type
dictionary
- pycbc.tmpltbank.lambda_mapping.generate_mapping(order)[source]¶
This function will take an order string and return a mapping between components in the metric and the various Lambda components. This must be used (and consistently used) when generating the metric and when transforming to/from the xi_i coordinates to the lambda_i coordinates.
NOTE: This is not a great way of doing this. It would be nice to clean this up. Hence pulling this function out. The valid PN orders are
zeroPN: Will only include the dominant term (proportional to chirp mass)
onePN: Will only the leading orbit term and first correction at 1PN
onePointFivePN: Will include orbit and spin terms to 1.5PN.
twoPN: Will include orbit and spin terms to 2PN.
twoPointFivePN: Will include orbit and spin terms to 2.5PN.
threePN: Will include orbit terms to 3PN and spin terms to 2.5PN.
threePointFivePN: Include orbit terms to 3.5PN and spin terms to 2.5PN
- Parameters
order (string) – A string containing a PN order. Valid values are given above.
- Returns
mapping – A mapping between the active Lambda terms and index in the metric
- Return type
dictionary
- pycbc.tmpltbank.lambda_mapping.get_chirp_params(mass1, mass2, spin1z, spin2z, f0, order, quadparam1=None, quadparam2=None, lambda1=None, lambda2=None)[source]¶
Take a set of masses and spins and convert to the various lambda coordinates that describe the orbital phase. Accepted PN orders are:
zeroPN: Will only include the dominant term (proportional to chirp mass)
onePN: Will only the leading orbit term and first correction at 1PN
onePointFivePN: Will include orbit and spin terms to 1.5PN.
twoPN: Will include orbit and spin terms to 2PN.
twoPointFivePN: Will include orbit and spin terms to 2.5PN.
threePN: Will include orbit terms to 3PN and spin terms to 2.5PN.
threePointFivePN: Include orbit terms to 3.5PN and spin terms to 2.5PN
- Parameters
mass1 (float or array) – Mass1 of input(s).
mass2 (float or array) – Mass2 of input(s).
spin1z (float or array) – Parallel spin component(s) of body 1.
spin2z (float or array) – Parallel spin component(s) of body 2.
f0 (float) – This is an arbitrary scaling factor introduced to avoid the potential for numerical overflow when calculating this. Generally the default value (70) is safe here. IMPORTANT, if you want to calculate the ethinca metric components later this MUST be set equal to f_low. This value must also be used consistently (ie. don’t change its value when calling different functions!).
order (string) – The Post-Newtonian order that is used to translate from masses and spins to the lambda_i coordinate system. Valid orders given above.
- Returns
lambdas – The lambda coordinates for the input system(s)
- Return type
list of floats or numpy.arrays
pycbc.tmpltbank.lattice_utils module¶
- pycbc.tmpltbank.lattice_utils.generate_anstar_3d_lattice(maxv1, minv1, maxv2, minv2, maxv3, minv3, mindist)[source]¶
This function calls into LAL routines to generate a 3-dimensional array of points using the An^* lattice.
- Parameters
maxv1 (float) – Largest value in the 1st dimension to cover
minv1 (float) – Smallest value in the 1st dimension to cover
maxv2 (float) – Largest value in the 2nd dimension to cover
minv2 (float) – Smallest value in the 2nd dimension to cover
maxv3 (float) – Largest value in the 3rd dimension to cover
minv3 (float) – Smallest value in the 3rd dimension to cover
mindist (float) – Maximum allowed mismatch between a point in the parameter space and the generated bank of points.
- Returns
v1s (numpy.array) – Array of positions in the first dimension
v2s (numpy.array) – Array of positions in the second dimension
v3s (numpy.array) – Array of positions in the second dimension
- pycbc.tmpltbank.lattice_utils.generate_hexagonal_lattice(maxv1, minv1, maxv2, minv2, mindist)[source]¶
This function generates a 2-dimensional lattice of points using a hexagonal lattice.
- Parameters
maxv1 (float) – Largest value in the 1st dimension to cover
minv1 (float) – Smallest value in the 1st dimension to cover
maxv2 (float) – Largest value in the 2nd dimension to cover
minv2 (float) – Smallest value in the 2nd dimension to cover
mindist (float) – Maximum allowed mismatch between a point in the parameter space and the generated bank of points.
- Returns
v1s (numpy.array) – Array of positions in the first dimension
v2s (numpy.array) – Array of positions in the second dimension
pycbc.tmpltbank.option_utils module¶
- class pycbc.tmpltbank.option_utils.IndentedHelpFormatterWithNL(prog, indent_increment=2, max_help_position=24, width=None)[source]¶
Bases:
argparse.ArgumentDefaultsHelpFormatter
This class taken from https://groups.google.com/forum/#!topic/comp.lang.python/bfbmtUGhW8I and is used to format the argparse help messages to deal with line breaking nicer. Specfically the pn-order help is large and looks crappy without this. This function is (C) Tim Chase
- pycbc.tmpltbank.option_utils.check_ethinca_against_bank_params(ethincaParams, metricParams)[source]¶
Cross-check the ethinca and bank layout metric calculation parameters and set the ethinca metric PN order equal to the bank PN order if not previously set.
- Parameters
ethincaParams (instance of ethincaParameters) –
metricParams (instance of metricParameters) –
- class pycbc.tmpltbank.option_utils.ethincaParameters(pnOrder, cutoff, freqStep, fLow=None, full_ethinca=False, time_ethinca=False)[source]¶
Bases:
object
This class holds all of the options that are parsed in the function insert_ethinca_metric_options and all products produced using these options. It can also be initialized from the __init__ function, providing directly the options normally provided on the command line
- pycbc.tmpltbank.option_utils.get_options_from_group(option_group)[source]¶
Take an option group and return all the options that are defined in that group.
- pycbc.tmpltbank.option_utils.insert_base_bank_options(parser)[source]¶
Adds essential common options for template bank generation to an ArgumentParser instance.
- pycbc.tmpltbank.option_utils.insert_ethinca_metric_options(parser)[source]¶
Adds the options used to calculate the ethinca metric, if required.
- Parameters
parser (object) – OptionParser instance.
- pycbc.tmpltbank.option_utils.insert_mass_range_option_group(parser, nonSpin=False)[source]¶
Adds the options used to specify mass ranges in the bank generation codes to an argparser as an OptionGroup. This should be used if you want to use these options in your code.
- Parameters
parser (object) – OptionParser instance.
nonSpin (boolean, optional (default=False)) – If this is provided the spin-related options will not be added.
- pycbc.tmpltbank.option_utils.insert_metric_calculation_options(parser)[source]¶
Adds the options used to obtain a metric in the bank generation codes to an argparser as an OptionGroup. This should be used if you want to use these options in your code.
- class pycbc.tmpltbank.option_utils.massRangeParameters(minMass1, maxMass1, minMass2, maxMass2, maxNSSpinMag=0, maxBHSpinMag=0, maxTotMass=None, minTotMass=None, maxEta=None, minEta=0, max_chirp_mass=None, min_chirp_mass=None, ns_bh_boundary_mass=None, nsbhFlag=False, remnant_mass_threshold=None, ns_eos=None, use_eos_max_ns_mass=False, delta_bh_spin=None, delta_ns_mass=None)[source]¶
Bases:
object
This class holds all of the options that are parsed in the function insert_mass_range_option_group and all products produced using these options. It can also be initialized from the __init__ function providing directly the options normally provided on the command line
- default_delta_bh_spin = 0.1¶
- default_delta_ns_mass = 0.1¶
- default_ns_eos = '2H'¶
- default_nsbh_boundary_mass = 3.0¶
- class pycbc.tmpltbank.option_utils.metricParameters(pnOrder, fLow, fUpper, deltaF, f0=70, write_metric=False)[source]¶
Bases:
object
This class holds all of the options that are parsed in the function insert_metric_calculation_options and all products produced using these options. It can also be initialized from the __init__ function, providing directly the options normally provided on the command line.
- property evals¶
The eigenvalues of the parameter space. This is a Dictionary of numpy.array Each entry in the dictionary corresponds to the different frequency ranges described in vary_fmax. If vary_fmax = False, the only entry will be f_upper, this corresponds to integrals in [f_low,f_upper). This entry is always present. Each other entry will use floats as keys to the dictionary. These floats give the upper frequency cutoff when it is varying. Each numpy.array contains the eigenvalues which, with the eigenvectors in evecs, are needed to rotate the coordinate system to one in which the metric is the identity matrix.
- property evecs¶
The eigenvectors of the parameter space. This is a Dictionary of numpy.matrix Each entry in the dictionary is as described under evals. Each numpy.matrix contains the eigenvectors which, with the eigenvalues in evals, are needed to rotate the coordinate system to one in which the metric is the identity matrix.
- property evecsCV¶
The eigenvectors of the principal directions of the mu space. This is a Dictionary of numpy.matrix Each entry in the dictionary is as described under evals. Each numpy.matrix contains the eigenvectors which, with the eigenvalues in evals, are needed to rotate the coordinate system to one in which the metric is the identity matrix.
- classmethod from_argparse(opts)[source]¶
Initialize an instance of the metricParameters class from an argparse.OptionParser instance. This assumes that insert_metric_calculation_options and verify_metric_calculation_options have already been called before initializing the class.
- property metric¶
The metric of the parameter space. This is a Dictionary of numpy.matrix Each entry in the dictionary is as described under evals. Each numpy.matrix contains the metric of the parameter space in the Lambda_i coordinate system.
- property moments¶
Moments structure This contains the result of all the integrals used in computing the metrics above. It can be used for the ethinca components calculation, or other similar calculations. This is composed of two compound dictionaries. The first entry indicates which moment is being calculated and the second entry indicates the upper frequency cutoff that was used.
In all cases x = f/f0.
For the first entries the options are:
moments[‘J%d’ %(i)][f_cutoff] This stores the integral of x**((-i)/3.) * delta X / PSD(x)
moments[‘log%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.))) x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**2 x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**3 x**((-i)/3.) * delta X / PSD(x)
moments[‘loglog%d’ %(i)][f_cutoff] This stores the integral of (numpy.log(x**(1./3.)))**4 x**((-i)/3.) * delta X / PSD(x)
The second entry stores the frequency cutoff that was used when computing the integral.
- property psd¶
A pyCBC FrequencySeries holding the appropriate PSD. Return the PSD used in the metric calculation.
- property time_unprojected_metric¶
The metric of the parameter space with the time dimension unprojected. This is a Dictionary of numpy.matrix Each entry in the dictionary is as described under evals. Each numpy.matrix contains the metric of the parameter space in the Lambda_i, t coordinate system. The time components are always in the last [-1] position in the matrix.
- pycbc.tmpltbank.option_utils.verify_ethinca_metric_options(opts, parser)[source]¶
Checks that the necessary options are given for the ethinca metric calculation.
- Parameters
opts (argparse.Values instance) – Result of parsing the input options with OptionParser
parser (object) – The OptionParser instance.
- pycbc.tmpltbank.option_utils.verify_mass_range_options(opts, parser, nonSpin=False)[source]¶
Parses the metric calculation options given and verifies that they are correct.
- Parameters
opts (argparse.Values instance) – Result of parsing the input options with OptionParser
parser (object) – The OptionParser instance.
nonSpin (boolean, optional (default=False)) – If this is provided the spin-related options will not be checked.
- pycbc.tmpltbank.option_utils.verify_metric_calculation_options(opts, parser)[source]¶
Parses the metric calculation options given and verifies that they are correct.
- Parameters
opts (argparse.Values instance) – Result of parsing the input options with OptionParser
parser (object) – The OptionParser instance.
pycbc.tmpltbank.partitioned_bank module¶
- class pycbc.tmpltbank.partitioned_bank.PartitionedTmpltbank(mass_range_params, metric_params, ref_freq, bin_spacing, bin_range_check=1)[source]¶
Bases:
object
This class is used to hold a template bank partitioned into numerous bins based on position in the Cartesian parameter space where the axes are the principal components. It can also be used to hold intermediary products used while constructing (e.g.) a stochastic template bank.
- add_point_by_chi_coords(chi_coords, mass1, mass2, spin1z, spin2z, point_fupper=None, mus=None)[source]¶
Add a point to the partitioned template bank. The point_fupper and mus kwargs must be provided for all templates if the vary fupper capability is desired. This requires that the chi_coords, as well as mus and point_fupper if needed, to be precalculated. If you just have the masses and don’t want to worry about translations see add_point_by_masses, which will do translations and then call this.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
mass1 (float) – The heavier mass of the point to add.
mass2 (float) – The lighter mass of the point to add.
spin1z (float) – The [aligned] spin on the heavier body.
spin2z (float) – The [aligned] spin on the lighter body. The upper frequency cutoff to use for this point. This value must be one of the ones already calculated in the metric.
mus (numpy.array) – A 2D array where idx 0 holds the upper frequency cutoff and idx 1 holds the coordinates in the [not covaried] mu parameter space for each value of the upper frequency cutoff.
- add_point_by_masses(mass1, mass2, spin1z, spin2z, vary_fupper=False)[source]¶
Add a point to the template bank. This differs from add point to bank as it assumes that the chi coordinates and the products needed to use vary_fupper have not already been calculated. This function calculates these products and then calls add_point_by_chi_coords. This function also carries out a number of sanity checks (eg. is the point within the ranges given by mass_range_params) that add_point_by_chi_coords does not do for speed concerns.
- add_tmpltbank_from_hdf_file(hdf_fp, vary_fupper=False)[source]¶
This function will take a pointer to an open HDF File object containing a list of templates and add them into the partitioned template bank object.
- Parameters
hdf_fp (h5py.File object) – The template bank in HDF5 format.
vary_fupper (False) – If given also include the additional information needed to compute distances with a varying upper frequency cutoff.
- add_tmpltbank_from_xml_table(sngl_table, vary_fupper=False)[source]¶
This function will take a sngl_inspiral_table of templates and add them into the partitioned template bank object.
- Parameters
sngl_table (sngl_inspiral_table) – List of sngl_inspiral templates.
vary_fupper (False) – If given also include the additional information needed to compute distances with a varying upper frequency cutoff.
- calc_point_distance(chi_coords)[source]¶
Calculate distance between point and the bank. Return the closest distance.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
- Returns
min_dist (float) – The smallest SQUARED metric distance between the test point and the bank.
indexes (The chi1_bin, chi2_bin and position within that bin at which) – the closest matching point lies.
- calc_point_distance_vary(chi_coords, point_fupper, mus)[source]¶
Calculate distance between point and the bank allowing the metric to vary based on varying upper frequency cutoff. Slower than calc_point_distance, but more reliable when upper frequency cutoff can change a lot.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
point_fupper (float) – The upper frequency cutoff to use for this point. This value must be one of the ones already calculated in the metric.
mus (numpy.array) – A 2D array where idx 0 holds the upper frequency cutoff and idx 1 holds the coordinates in the [not covaried] mu parameter space for each value of the upper frequency cutoff.
- Returns
min_dist (float) – The smallest SQUARED metric distance between the test point and the bank.
indexes (The chi1_bin, chi2_bin and position within that bin at which) – the closest matching point lies.
- check_bin_existence(chi1_bin, chi2_bin)[source]¶
Given indices for bins in chi1 and chi2 space check that the bin exists in the object. If not add it. Also check for the existence of all bins within +/- self.bin_range_check and add if not present.
- find_point_bin(chi_coords)[source]¶
Given a set of coordinates in the chi parameter space, identify the indices of the chi1 and chi2 bins that the point occurs in. Returns these indices.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
- Returns
chi1_bin (int) – Index of the chi_1 bin.
chi2_bin (int) – Index of the chi_2 bin.
- get_freq_map_and_normalizations(frequency_list, upper_freq_formula)[source]¶
If using the –vary-fupper capability we need to store the mapping between index and frequencies in the list. We also precalculate the normalization factor at every frequency, which is used when estimating overlaps to account for abrupt changes in termination frequency.
- Parameters
frequency_list (array of floats) – The frequencies for which the metric has been computed and lie within the parameter space being considered.
upper_freq_formula (string) –
- get_point_from_bins_and_idx(chi1_bin, chi2_bin, idx)[source]¶
Find masses and spins given bin numbers and index.
Given the chi1 bin, chi2 bin and an index, return the masses and spins of the point at that index. Will fail if no point exists there.
- Parameters
- Returns
mass1 (float) – Mass of heavier body.
mass2 (float) – Mass of lighter body.
spin1z (float) – Spin of heavier body.
spin2z (float) – Spin of lighter body.
- output_all_points()[source]¶
Return all points in the bank.
Return all points in the bank as lists of m1, m2, spin1z, spin2z.
- Returns
mass1 (list) – List of mass1 values.
mass2 (list) – List of mass2 values.
spin1z (list) – List of spin1z values.
spin2z (list) – List of spin2z values.
- test_point_distance(chi_coords, distance_threshold)[source]¶
Test if the distance between the supplied point and the bank is less than the supplied distance theshold.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
distance_threshold (float) – The SQUARE ROOT of the metric distance to test as threshold. E.g. if you want to test to a minimal match of 0.97 you would use 1 - 0.97 = 0.03 for this value.
- Returns
True if point is within the distance threshold. False if not.
- Return type
Boolean
- test_point_distance_vary(chi_coords, point_fupper, mus, distance_threshold)[source]¶
Test if distance between point and the bank is greater than distance threshold while allowing the metric to vary based on varying upper frequency cutoff. Slower than test_point_distance, but more reliable when upper frequency cutoff can change a lot.
- Parameters
chi_coords (numpy.array) – The position of the point in the chi coordinates.
point_fupper (float) – The upper frequency cutoff to use for this point. This value must be one of the ones already calculated in the metric.
mus (numpy.array) – A 2D array where idx 0 holds the upper frequency cutoff and idx 1 holds the coordinates in the [not covaried] mu parameter space for each value of the upper frequency cutoff.
distance_threshold (float) – The SQUARE ROOT of the metric distance to test as threshold. E.g. if you want to test to a minimal match of 0.97 you would use 1 - 0.97 = 0.03 for this value.
- Returns
True if point is within the distance threshold. False if not.
- Return type
Boolean