pycbc.strain package
Submodules
pycbc.strain.calibration module
Functions for adding calibration factors to waveform templates.
- class pycbc.strain.calibration.CubicSpline(minimum_frequency, maximum_frequency, n_points, ifo_name)[source]
Bases:
Recalibrate
- apply_calibration(strain)[source]
Apply calibration model
This applies cubic spline calibration to the strain.
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- name = 'cubic_spline'
- class pycbc.strain.calibration.Recalibrate(ifo_name)[source]
Bases:
object
- abstract apply_calibration(strain)[source]
Apply calibration model
This method should be overwritten by subclasses
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- classmethod from_config(cp, ifo, section)[source]
Read a config file to get calibration options and transfer functions which will be used to intialize the model.
- Parameters:
cp (WorkflowConfigParser) – An open config file.
ifo (string) – The detector (H1, L1) for which the calibration model will be loaded.
section (string) – The section name in the config file from which to retrieve the calibration options.
- Returns:
An instance of the class.
- Return type:
instance
- map_to_adjust(strain, prefix='recalib_', **params)[source]
Map an input dictionary of sampling parameters to the adjust_strain function by filtering the dictionary for the calibration parameters, then calling adjust_strain.
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
prefix (str) – Prefix for calibration parameter names
params (dict) – Dictionary of sampling parameters which includes calibration parameters.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- name = None
pycbc.strain.gate module
Functions for applying gates to data.
- pycbc.strain.gate.add_gate_option_group(parser)[source]
Adds the options needed to apply gates to data.
- Parameters:
parser (object) – ArgumentParser instance.
- pycbc.strain.gate.apply_gates_to_fd(stilde_dict, gates)[source]
Applies the given dictionary of gates to the given dictionary of strain in the frequency domain.
Gates are applied by IFFT-ing the strain data to the time domain, applying the gate, then FFT-ing back to the frequency domain.
- Parameters:
- Returns:
Dictionary of frequency-domain strain with the gates applied.
- Return type:
- pycbc.strain.gate.apply_gates_to_td(strain_dict, gates)[source]
Applies the given dictionary of gates to the given dictionary of strain.
- Parameters:
- Returns:
Dictionary of time-domain strain with the gates applied.
- Return type:
- pycbc.strain.gate.gate_and_paint(data, lindex, rindex, invpsd, copy=True)[source]
Gates and in-paints data.
- Parameters:
data (TimeSeries) – The data to gate.
lindex (int) – The start index of the gate.
rindex (int) – The end index of the gate.
invpsd (FrequencySeries) – The inverse of the PSD.
copy (bool, optional) – Copy the data before applying the gate. Otherwise, the gate will be applied in-place. Default is True.
- Returns:
The gated and in-painted time series.
- Return type:
pycbc.strain.lines module
Functions for removing frequency lines from real data.
- pycbc.strain.lines.avg_inner_product(data1, data2, bin_size)[source]
Calculate the time-domain inner product averaged over bins.
- Parameters:
data1 (pycbc.types.TimeSeries) – First data set.
data2 (pycbc.types.TimeSeries) – Second data set, with same duration and sample rate as data1.
bin_size (float) – Duration of the bins the data will be divided into to calculate the inner product.
- Returns:
inner_prod (list) – The (complex) inner product of data1 and data2 obtained in each bin.
amp (float) – The absolute value of the median of the inner product.
phi (float) – The angle of the median of the inner product.
- pycbc.strain.lines.calibration_lines(freqs, data, tref=None)[source]
Extract the calibration lines from strain data.
- Parameters:
freqs (list) – List containing the frequencies of the calibration lines.
data (pycbc.types.TimeSeries) – Strain data to extract the calibration lines from.
tref ({None, float}, optional) – Reference time for the line. If None, will use data.start_time.
- Returns:
data – The strain data with the calibration lines removed.
- Return type:
pycbc.types.TimeSeries
- pycbc.strain.lines.clean_data(freqs, data, chunk, avg_bin)[source]
Extract time-varying (wandering) lines from strain data.
- Parameters:
freqs (list) – List containing the frequencies of the wandering lines.
data (pycbc.types.TimeSeries) – Strain data to extract the wandering lines from.
chunk (float) – Duration of the chunks the data will be divided into to account for the time variation of the wandering lines. Should be smaller than data.duration, and allow for at least a few chunks.
avg_bin (float) – Duration of the bins each chunk will be divided into for averaging the inner product when measuring the parameters of the line. Should be smaller than chunk.
- Returns:
data – The strain data with the wandering lines removed.
- Return type:
pycbc.types.TimeSeries
- pycbc.strain.lines.complex_median(complex_list)[source]
Get the median value of a list of complex numbers.
- Parameters:
complex_list (list) – List of complex numbers to calculate the median.
- Returns:
a + 1.j*b – The median of the real and imaginary parts.
- Return type:
complex number
- pycbc.strain.lines.line_model(freq, data, tref, amp=1, phi=0)[source]
Simple time-domain model for a frequency line.
- Parameters:
freq (float) – Frequency of the line.
data (pycbc.types.TimeSeries) – Reference data, to get delta_t, start_time, duration and sample_times.
tref (float) – Reference time for the line model.
amp ({1., float}, optional) – Amplitude of the frequency line.
phi ({0. float}, optional) – Phase of the frequency line (radians).
- Returns:
freq_line – A timeseries of the line model with frequency ‘freq’. The returned data are complex to allow measuring the amplitude and phase of the corresponding frequency line in the strain data. For extraction, use only the real part of the data.
- Return type:
pycbc.types.TimeSeries
- pycbc.strain.lines.matching_line(freq, data, tref, bin_size=1)[source]
Find the parameter of the line with frequency ‘freq’ in the data.
- Parameters:
- Returns:
line_model – A timeseries containing the frequency line with the amplitude and phase measured from the data.
- Return type:
pycbc.types.TimeSeries
pycbc.strain.recalibrate module
Classes and functions for adjusting strain data.
- class pycbc.strain.recalibrate.CubicSpline(minimum_frequency, maximum_frequency, n_points, ifo_name)[source]
Bases:
Recalibrate
Cubic spline recalibration
see https://dcc.ligo.org/LIGO-T1400682/public
This assumes the spline points follow np.logspace(np.log(minimum_frequency), np.log(maximum_frequency), n_points)
- Parameters:
- apply_calibration(strain)[source]
Apply calibration model
This applies cubic spline calibration to the strain.
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- name = 'cubic_spline'
- class pycbc.strain.recalibrate.PhysicalModel(freq=None, fc0=None, c0=None, d0=None, a_tst0=None, a_pu0=None, fs0=None, qinv0=None)[source]
Bases:
object
Class for adjusting time-varying calibration parameters of given strain data.
- Parameters:
strain (FrequencySeries) – The strain to be adjusted.
freq (array) – The frequencies corresponding to the values of c0, d0, a0 in Hertz.
fc0 (float) – Coupled-cavity (CC) pole at time t0, when c0=c(t0) and a0=a(t0) are measured.
c0 (array) – Initial sensing function at t0 for the frequencies.
d0 (array) – Digital filter for the frequencies.
a_tst0 (array) – Initial actuation function for the test mass at t0 for the frequencies.
a_pu0 (array) – Initial actuation function for the penultimate mass at t0 for the frequencies.
fs0 (float) – Initial spring frequency at t0 for the signal recycling cavity.
qinv0 (float) – Initial inverse quality factor at t0 for the signal recycling cavity.
- adjust_strain(strain, delta_fs=None, delta_qinv=None, delta_fc=None, kappa_c=1.0, kappa_tst_re=1.0, kappa_tst_im=0.0, kappa_pu_re=1.0, kappa_pu_im=0.0)[source]
Adjust the FrequencySeries strain by changing the time-dependent calibration parameters kappa_c(t), kappa_a(t), f_c(t), fs, and qinv.
- Parameters:
strain (FrequencySeries) – The strain data to be adjusted.
delta_fc (float) – Change in coupled-cavity (CC) pole at time t.
kappa_c (float) – Scalar correction factor for sensing function c0 at time t.
kappa_tst_re (float) – Real part of scalar correction factor for actuation function A_{tst0} at time t.
kappa_tst_im (float) – Imaginary part of scalar correction factor for actuation function A_tst0 at time t.
kappa_pu_re (float) – Real part of scalar correction factor for actuation function A_{pu0} at time t.
kappa_pu_im (float) – Imaginary part of scalar correction factor for actuation function A_{pu0} at time t.
fs (float) – Spring frequency for signal recycling cavity.
qinv (float) – Inverse quality factor for signal recycling cavity.
- Returns:
strain_adjusted – The adjusted strain.
- Return type:
- classmethod from_config(cp, ifo, section)[source]
Read a config file to get calibration options and transfer functions which will be used to intialize the model.
- Parameters:
cp (WorkflowConfigParser) – An open config file.
ifo (string) – The detector (H1, L1) for which the calibration model will be loaded.
section (string) – The section name in the config file from which to retrieve the calibration options.
- Returns:
An instance of the Recalibrate class.
- Return type:
instance
- map_to_adjust(strain, **params)[source]
Map an input dictionary of sampling parameters to the adjust_strain function by filtering the dictionary for the calibration parameters, then calling adjust_strain.
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
params (dict) – Dictionary of sampling parameters which includes calibration parameters.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- name = 'physical_model'
- classmethod tf_from_file(path, delimiter=' ')[source]
Convert the contents of a file with the columns [freq, real(h), imag(h)] to a numpy.array with columns [freq, real(h)+j*imag(h)].
- Parameters:
path (string)
delimiter ({" ", string})
- Return type:
numpy.array
- update_c(fs=None, qinv=None, fc=None, kappa_c=1.0)[source]
Calculate the sensing function c(f,t) given the new parameters kappa_c(t), kappa_a(t), f_c(t), fs, and qinv.
- Parameters:
- Returns:
c – The new sensing function c(f,t).
- Return type:
numpy.array
- update_g(fs=None, qinv=None, fc=None, kappa_tst_re=1.0, kappa_tst_im=0.0, kappa_pu_re=1.0, kappa_pu_im=0.0, kappa_c=1.0)[source]
Calculate the open loop gain g(f,t) given the new parameters kappa_c(t), kappa_a(t), f_c(t), fs, and qinv.
- Parameters:
fc (float) – Coupled-cavity (CC) pole at time t.
kappa_c (float) – Scalar correction factor for sensing function c at time t.
kappa_tst_re (float) – Real part of scalar correction factor for actuation function a_tst0 at time t.
kappa_pu_re (float) – Real part of scalar correction factor for actuation function a_pu0 at time t.
kappa_tst_im (float) – Imaginary part of scalar correction factor for actuation function a_tst0 at time t.
kappa_pu_im (float) – Imaginary part of scalar correction factor for actuation function a_pu0 at time t.
fs (float) – Spring frequency for signal recycling cavity.
qinv (float) – Inverse quality factor for signal recycling cavity.
- Returns:
g – The new open loop gain g(f,t).
- Return type:
numpy.array
- update_r(fs=None, qinv=None, fc=None, kappa_c=1.0, kappa_tst_re=1.0, kappa_tst_im=0.0, kappa_pu_re=1.0, kappa_pu_im=0.0)[source]
Calculate the response function R(f,t) given the new parameters kappa_c(t), kappa_a(t), f_c(t), fs, and qinv.
- Parameters:
fc (float) – Coupled-cavity (CC) pole at time t.
kappa_c (float) – Scalar correction factor for sensing function c at time t.
kappa_tst_re (float) – Real part of scalar correction factor for actuation function a_tst0 at time t.
kappa_pu_re (float) – Real part of scalar correction factor for actuation function a_pu0 at time t.
kappa_tst_im (float) – Imaginary part of scalar correction factor for actuation function a_tst0 at time t.
kappa_pu_im (float) – Imaginary part of scalar correction factor for actuation function a_pu0 at time t.
fs (float) – Spring frequency for signal recycling cavity.
qinv (float) – Inverse quality factor for signal recycling cavity.
- Returns:
r – The new response function r(f,t).
- Return type:
numpy.array
- class pycbc.strain.recalibrate.Recalibrate(ifo_name)[source]
Bases:
object
Base class for modifying calibration
- abstract apply_calibration(strain)[source]
Apply calibration model
This method should be overwritten by subclasses
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- classmethod from_config(cp, ifo, section)[source]
Read a config file to get calibration options and transfer functions which will be used to intialize the model.
- Parameters:
cp (WorkflowConfigParser) – An open config file.
ifo (string) – The detector (H1, L1) for which the calibration model will be loaded.
section (string) – The section name in the config file from which to retrieve the calibration options.
- Returns:
An instance of the class.
- Return type:
instance
- map_to_adjust(strain, prefix='recalib_', **params)[source]
Map an input dictionary of sampling parameters to the adjust_strain function by filtering the dictionary for the calibration parameters, then calling adjust_strain.
- Parameters:
strain (FrequencySeries) – The strain to be recalibrated.
prefix (str) – Prefix for calibration parameter names
params (dict) – Dictionary of sampling parameters which includes calibration parameters.
- Returns:
strain_adjusted – The recalibrated strain.
- Return type:
- name = None
pycbc.strain.strain module
This modules contains functions reading, generating, and segmenting strain data
- class pycbc.strain.strain.StrainBuffer(frame_src, channel_name, start_time, max_buffer, sample_rate, low_frequency_cutoff=20, highpass_frequency=15.0, highpass_reduction=200.0, highpass_bandwidth=5.0, psd_samples=30, psd_segment_length=4, psd_inverse_length=3.5, trim_padding=0.25, autogating_threshold=None, autogating_cluster=None, autogating_pad=None, autogating_width=None, autogating_taper=None, autogating_duration=None, autogating_psd_segment_length=None, autogating_psd_stride=None, state_channel=None, data_quality_channel=None, idq_channel=None, idq_state_channel=None, idq_threshold=None, dyn_range_fac=5.902958103587057e+20, psd_abort_difference=None, psd_recalculate_difference=None, force_update_cache=True, increment_update_cache=None, analyze_flags=None, data_quality_flags=None, dq_padding=0)[source]
Bases:
DataBuffer
- add_hard_count()[source]
Reset the countdown timer, so that we don’t analyze data long enough to generate a new PSD.
- advance(blocksize, timeout=10)[source]
Advanced buffer blocksize seconds.
Add blocksize seconds more to the buffer, push blocksize seconds from the beginning.
- Parameters:
blocksize (int) – The number of seconds to attempt to read from the channel
- Returns:
status – Returns True if this block is analyzable.
- Return type:
boolean
- check_psd_dist(min_dist, max_dist)[source]
Check that the horizon distance of a detector is within a required range. If so, return True, otherwise log a warning and return False.
- property end_time
Return the end time of the current valid segment of data
- classmethod from_cli(ifo, args)[source]
Initialize a StrainBuffer object (data reader) for a particular detector.
- invalidate_psd()[source]
Make the current PSD invalid. A new one will be generated when it is next required
- near_hwinj()[source]
Check that the current set of triggers could be influenced by a hardware injection.
- null_advance_strain(blocksize)[source]
Advance and insert zeros
- Parameters:
blocksize (int) – The number of seconds to attempt to read from the channel
- overwhitened_data(delta_f)[source]
Return overwhitened data
- Parameters:
delta_f (float) – The sample step to generate overwhitened frequency domain data for
- Returns:
htilde – Overwhited strain data
- Return type:
- property start_time
Return the start time of the current valid segment of data
- class pycbc.strain.strain.StrainSegments(strain, segment_length=None, segment_start_pad=0, segment_end_pad=0, trigger_start=None, trigger_end=None, filter_inj_only=False, injection_window=None, allow_zero_padding=False)[source]
Bases:
object
Class for managing manipulation of strain data for the purpose of matched filtering. This includes methods for segmenting and conditioning.
- fourier_segments()[source]
Return a list of the FFT’d segments. Return the list of FrequencySeries. Additional properties are added that describe the strain segment. The property ‘analyze’ is a slice corresponding to the portion of the time domain equivalent of the segment to analyze for triggers. The value ‘cumulative_index’ indexes from the beginning of the original strain series.
- classmethod from_cli(opt, strain)[source]
Calculate the segmentation of the strain data for analysis from the command line options.
- classmethod from_cli_multi_ifos(opt, strain_dict, ifos)[source]
Calculate the segmentation of the strain data for analysis from the command line options.
- classmethod from_cli_single_ifo(opt, strain, ifo)[source]
Calculate the segmentation of the strain data for analysis from the command line options.
- required_opts_list = ['--segment-length', '--segment-start-pad', '--segment-end-pad']
- pycbc.strain.strain.create_memory_and_engine_for_class_based_fft(npoints_time, dtype, delta_t=1, ifft=False, uid=0)[source]
Create memory and engine for class-based FFT/IFFT
Currently only supports R2C FFT / C2R IFFTs, but this could be expanded if use-cases arise.
- Parameters:
npoints_time (int) – Number of time samples of the real input vector (or real output vector if doing an IFFT).
dtype (np.dtype) – The dtype for the real input vector (or real output vector if doing an IFFT). np.float32 or np.float64 I think in all cases.
delta_t (float (default: 1)) – delta_t of the real vector. If not given this will be set to 1, and we will assume it is not needed in the returned TimeSeries/FrequencySeries
ifft (boolean (default: False)) – By default will use the FFT class, set to true to use IFFT.
uid (int (default: 0)) – Provide a unique identifier. This is used to provide a separate set of memory in the cache, for instance if calling this from different codes.
- pycbc.strain.strain.detect_loud_glitches(strain, psd_duration=4.0, psd_stride=2.0, psd_avg_method='median', low_freq_cutoff=30.0, threshold=50.0, cluster_window=5.0, corrupt_time=4.0, high_freq_cutoff=None, output_intermediates=False)[source]
Automatic identification of loud transients for gating purposes.
This function first estimates the PSD of the input time series using the FindChirp Welch method. Then it whitens the time series using that estimate. Finally, it computes the magnitude of the whitened series, thresholds it and applies the FindChirp clustering over time to the surviving samples.
- Parameters:
strain (TimeSeries) – Input strain time series to detect glitches over.
psd_duration ({float, 4}) – Duration of the segments for PSD estimation in seconds.
psd_stride ({float, 2}) – Separation between PSD estimation segments in seconds.
psd_avg_method ({string, 'median'}) – Method for averaging PSD estimation segments.
low_freq_cutoff ({float, 30}) – Minimum frequency to include in the whitened strain.
threshold ({float, 50}) – Minimum magnitude of whitened strain for considering a transient to be present.
cluster_window ({float, 5}) – Length of time window to cluster surviving samples over, in seconds.
corrupt_time ({float, 4}) – Amount of time to be discarded at the beginning and end of the input time series.
high_frequency_cutoff ({float, None}) – Maximum frequency to include in the whitened strain. If given, the input series is downsampled accordingly. If omitted, the Nyquist frequency is used.
output_intermediates ({bool, False}) – Save intermediate time series for debugging.
- pycbc.strain.strain.execute_cached_fft(invec_data, normalize_by_rate=True, ifft=False, copy_output=True, uid=0)[source]
Executes a cached FFT
- Parameters:
invec_data (Array) – Array which will be used as input when fft_class is executed.
normalize_by_rate (boolean (optional, default:False)) – If True, then normalize by delta_t (for an FFT) or delta_f (for an IFFT).
ifft (boolean (optional, default:False)) – If true assume this is an IFFT and multiply by delta_f not delta_t. Will do nothing if normalize_by_rate is False.
copy_output (boolean (optional, default:True)) – If True we will copy the output into a new array. This avoids the issue that calling this function again might overwrite output. However, if you know that the output array will not be used before this function might be called again with the same length, then setting this to False will provide some increase in efficiency. The uid can also be used to help ensure that data doesn’t get unintentionally overwritten!
uid (int (default: 0)) – Provide a unique identifier. This is used to provide a separate set of memory in the cache, for instance if calling this from different codes.
- pycbc.strain.strain.execute_cached_ifft(*args, **kwargs)[source]
Executes a cached IFFT
- Parameters:
invec_data (Array) – Array which will be used as input when fft_class is executed.
normalize_by_rate (boolean (optional, default:False)) – If True, then normalize by delta_t (for an FFT) or delta_f (for an IFFT).
copy_output (boolean (optional, default:True)) – If True we will copy the output into a new array. This avoids the issue that calling this function again might overwrite output. However, if you know that the output array will not be used before this function might be called again with the same length, then setting this to False will provide some increase in efficiency. The uid can also be used to help ensure that data doesn’t get unintentionally overwritten!
uid (int (default: 0)) – Provide a unique identifier. This is used to provide a separate set of memory in the cache, for instance if calling this from different codes.
- pycbc.strain.strain.from_cli(opt, dyn_range_fac=1, precision='single', inj_filter_rejector=None)[source]
Parses the CLI options related to strain data reading and conditioning.
- Parameters:
opt (object) – Result of parsing the CLI with OptionParser, or any object with the required attributes (gps-start-time, gps-end-time, strain-high-pass, pad-data, sample-rate, (frame-cache or frame-files), channel-name, fake-strain, fake-strain-seed, fake-strain-from-file, gating_file).
dyn_range_fac ({float, 1}, optional) – A large constant to reduce the dynamic range of the strain.
precision (string) – Precision of the returned strain (‘single’ or ‘double’).
inj_filter_rejector (InjFilterRejector instance; optional, default=None) – If given send the InjFilterRejector instance to the inject module so that it can store a reduced representation of injections if necessary.
- Returns:
strain – The time series containing the conditioned strain data.
- Return type:
- pycbc.strain.strain.from_cli_multi_ifos(opt, ifos, inj_filter_rejector_dict=None, **kwargs)[source]
Get the strain for all ifos when using the multi-detector CLI
- pycbc.strain.strain.from_cli_single_ifo(opt, ifo, inj_filter_rejector=None, **kwargs)[source]
Get the strain for a single ifo when using the multi-detector CLI
- pycbc.strain.strain.gate_data(data, gate_params)[source]
Apply a set of gating windows to a time series.
Each gating window is defined by a central time, a given duration (centered on the given time) to zero out, and a given duration of smooth tapering on each side of the window. The window function used for tapering is a Tukey window.
- Parameters:
data (TimeSeries) – The time series to be gated.
gate_params (list) – List of parameters for the gating windows. Each element should be a list or tuple with 3 elements: the central time of the gating window, the half-duration of the portion to zero out, and the duration of the Tukey tapering on each side. All times in seconds. The total duration of the data affected by one gating window is thus twice the second parameter plus twice the third parameter.
- Returns:
data – The gated time series.
- Return type:
- pycbc.strain.strain.insert_strain_option_group(parser, gps_times=True)[source]
Add strain-related options to the optparser object.
Adds the options used to call the pycbc.strain.from_cli function to an optparser as an OptionGroup. This should be used if you want to use these options in your code.
- pycbc.strain.strain.insert_strain_option_group_multi_ifo(parser, gps_times=True)[source]
Adds the options used to call the pycbc.strain.from_cli function to an optparser as an OptionGroup. This should be used if you want to use these options in your code.
- pycbc.strain.strain.next_power_of_2(n)[source]
Return the smallest integer power of 2 larger than the argument.
- pycbc.strain.strain.verify_strain_options(opts, parser)[source]
Sanity check provided strain arguments.
Parses the strain data CLI options and verifies that they are consistent and reasonable.
- pycbc.strain.strain.verify_strain_options_multi_ifo(opts, parser, ifos)[source]
Sanity check provided strain arguments.
Parses the strain data CLI options and verifies that they are consistent and reasonable.
- Parameters:
opt (object) – Result of parsing the CLI with OptionParser, or any object with the required attributes (gps-start-time, gps-end-time, strain-high-pass, pad-data, sample-rate, frame-cache, channel-name, fake-strain, fake-strain-seed).
parser (object) – OptionParser instance.
ifos (list of strings) – List of ifos for which to verify options for
Module contents
- pycbc.strain.read_model_from_config(cp, ifo, section='calibration')[source]
Returns an instance of the calibration model specified in the given configuration file.
- Parameters:
cp (WorflowConfigParser) – An open config file to read.
ifo (string) – The detector (H1, L1) whose model will be loaded.
section ({"calibration", string}) – Section name from which to retrieve the model.
- Returns:
An instance of the calibration model class.
- Return type:
instance