>
Contents

minke package

Submodules

minke.antenna module

This module provides functions to calculate antenna factors for a given time, a given sky location and a given detector.

Adapted from the implementation in pylal

minke.antenna.response(gpsTime, rightAscension, declination, inclination, polarization, unit, detector)[source]

Calculates the antenna factors for a detector ‘detector’ (e.g. ‘H1’) at a given gps time (as integer) for a given sky location (rightAscension, declination) in some unit (degree/radians). This computation also takes into account a specific inclination and polarization.

The returned values are: (f-plus, f-cross, f-average, q-value).

Example: antenna.response( 854378604.780, 11.089, 42.308, 0, 0, ‘radians’, ‘H1’ )

minke.antenna.timeDelay(gpsTime, rightAscension, declination, unit, det1, det2)[source]

Calculates the time delay in seconds between the detectors ‘det1’ and ‘det2’ (e.g. ‘H1’) for a sky location at (rightAscension and declination) which must be given in certain units (‘radians’ or ‘degree’). The time is passes as GPS time. A positive time delay means the GW arrives first at ‘det2’, then at ‘det1’.

Example: antenna.timeDelay( 877320548.000, 355.084,31.757, ‘degree’,’H1’,’L1’) 0.0011604683260994519 Given these values, the signal arrives first at detector L1, and 1.16 ms later at H2

minke.distribution module

minke.distribution.burst_dist(minimum, maximum, size=1)[source]

Generate an hrss drawn from the distance distribution [ r + 50/r ] used for Burst allsky mock data challenges.

Parameters
minimumfloat

The lowest hrss to be produced.

maximumfloat

The largest hrss to be produced.

sizeint

The number of draws to produce. Defaults to 1.

minke.distribution.even_time(start, stop, rate, jitter=0)[source]

Produce an evenly-distributed set of times which has some jitter.

Parameters
startint

The gps start time of the set.

stopint

The gps end time of the set.

rateint

The rate of events per year.

jitterfloat

The number of seconds around which the time can “jitter”. This has the effect of adding a small random time on to each time returned in the distribution, up to a maximum of half the jitter value, or taking away up to half the jitter value.

Returns
timesndarray

A numpy array of gps times where injections should be made.

minke.distribution.favorable_sky(net, time)[source]

Wander through the skies, searching for a most favorable location. — draw extrinsic parameters as if the network antenna pattern magnitude were the PDF.

minke.distribution.log_uniform(lower, upper, number)[source]

Draw uniformly in the log of a predefined range.

Parameters
lowerfloat

The lower hrss value (n.b. not the lower log(hrss) value)

upperfloat

The upper hrss value (n.b. not the upper log(hrss) value)

numberint

The number of samples to be drawn from the distribution.

Returns
samplefloat

An sample value drawn from the log uniform distribution.

minke.distribution.sky_params(net, time, ra, dec, psi, incl=3.141592653589793)[source]
minke.distribution.supernova_angle(num, divisions=10)[source]

Draw from a discrete distribution of angles.

minke.distribution.uniform_dec(num)[source]

Declination distribution: uniform in sin(dec), which leads to a uniform distribution across all declinations.

Parameters
numint

The number of random declinations to produce.

minke.distribution.uniform_interval(interval, num)[source]

Return a number, or a list of numbers which are sampled from a uniform distribution.

Parameters
intervaltuple (lower, upper)

The interval which the uniform distribution covers.

numint

The number of samples which should be drawn from the distribution.

Returns
samplefloat or list

A sample, or a list of samples.

Notes

http://docs.scipy.org/doc/numpy/reference/generated/numpy.random.random_sample.html#numpy.random.random_sample

minke.distribution.uniform_phi(num)[source]

Uniform in (0, 2pi) distribution. num controls the number of draws.

minke.distribution.uniform_sky(number=1)[source]

Get a set of (RA, declination, polarization) drawn from an isotropic distribution over the whole sky.

Parameters
numberint

The number of random sky locations and polarisations to be produced.

Returns
raarray of float

Randomly drawn right ascensions.

decarray of float

Randomly drawn declinations.

polarray of float

Randomly drawn polarisations.

minke.distribution.uniform_theta(num)[source]

Uniform in cos distribution. num controls the number of draws.

minke.distribution.uniform_time(tstart, tstop, number)[source]

Get a set of randomized (integer) event times.

minke.mattermost module

This is an experminetal module to allow Minke production pipelines to communicate wasily with the Mattermost chat server which is used by LIGO.

minke.mattermost.send_mattermost(message, hook)[source]

minke.mdctools module

88b d88 88 88 888b d888 “” 88 88`8b d8’88 88 88 8b d8’ 88 88 8b,dPPYba, 88 ,d8 ,adPPYba, 88 `8b d8’ 88 88 88P’ `”8a 88 ,a8” a8P_____88 88 `8b d8’ 88 88 88 88 8888[ 8PP” 88 `888’ 88 88 88 88 88”Yba, “8b, ,aa 88 `8’ 88 88 88 88 88 `Y8a `”Ybbd8”’


This file is a part of Minke, a tool for generating simulated gravitational wave signals, used for characterising and training search algorithms.

Minke was created by Daniel Williams, based on work started by Chris Pankow and others, and is built around the LALSimulation library.

class minke.mdctools.Frame(start, duration, ifo, number=- 1)[source]

Bases: object

Represents a frame, in order to prepare the injection frames

Methods

generate_gwf(mdc, directory[, project, …])

Produce the gwf file which corresponds to the MDC set over the period of this frame.

get_rowlist(mdcs)

Return the rows from an MDC set which correspond to this frame.

calculate_n_injections

generate_log

calculate_n_injections(mdcs)[source]
generate_gwf(mdc, directory, project='Minke', channel='SCIENCE', force=False, rate=16384.0)[source]

Produce the gwf file which corresponds to the MDC set over the period of this frame.

Parameters
mdcMDCSet object

The MDC set which should be used to produce this frame.

directorystr

The root directory where all of the frames are to be stored, for example “/home/albert.einstein/data/mdc/frames/” would cause the SineGaussian injections to be made in the directories under “/home/albert.einstein/data/mdc/frames/sg”

projectstr

The name of the project which this frame is a part of. Defaults to ‘Minke’.

channelstr

The name of the channel which the injections should be made into. This is prepended by the initials for each interferometer, so there will be a channel for each interferometer in the gwf.

forcebool

If true this forces the recreation of a GWF file even if it already exists.

generate_log(mdc)[source]
get_rowlist(mdcs)[source]

Return the rows from an MDC set which correspond to this frame.

Parameters
mdcsMDCSet object

The set of MDCs from which the rows are to be found.

class minke.mdctools.FrameSet(frame_list)[source]

Bases: object

Methods

full_frameset(mdc, directory[, channel, force])

Produce the gwf files which corresponds to the MDC set over the period of the frames in this collection.

full_logfile(mdc, location)

Produce a log file for the entire frame set

full_frameset(mdc, directory, channel='SCIENCE', force=False)[source]

Produce the gwf files which corresponds to the MDC set over the period of the frames in this collection.

Parameters
mdcMDCSet object

The MDC set which should be used to produce this frame.

directorystr

The root directory where all of the frames are to be stored, for example “/home/albert.einstein/data/mdc/frames/” would cause the SineGaussian injections to be made in the directories under “/home/albert.einstein/data/mdc/frames/sg”

channelstr

The name of the channel which the injections should be made into. This is prepended by the initials for each interferometer, so there will be a channel for each interferometer in the gwf.

forcebool

If true this forces the recreation of a GWF file even if it already exists.

full_logfile(mdc, location)[source]

Produce a log file for the entire frame set

class minke.mdctools.HWFrameSet(ifos=['H1', 'L1'])[source]

Bases: object

Methods

full_frameset(mdc, directory[, force])

Produce the gwf files which corresponds to the MDC set over the period of the frames in this collection.

full_frameset(mdc, directory, force=False)[source]

Produce the gwf files which corresponds to the MDC set over the period of the frames in this collection.

Parameters
mdcMDCSet object

The MDC set which should be used to produce this frame.

directorystr

The root directory where all of the frames are to be stored, for example “/home/albert.einstein/data/mdc/frames/” would cause the SineGaussian injections to be made in the directories under “/home/albert.einstein/data/mdc/frames/sg”

forcebool

If true this forces the recreation of a GWF file even if it already exists.

class minke.mdctools.HWInj(ifos)[source]

Bases: minke.mdctools.Frame

Represents a hardware injection frame.

Injection frames must be an ASCII file of the hoft sampled at the antenna sampling rate, appropriately convolved with an antenna response function.

As a result of the simplicity of this specific output format we do not need information such as start-time in the file itself, however we should have a sensible naming scheme for the ASCII files since they will need to be produced as sidecars for an xml file.

Methods

generate_gwf(mdc, directory[, project, …])

Produce the gwf file which corresponds to the MDC set over the period of this frame.

generate_pcal(mdc, directory[, force, rate])

Produce the PCAL-ready hardware injection files as an ASCII list sampled at the detector’s sample rate.

get_rowlist(mdcs)

Return the rows from an MDC set which correspond to this frame.

calculate_n_injections

generate_log

generate_pcal(mdc, directory, force=False, rate=16384)[source]

Produce the PCAL-ready hardware injection files as an ASCII list sampled at the detector’s sample rate.

Parameters
mdcMDCSet object

The signal set which should be used to generate the frame.

directorystr

The root directory where all of the frames are to be stored, for example “/home/albert.einstein/data/mdc/frames/” would cause the SineGaussian injections to be made in the directories under “/home/albert.einstein/data/mdc/frames/sg”

forcebool

If true this forces the regeneration of the file, even if it already exists.

class minke.mdctools.MDCSet(detectors, name='MDC Set', table_type='burst')[source]

Bases: object

Methods

directory_path()

Generate the directory where the frames from this MDC should be stored, so, e.g.

gravEn_row(row, frame)

Produces a gravEn-style log row for a row of the simBurstTable.

load_xml(filename[, full, start, stop])

Load the MDC Set from an XML file containing the SimBurstTable.

plot_hist(parameter)

Plot a histogram of a waveform parameter.

plot_skymap()

Plot a skymap of the injections distribution in RA and DEC on a Hammer projection.

save_xml(filename)

Save the MDC set as an XML SimBurstTable.

directory_path()[source]

Generate the directory where the frames from this MDC should be stored, so, e.g. Gaussians 0d100 would go in “ga/ga0d100/”

Returns
str

the folder structure

gravEn_row(row, frame)[source]

Produces a gravEn-style log row for a row of the simBurstTable.

Parameters
rowint

The row number of the waveforms to be measured

Returns
str

A string in the gravEn format which describes the injection.

hist_parameters = {'BTLWNB': ['hrss', 'ra', 'dec'], 'Dimmelmeier+08': ['hrss', 'ra', 'dec'], 'Gaussian': ['hrss', 'psi', 'ra', 'dec'], 'SineGaussian': ['hrss', 'psi', 'ra', 'dec'], 'StringCusp': ['amplitude', 'ra', 'dec']}
inj_families_abb = {'ADI': 'adi', 'BBHRingdown': 'rng', 'BTLWNB': 'wnb', 'Dimmelmeier+08': 'd08', 'Gaussian': 'ga', 'GenericRingdown': 'gng', 'Mueller+12': 'm12', 'Ott+13': 'o13', 'Scheidegger+10': 's10', 'SineGaussian': 'sg', 'StringCusp': 'sc', 'Yakunin+10': 'y10'}
inj_families_names = {'adi': 'ADI', 'd08': 'Dimmelmeier+08', 'ga': 'Gaussian', 'gng': 'GenericRingdown', 'm12': 'Mueller+12', 'o13': 'Ott+13', 'rng': 'BBHRingdown', 's10': 'Scheidegger+10', 'sc': 'StringCusp', 'sg': 'SineGaussian', 'wnb': 'BTLWNB', 'y10': 'Yakunin+10'}
load_xml(filename, full=True, start=None, stop=None)[source]

Load the MDC Set from an XML file containing the SimBurstTable.

Parameters
filenamestr

The filename of the XML file.

fullbool

If this is true (which is the default) then all of the calculated parameters are computed from the waveform definintion.

startfloat

The time at which the xml read-in should start. The default is “None”, in which case the xml file will be read-in from the start.

endfloat

The last time to be read from the xml file. The default is None, which causes the xml to be read right-up to the last time in the file.

plot_hist(parameter)[source]

Plot a histogram of a waveform parameter.

Parameters
parameterstr

The name of the simburst table parameter which is desired for the plot.

Returns
matplotlib figure
plot_skymap()[source]

Plot a skymap of the injections distribution in RA and DEC on a Hammer projection.

Returns
matplotlib figure
save_xml(filename)[source]

Save the MDC set as an XML SimBurstTable.

Parameters
filenamestr

The location to save the xml file. The output is gzipped, so ending it with a “.gz” would stick with convention.

waveforms = []
exception minke.mdctools.TableTypeError[source]

Bases: Exception

minke.mdctools.mkdir(path)[source]

Make all of the tree of directories in a given path if they don’t already exist.

Parameters
pathstr

The path to the desired directory.

minke.mdctools.source_from_dict(params)[source]
minke.mdctools.source_from_row(row)[source]

minke.noise module

class minke.noise.NoiseTimeseries(psd, length, epoch, rate=16384.0, seed=0)[source]

Bases: object

Create timeseries of noise coloured by a PSD.

class minke.noise.PSD(filename, fmin=20, fmax=3000, rate=16384.0)[source]

Bases: object

Create a power spectral distribution from a file which specifies the PSD.

Methods

SNR(waveform, detectors[, fmin, fmax])

Calculate the SNR for a given waveform compared to this SNR.

plot()

Plot the PSD

SNR(waveform, detectors, fmin=0, fmax=- 1)[source]

Calculate the SNR for a given waveform compared to this SNR.

Parameters
waveformminke source object

The waveform which should have its snr measured

detectorslist

A list of detector names in the format “L1” for the Livingston 4km detector, etc.

fminfloat

The minimum frequency at which to evaluate the SNR. If 0 the lower frequency is unbounded.

fmaxfloat

The maximum frequency at which the SNR is evaluated. If negative the upper frequency is unbounded.

Returns
network snrfloat

The SNR over the whole network described by the detector list

SNRslist

A list of SNRs for each detector

plot()[source]

Plot the PSD

minke.sources module

class minke.sources.ADI(time, sky_dist=<function uniform_sky>, filepath='stamp_adi_a_tapered.mat', decomposed_path=None)[source]

Bases: minke.sources.LongDuration

Accretion disk instability waveforms which are generated using the method described in LIGO-T1100093, at https://dcc.ligo.org/LIGO-T1100093. The waveforms are based off a model by MH van Putten,

M. H. van Putten, A. Levinson, H. K. Lee, T. Regimbau, M. Punturo, and G. M. Harry. Phys. Rev. D., 69(4), 044007, 2004. M. H. van Putten. Phys. Rev. Lett., 87(9), 091101, 2001.

The waveforms are stored in .mat binary files which can be read-in by SciPy.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of the ADI waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

decompose(numrel_file, sample_rate=16384.0, step_back=0.01)[source]

Produce the spherial harmonic decompositions of the ADI waveform. This is a special case since it is axisymmetric.

Parameters
numrel_filestr

The location of the numerical relativity waveform file.

sample_ratefloat

The sample rate of the output NR file. Defaults to 16384.0 Hz, and should be the same as the data rate of the detector.

step_backfloat

The amount of time, in seconds, of the data which should be included before the peak amplitude. Defaults to 0.01 sec.

Returns
decompositionndarray

The re-interpolated file at the desired sample rate which is in the <time hp hc> format which can be accepted by LALSimulation.

waveform = 'ADI'
class minke.sources.BBHRingdown(time, phi0, mass, spin, massloss, distance, inclination, sky_dist=<function uniform_sky>)[source]

Bases: minke.sources.Ringdown

A class to represent BBH ringdowns.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'BBHRingdown'
class minke.sources.Dimmelmeier08(time, distance=0.01, sky_dist=<function uniform_sky>, filepath='signal_s15a2o05_ls.dat', decomposed_path=None)[source]

Bases: minke.sources.Supernova

The Dimmelmeier08 waveform.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of the Dimmelmeier numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

decompose(numrel_file, sample_rate=16384.0, step_back=0.01)[source]

Produce the spherial harmonic decompositions of the Dimmelmeier numerical waveform. This is a special case since it is axisymmetric.

Parameters
numrel_filestr

The location of the numerical relativity waveform file.

sample_ratefloat

The sample rate of the NR file. Defaults to 16384.0 Hz.

step_backfloat

The amount of time, in seconds, of the data which should be included before the peak amplitude. Defaults to 0.01 sec.

Returns
decompositionndarray

The l=2 mode spherical decompositions of the waveform.

waveform = 'Dimmelmeier+08'
class minke.sources.Gaussian(duration, hrss, time, sky_dist=<function uniform_sky>, seed=0)[source]

Bases: minke.sources.Waveform

A class to represent a Gaussian injection.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'Gaussian'
class minke.sources.Hyperbolic(datafile, total_mass, distance, extraction, sky_dist=<function uniform_sky>, **kwargs)[source]

Bases: minke.sources.Numerical2Column

Attributes
extraction

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

class minke.sources.LongDuration[source]

Bases: minke.sources.Supernova

A superclass to handle the spherial harmonic decompositions which long duration numerical relativity bursts may require.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of a numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

supernova = True
waveform = 'LongDuration'
class minke.sources.Mueller2012(theta, phi, time, distance=0.01, sky_dist=<function uniform_sky>, filepath=None, family='L15-3', decomposed_path=None)[source]

Bases: minke.sources.Supernova

The Mueller2012 waveform.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of a numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

decompose(numrel_file, sample_rate=16384.0, step_back=0.01)[source]

Produce the spherial harmonic decompositions of a numerical waveform.

Parameters
numrel_filestr

The location of the numerical relativity waveform file.

sample_ratefloat

The sample rate of the NR file. Defaults to 16384.0 Hz.

step_backfloat

The amount of time, in seconds, of the data which should be included before the peak amplitude. Defaults to 0.01 sec.

distancefloat

The distance, in megaparsecs, from the observer at which the NR waveforms were simulated. Defaults to 10 kpc (i.e. 10e-3 Mpc).

Returns
decompositionndarray

The l=2 mode spherical decompositions of the waveform.

has_memory = True
waveform = 'Mueller+12'
class minke.sources.Numerical2Column[source]

Bases: minke.sources.Waveform

A superclass to handle ninja-based numerical relativity waveforms.

Attributes
extraction

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

extraction = None
interpolate(x_old, y_old, x_new)[source]

Convenience funtion to avoid repeated code

supernova = False
waveform = 'Numerical'
class minke.sources.Ott2013(theta, phi, time, sky_dist=<function uniform_sky>, distance=0.01, filepath=None, family='s27fheat1p05', decomposed_path=None)[source]

Bases: minke.sources.Supernova

The Ott+2013 supernova waveform

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of a numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

has_memory = True
waveform = 'Ott+13'
class minke.sources.Ringdown[source]

Bases: minke.sources.Waveform

A class to handle Rindown waveforms.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

table_type

alias of glue.ligolw.lsctables.SimRingdownTable

waveform = 'GenericRingdown'
class minke.sources.Scheidegger2010(theta, phi, time, distance=0.01, sky_dist=<function uniform_sky>, filepath=None, family='R1E1CA_L', decomposed_path=None)[source]

Bases: minke.sources.Supernova

The Scheidegger2010 waveform.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of a numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'Scheidegger+10'
class minke.sources.SineGaussian(q, frequency, hrss, polarisation, time, sky_dist=<function uniform_sky>, seed=0)[source]

Bases: minke.sources.Waveform

A class to represent a SineGaussian injection.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'SineGaussian'
class minke.sources.StringCusp(amplitude, f_max, time, sky_dist=<function uniform_sky>)[source]

Bases: minke.sources.Waveform

A class to represent a StringCusp injection.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'StringCusp'
class minke.sources.Supernova[source]

Bases: minke.sources.Waveform

A superclass to handle the spherial harmonic decompositions which all supernova waveforms require.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of a numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz, l=2, m=2)[source]

Construct the expansion parameters Hlm from T1000553. Returns the expansion parameters for l=2, m=m

decompose(numrel_file, sample_rate=16384.0, step_back=0.01)[source]

Produce the spherial harmonic decompositions of a numerical waveform.

Parameters
numrel_filestr

The location of the numerical relativity waveform file.

sample_ratefloat

The sample rate of the NR file. Defaults to 16384.0 Hz.

step_backfloat

The amount of time, in seconds, of the data which should be included before the peak amplitude. Defaults to 0.01 sec.

Returns
decompositionndarray

The l=2 mode spherical decompositions of the waveform.

file_distance = 0.01
generate_tail(sampling=16384.0, length=1, h_max=1e-23, h_min=0)[source]

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

This code was adapted from an iPython notebook provided by Marek Szczepanczyk.

The tail needs to be added to the waveform after all of the other corrections have been applied (DW: I think)

Parameters
samplingfloat

The sample rate of the injection data. By default this is 16384 Hz, which is the standard advanced LIGO sampling rate.

lengthfloat

The length of the tail to be added, in seconds; defaults to 1.

h_maxfloat

The strain at the beginning of the tail – the strain at the end of the NR data.

Notes

  • TODO Confirm that the tail is added-on after the waveform is

convolved with the antenna pattern.

has_memory = False
interpolate(x_old, y_old, x_new)[source]

Convenience funtion to avoid repeated code

supernova = True
waveform = 'Supernova'
class minke.sources.Waveform[source]

Bases: object

Generic container for different source types. Currently, it checks for the waveform type and initializes itself appropriately. In the future, different sources should subclass this and override the generation routines.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

expnum = 1
generate_tail(sampling=16384.0, length=1, h_max=1e-23, h_min=0)[source]

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

This code was adapted from an iPython notebook provided by Marek Szczepanczyk.

The tail needs to be added to the waveform after all of the other corrections have been applied (DW: I think)

Parameters
samplingfloat

The sample rate of the injection data. By default this is 16384 Hz, which is the standard advanced LIGO sampling rate.

lengthfloat

The length of the tail to be added, in seconds; defaults to 1.

h_maxfloat

The strain at the beginning of the tail – the strain at the end of the NR data.

Notes

  • TODO Confirm that the tail is added-on after the waveform is

convolved with the antenna pattern.

interpolate(x_old, y_old, x_new, method='linear')[source]

Convenience funtion to avoid repeated code

numrel_data = []
params = {}
parse_polarisation(polarisation)[source]

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

Parameters
polarisationstr, {‘linear’, ‘circular’, ‘elliptical’, ‘inclination’}

The description of the polarisation, in words.

plot(figsize=(10, 5))[source]

Produce a plot of the injection.

sim = []
table_type

alias of glue.ligolw.lsctables.SimBurstTable

waveform = 'Generic'
class minke.sources.WhiteNoiseBurst(duration, bandwidth, frequency, time, hrss=None, egw=None, sky_dist=<function uniform_sky>, seed=0)[source]

Bases: minke.sources.Waveform

A class to represent a WNB injection.

Methods

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new[, method])

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

waveform = 'BTLWNB'
class minke.sources.Yakunin10(time, distance=0.01, sky_dist=<function uniform_sky>, filepath='Yakunin2010/hplus-B12-WH07_tail.txt', decomposed_path=None)[source]

Bases: minke.sources.Supernova

The Yakunin10 waveform.

Methods

construct_Hlm(Ixx, Ixy, Ixz, Iyy, Iyz, Izz)

Construct the expansion parameters Hlm from T1000553.

decompose(numrel_file[, sample_rate, step_back])

Produce the spherial harmonic decompositions of the Dimmelmeier numerical waveform.

generate_tail([sampling, length, h_max, h_min])

Generate a “low frequency tail” to append to the end of the waveform to overcome problems related to memory in the waveform.

interpolate(x_old, y_old, x_new)

Convenience funtion to avoid repeated code

parse_polarisation(polarisation)

Convert a string description of a polarisation to an ellipse eccentricity and an ellipse angle.

plot([figsize])

Produce a plot of the injection.

table_type

decompose(numrel_file, sample_rate=16384.0, step_back=0.01)[source]

Produce the spherial harmonic decompositions of the Dimmelmeier numerical waveform. This is a special case since it is axisymmetric.

Parameters
numrel_filestr

The location of the numerical relativity waveform file.

sample_ratefloat

The sample rate of the NR file. Defaults to 16384.0 Hz.

step_backfloat

The amount of time, in seconds, of the data which should be included before the peak amplitude. Defaults to 0.01 sec.

Returns
decompositionndarray

The l=2 mode spherical decompositions of the waveform.

waveform = 'Yakunin+10'

Module contents