Omicron  4.2.0
An algorithm to detect and characterize transient events in gravitational-wave detectors
OxEvent Class Reference

Ox event management. More...

#include <OxEvent.h>

Inheritance diagram for OxEvent:
Collaboration diagram for OxEvent:

Public Member Functions

bool CreateOutputFile (const unsigned int aTime)
 Creates an ouput file. More...
 
Long64_t GetEntries (void)
 Returns the number of events (read-only). More...
 
int GetEntry (const Long64_t aEntry)
 Loads an entry in the event TTree (read-only). More...
 
UChar_t GetEventAnalysisIndex (void)
 Returns the analysis index for the current event. More...
 
double GetEventFrequency (void)
 Returns the frequency for the current event [Hz]. More...
 
double GetEventRank (void)
 Returns the rank for the current event. More...
 
unsigned int GetEventSliceIndex (void)
 Returns the slice index for the current event. More...
 
double GetSliceHalfDuration (void)
 Returns the half duration of one slice [s]. More...
 
unsigned int GetSliceIndexCenter (void)
 Returns the index of the slice at the center of the cross-correlation map. More...
 
unsigned int GetSliceN (void)
 Returns the number of slices in one cross-correlation map. More...
 
double GetSliceTimeDelay (const unsigned int aSliceIndex)
 Returns the time delay at the center of a slice [s]. More...
 
unsigned int GetSliceTimeDelayBinsN (void)
 Returns the number of time-delay bins in one slice. More...
 
unsigned int GetSliceUnusedTimeDelayBinsN (void)
 Returns the number of unsed time-delay bins after slicing. More...
 
bool SaveEvents (const UChar_t aAnalysisIndex=0)
 Save events in the TTree. More...
 
Constructors and destructors
 OxEvent (const string aOptionFile, const unsigned int aGpsRef)
 Constructor of the OxEvent class. More...
 
virtual ~OxEvent (void)
 Destructor of the OxEvent class. More...
 
- Public Member Functions inherited from OxCorr
double GetTimeDelay (const unsigned int aTimeDelayBinIndex)
 Returns the time delay of a bin in the cross-correlation map [s]. More...
 
unsigned int GetTimeDelayBinIndex (const double aTimeDelay)
 Returns the time-delay bin index for a given time delay. More...
 
unsigned int GetTimeDelayBinsN (void)
 Returns the number of time-delay bins. More...
 
double GetTimeDelayResolution (void)
 Returns the time-delay resolution [s]. More...
 
double GetXi (const unsigned int aQindex, const unsigned int aFrequencyIndex, const unsigned int aTimeDelayIndex)
 Returns the cross-correlation coefficient \(\xi(Q,f,dt)\). More...
 
bool Process (void)
 Runs the cross-correlation analysis. More...
 
bool Transform (const unsigned int aDetectorIndex)
 Transforms the Omicron spectrogram in the Fourier domain. More...
 
 OxCorr (const string aOptionFile, const unsigned int aGpsRef)
 Constructor of the OxCorr class. More...
 
virtual ~OxCorr (void)
 Destructor of the OxCorr class. More...
 
- Public Member Functions inherited from OxOmicron
bool CreateOutputFile (const unsigned int aTime)
 Creates an ouput file. More...
 
unsigned int GetAnalysisDuration (void)
 Returns the Ox analysis duration [s]. More...
 
Long64_t GetEntries (const unsigned int aDetectorIndex)
 Returns the number of Omicron results. More...
 
int GetEntry (const unsigned int aDetectorIndex, Long64_t aEntry)
 Loads an entry in the Omicron results. More...
 
double GetLightTravelTime (void)
 Returns the light travel time between the two detectors. More...
 
string GetOmicronChannelName (const unsigned int aDetectorIndex)
 Returns the Omicron channel name. More...
 
unsigned int GetOmicronChunkDuration (void)
 Returns the Omicron chunk duration [s]. More...
 
unsigned int GetOmicronEndTime (const unsigned int aDetectorIndex)
 Returns the end time of the current Omicron entry. More...
 
unsigned int GetOmicronOverlapDuration (void)
 Returns the Omicron overlap duration [s]. More...
 
unsigned int GetOmicronQ (const unsigned int aQindex)
 Returns Q value of a given Q plane. More...
 
unsigned int GetOmicronQN (void)
 Returns the number of Omicron Q planes. More...
 
unsigned int GetOmicronStartTime (const unsigned int aDetectorIndex)
 Returns the start time of the current Omicron entry. More...
 
bool Process (const unsigned int aDetectorIndex, const unsigned int aTimeStart, const bool aResetPsd)
 Runs the Omicron analysis. More...
 
 OxOmicron (const string aOptionFile, const unsigned int aGpsRef)
 Constructor of the OxOmicron class. More...
 
virtual ~OxOmicron (void)
 Destructor of the OxOmicron class. More...
 
- Public Member Functions inherited from OxInit
void AttachTree (TTree *aTree)
 Attach a TTree to the output file. More...
 
void CloseOutputFile (void)
 Closes the current ouput file (if any). More...
 
bool CreateOutputFile (const unsigned int aTime)
 Creates an ouput file. More...
 
string GetName (void)
 Returns the object name. More...
 
string GetOutputDirectory (void)
 Returns the path to the output directory. More...
 
bool GetStatus (void)
 Returns the class status. More...
 
 OxInit (const string aOptionFile)
 Constructor of the OxInit class. More...
 
virtual ~OxInit (void)
 Destructor of the OxInit class. More...
 

Private Member Functions

void ComputeEvent (const unsigned int aSliceIndex)
 Computes the event parameters. More...
 
void MakeSlices (void)
 Slices a cross-correlation map. More...
 

Private Attributes

UChar_t oxe_aindex
 Event tree: analysis index. More...
 
double oxe_dt
 Event tree: time delay [s]. More...
 
double * oxe_dt_var_q
 Event tree: time-delay variance [s^2] / Q. More...
 
double * oxe_frequency_q
 Event tree: frequency [Hz] / Q. More...
 
unsigned int * oxe_n_sveto_q
 Event tree: number of frequency rows rejected by the slice veto / Q. More...
 
Long64_t oxe_oxo_0
 Event tree: Omicron TTree entry (det 0). More...
 
Long64_t oxe_oxo_1
 Event tree: Omicron TTree entry (det 1). More...
 
unsigned int oxe_sindex
 Event tree: slice index. More...
 
TTree * oxe_tree
 Event tree. Do not delete as it is owned by the output file. More...
 
TChain * oxe_tree_read
 TChain: Ox events (read-only). More...
 
double * oxe_Xi_q
 Event tree: \(\Xi\) /Q. More...
 
int ** slice_ctn_binmax
 Slice container: time-delay bin where xi is max / f row / Q. More...
 
unsigned int * slice_firstbin
 First time-delay bin index of each slice. More...
 
unsigned int slice_n
 Number of slices in a cross-correlation map. More...
 
unsigned int slice_nbins
 Number of time-delay bins in a slice. More...
 

Additional Inherited Members

- Protected Member Functions inherited from OxInit
void AddOptions (void)
 Adds options from the text file. More...
 
void DefineOption (const string aTag, const string aKey, const double aDefaultValue, const unsigned int aSize)
 Defines an option (double). More...
 
void DefineOption (const string aTag, const string aKey, const int aDefaultValue, const unsigned int aSize)
 Defines an option (integer). More...
 
void DefineOption (const string aTag, const string aKey, const string aDefaultValue, const unsigned int aSize)
 Defines an option (string). More...
 
void DefineOption (const string aTag, const string aKey, const unsigned int aDefaultValue, const unsigned int aSize)
 Defines an option (unsigned integer). More...
 
void OverloadOption (const string aTag, const string aKey)
 Overloads an option using the text file as a reference. More...
 
- Protected Attributes inherited from OxOmicron
Omicronomicron [2]
 Omicron objects (detector 1 and 2) - do not delete. More...
 
TTree * oxo_tree [2]
 TTree: Omicron results. Do not delete as it is owned by the output file. More...
 
TChain * oxo_tree_read [2]
 TChain: Omicron results (read-only). More...
 
- Protected Attributes inherited from OxInit
GwollumOptions * ox_opt
 List of options. More...
 
bool status
 Class status. More...
 

Detailed Description

Ox event management.

An Ox event is defined as a time-delay slice in a cross-correlation map built with OxCorr. Each event is parameterized by (see ComputeEvent()):

  • An analysis index: 8-bit integer to associate an event to a given analysis.
  • A slice index: see MakeSlices() for a definition.
  • A time delay: this is the time shift between detector 1 and 2 which maximizes cross-correlation.
  • For each Q plane, the integral of \(xi\), noted \(Xi\), over all frequencies, calculated for the event time delay.
  • For each Q plane, an average frequency.
  • For each Q plane, the number of frequency rows which are rejected by the slice veto.
  • For each Q plane, the time-delay variance around the event time delay.

This class can be used in two different ways:

  1. To build a list of Ox events. After creating an Ox file with CreateOutputFile(), use the OxCorr engine to build cross-correlation maps. Then call SaveEvents to record the events in the file.
  2. To process a list of Ox events. In this case, a list of Ox files (generated in 1.) must be provided in the option file (OXI/OXFILELIST). Then, call the events with GetEntry().

Constructor & Destructor Documentation

◆ OxEvent()

OxEvent::OxEvent ( const string  aOptionFile,
const unsigned int  aGpsRef 
)

Constructor of the OxEvent class.

An option file must be provided to configure the Ox event management:

  • OXE SLICEDELTA Size of slices. Two time delays [s] are expected: the calibration time uncertainty, noted \(\delta_\mathrm{cal}\), and the slice veto time delay \(\delta_\mathrm{veto}\).

The cross-correlation map is sliced with MakeSlices().

Parameters
[in]aOptionFilePath to the option file.
[in]aGpsRefReference GPS time to intiate Omicron objects.

◆ ~OxEvent()

OxEvent::~OxEvent ( void  )
virtual

Destructor of the OxEvent class.

Member Function Documentation

◆ ComputeEvent()

void OxEvent::ComputeEvent ( const unsigned int  aSliceIndex)
private

Computes the event parameters.

TBC.

◆ CreateOutputFile()

bool OxEvent::CreateOutputFile ( const unsigned int  aTime)

Creates an ouput file.

A new ROOT file is opened in the output directory to save the Ox analysis results (OxInit::CreateOutputFile()). The Ox event TTree "oxe" is initialized in the ouput file.

See also
OxInit::CreateOutputFile() for the file management.
Parameters
[in]aTimeTime used to name the ouput file.

◆ GetEntries()

Long64_t OxEvent::GetEntries ( void  )
inline

Returns the number of events (read-only).

◆ GetEntry()

int OxEvent::GetEntry ( const Long64_t  aEntry)

Loads an entry in the event TTree (read-only).

The corresponding Omicron result is also loaded.

Returns
The function returns the number of bytes read from the input buffer. If entry does not exist the function returns 0. If an I/O error occurs, the function returns -1.
Parameters
[in]aEntryEntry number.

◆ GetEventAnalysisIndex()

UChar_t OxEvent::GetEventAnalysisIndex ( void  )
inline

Returns the analysis index for the current event.

◆ GetEventFrequency()

double OxEvent::GetEventFrequency ( void  )

Returns the frequency for the current event [Hz].

The event frequency is averaged over all the Omicron Q planes.

◆ GetEventRank()

double OxEvent::GetEventRank ( void  )

Returns the rank for the current event.

The event rank is the average \(\Xi\) over all the Omicron Q planes.

◆ GetEventSliceIndex()

unsigned int OxEvent::GetEventSliceIndex ( void  )
inline

Returns the slice index for the current event.

◆ GetSliceHalfDuration()

double OxEvent::GetSliceHalfDuration ( void  )
inline

Returns the half duration of one slice [s].

See also
MakeSlices()

◆ GetSliceIndexCenter()

unsigned int OxEvent::GetSliceIndexCenter ( void  )
inline

Returns the index of the slice at the center of the cross-correlation map.

This is the slice indexed by \(s= n_s = (N_s-1)/2\).

◆ GetSliceN()

unsigned int OxEvent::GetSliceN ( void  )
inline

Returns the number of slices in one cross-correlation map.

See also
MakeSlices()

◆ GetSliceTimeDelay()

double OxEvent::GetSliceTimeDelay ( const unsigned int  aSliceIndex)
inline

Returns the time delay at the center of a slice [s].

Parameters
[in]aSliceIndexSlice index.
Precondition
The slice index must be valid.
See also
MakeSlices()

◆ GetSliceTimeDelayBinsN()

unsigned int OxEvent::GetSliceTimeDelayBinsN ( void  )
inline

Returns the number of time-delay bins in one slice.

See also
MakeSlices()

◆ GetSliceUnusedTimeDelayBinsN()

unsigned int OxEvent::GetSliceUnusedTimeDelayBinsN ( void  )
inline

Returns the number of unsed time-delay bins after slicing.

Note
It is always an even number (by construction).
See also
MakeSlices()

◆ MakeSlices()

void OxEvent::MakeSlices ( void  )
private

Slices a cross-correlation map.

A slice in the cross-correlation map is a time-delay window the size of which is \(\pm (\delta_\mathrm{det} + \delta_\mathrm{cal} + \delta_\mathrm{veto})\). The maximum light travel time between the two detectors, \(\delta_\mathrm{det}\), is automatically determined with the GetLightTravelTime() function from GWOLLUM. The calibration time uncertainty is noted \(\delta_\mathrm{cal}\). An additional time delay \(\delta_\mathrm{veto}\) can be provided. The slice width includes all correlated signals between the 2 detectors.

A sequence of slices is constructed with the following conditions:

  • Slices are contiguous.
  • The number of slices is maximal given the size of the cross-correlation map.
  • There must be an odd number of slices in a map: \(N_s = 2n_s+1\).
  • Slices in a map are indexed by \(s\) running from \(0\) to \(N_s-1\). This is called the "slice index".
  • The slice \(s = n_s = (N_s-1)/2\) is centered at the center of the cross-correlation map.
Note
Some part of the map may be left unused at both edges of the cross-correlation map.

The following image shows an example of slicing. There is \(N_s=7\) slices. The slices are positioned symmetrically with respect to the center of the map.

The cross-correlation map is sliced

◆ SaveEvents()

bool OxEvent::SaveEvents ( const UChar_t  aAnalysisIndex = 0)

Save events in the TTree.

After the cross-correlation analysis is finished (OxCorr::Process()), call this function to save all the events in the output TTree. Each slice is processed with ComputeEvent() to calculate all the event parameters.

An event can be tagged with an analysis index. This option is useful to manage multiple sets of events.

Note
This function can only be called after an ouput file is created with CreateOutputFile().
Parameters
[in]aAnalysisIndexAnalysis index.

Member Data Documentation

◆ oxe_aindex

UChar_t OxEvent::oxe_aindex
private

Event tree: analysis index.

◆ oxe_dt

double OxEvent::oxe_dt
private

Event tree: time delay [s].

◆ oxe_dt_var_q

double* OxEvent::oxe_dt_var_q
private

Event tree: time-delay variance [s^2] / Q.

◆ oxe_frequency_q

double* OxEvent::oxe_frequency_q
private

Event tree: frequency [Hz] / Q.

◆ oxe_n_sveto_q

unsigned int* OxEvent::oxe_n_sveto_q
private

Event tree: number of frequency rows rejected by the slice veto / Q.

◆ oxe_oxo_0

Long64_t OxEvent::oxe_oxo_0
private

Event tree: Omicron TTree entry (det 0).

◆ oxe_oxo_1

Long64_t OxEvent::oxe_oxo_1
private

Event tree: Omicron TTree entry (det 1).

◆ oxe_sindex

unsigned int OxEvent::oxe_sindex
private

Event tree: slice index.

◆ oxe_tree

TTree* OxEvent::oxe_tree
private

Event tree. Do not delete as it is owned by the output file.

◆ oxe_tree_read

TChain* OxEvent::oxe_tree_read
private

TChain: Ox events (read-only).

◆ oxe_Xi_q

double* OxEvent::oxe_Xi_q
private

Event tree: \(\Xi\) /Q.

◆ slice_ctn_binmax

int** OxEvent::slice_ctn_binmax
private

Slice container: time-delay bin where xi is max / f row / Q.

◆ slice_firstbin

unsigned int* OxEvent::slice_firstbin
private

First time-delay bin index of each slice.

◆ slice_n

unsigned int OxEvent::slice_n
private

Number of slices in a cross-correlation map.

◆ slice_nbins

unsigned int OxEvent::slice_nbins
private

Number of time-delay bins in a slice.


The documentation for this class was generated from the following files: