Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_Gating.C
Go to the documentation of this file.
1 #define XIFO 4
2 
3 #pragma GCC system_header
4 
5 #include "cwb.hh"
6 #include "cwb2G.hh"
7 #include "config.hh"
8 #include "network.hh"
9 #include "wavearray.hh"
10 #include "TString.h"
11 #include "TObjArray.h"
12 #include "TObjString.h"
13 #include "TRandom.h"
14 #include "Toolbox.hh"
15 
16 //!DATA_CONDITIONING
17 
18 // This plugin implements the time gating in the pixel's selection stage.
19 // Gating is a veto of the pixels belonging to a time interval.
20 // This plugin is used to exclude from the analysis the pixels
21 // in a time interval where there is a huge glitch.
22 // Huge glitches are harmful because they affect the correct estimation
23 // of the selection threshold and could mask the nearby events at lower SNR.
24 // Warning : events with high SNR can be rejected by this procedure (see SETHR)
25 
26 #define SETHR 1000000 // Is the threshold which define the time slices to be cutted
27  // Warning : this value depends on the frequency interval [fHigh:fLow]
28 
29 #define TWIN 0.5 // Is the window (sec) used to integrate the energies in time
30  // TWIN must be a multiple of the greatest time resolution used in the analysis
31 
32 #define TEDG 1.5 // Is the time window (sec) to be cutted when the time integrated energy is > SETHR
33 
34 void
36 
37 // this plugin is called in cwb2G.cc to the whitened data at the beginning of the COHERENT stage
38 
39  cout << endl;
40  cout << "-----> CWB_Plugin_Gating.C" << endl;
41  cout << "ifo " << ifo.Data() << endl;
42  cout << "type " << type << endl;
43  cout << endl;
44 
45  if(type==CWB_PLUGIN_ICOHERENCE) {
46 
47  // get data range
48  cwb2G* gCWB2G; IMPORT(cwb2G*,gCWB2G)
49  double Tb = gCWB2G->GetSegBegin()-cfg->segEdge;
50  double Te = gCWB2G->GetSegEnd()+cfg->segEdge;
51 
52 
53  //cout << "CWB_Plugin_Gating.C - " << "Tb : " << int(Tb) << " Te : " << int(Te) << endl;
54 
55  int nIFO = (gCWB2G->IsSingleDetector()) ? 1 : net->ifoListSize(); // number of detectors
56  int size = net->getifo(0)->getHoT()->size(); // number of time bins
57  double rate = net->getifo(0)->getHoT()->rate(); // rate
58  double dt = 1./rate; // time bin resolution (sec)
59 
60  //cout << "-----> CWB_Plugin_Gating.C - size : " << size << "\t rate : " << rate << "\t dt : " << dt << " length(sec) " << size*dt << endl;
61 
63  for(int n=0; n<nIFO; n++) hot[n] = net->getifo(n)->getHoT(); // whitened data
64 
65  // A new array 'SE' is obtained from the arrays 'hot whitened data'
66  // It contains the time integrated energies over the time window TWIN
67  // A time sample is marked as to be vetoed when the energy SE>SETHR
68  // When a time sample 'i' is vetoed also the nearby M time size are vetoed [i-M:i+M]
69  // where M is the number of samples contained in TEDG sec
70 
71  int X = int(cfg->segEdge/dt+0.5); // samples in segEdge : this the is scratch time
72  int M = int(TEDG/dt)+1; // samples in TEDG sec
73  int N = int(TWIN/dt); // number of time samples in TWIN
74  if(N<1) N=1;
75  vector<waveSegment> glist; // gating segment list
76  wavearray<double> SE(size);
77  wavearray<int> sVeto(size); // array which contains samples to be vetoed
78 
79  int gsize=0;
80  for(int n=0; n<nIFO; n++) { // loop over detectors
81  SE=0;
82  for(int i=N-1;i<size;i++) for(int j=0;j<N;j++) SE[i]+=pow(hot[n]->data[i-j],2);
83  sVeto=0;
84  for(int i=0;i<size;i++) {
85  if(SE[i]>SETHR) {
86  int is = i-M>X ? i-M : X;
87  int es = i+M<(size-X) ? i+M : size-X;
88  for(int k=is;k<es;k++) sVeto[k]=1;
89  }
90  }
91  // array derivative -> gating_start=1, gating_stop=-1
92  sVeto[0]=0;sVeto[size-1]=0;
93  for(int i=0;i<size-1;i++) sVeto[i]=sVeto[i+1]-sVeto[i];
94 
95  // build the gating segment list
96  waveSegment gseg; // gating segment
97  gseg.index = 0;
98  for(int i=0;i<size;i++) {
99  if(sVeto[i]== 1) gseg.start = Tb+int(i*dt); // round to nearest lower integer
100  if(sVeto[i]==-1) {
101  gseg.stop = Tb+int(i*dt+0.5); // round to nearest higher integer
102  gseg.index++;
103  glist.push_back(gseg);
104  }
105  }
106 
107  if(glist.size()>0) {
108  cout << endl;
109  cout.precision(10);
110  for(int j=gsize;j<glist.size();j++) {
111  cout << j << " -----> CWB_Plugin_Gating.C - " << net->getifo(n)->Name << " gating time : start="
112  << glist[j].start << " stop=" << glist[j].stop << " len=" << glist[j].stop-glist[j].start << endl;
113  }
114  cout << endl;
115  gsize=glist.size();
116  }
117  }
118 
119  double gating_time = 0; // total gating time
120  if(glist.size()) {
121  glist = CWB::Toolbox::unionSegments(glist); // Join & sort a waveSegment list
122  gating_time = CWB::Toolbox::getTimeSegList(glist); // get total gating time
123  glist = CWB::Toolbox::invertSegments(glist); // Invert waveSegment list
124  net->segList = CWB::Toolbox::mergeSegLists(glist,net->segList); // merge gating veto list to the net->segList
125  }
126 
127  cout<< "-----> CWB_Plugin_Gating.C - gating time : " << gating_time << endl << endl;
128 
129  // add infos to history
130  char info[256];
131  sprintf(info,"-GT:%g",gating_time);
132  gCWB2G->PrintAnalysisInfo(CWB_STAGE_COHERENCE,"cwb2G::Coherence",info,false);
133  }
134 
135  return;
136 }
detector * getifo(size_t n)
param: detector index
Definition: network.hh:418
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
virtual size_t size() const
Definition: wavearray.hh:127
double start
Definition: network.hh:37
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *net, WSeries< double > *x, TString ifo, int type)
COHERENCE.
virtual void rate(double r)
Definition: wavearray.hh:123
int n
Definition: cwb_net.C:10
TString("c")
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
double Te
Long_t size
#define M
Definition: UniqSLagsList.C:3
static double getTimeSegList(vector< waveSegment > list)
Definition: Toolbox.cc:592
int j
Definition: cwb_net.C:10
i drho i
static vector< waveSegment > invertSegments(vector< waveSegment > &ilist)
Definition: Toolbox.cc:229
#define N
#define TEDG
char ifo[NIFO_MAX][8]
Definition: cwb2G.hh:15
network ** net
NOISE_MDC_SIMULATION.
size_t ifoListSize()
Definition: network.hh:413
static vector< waveSegment > unionSegments(vector< waveSegment > &ilist)
Definition: Toolbox.cc:101
double segEdge
Definition: config.hh:146
#define nIFO
#define SETHR
DATA_CONDITIONING.
double GetSegEnd()
Definition: cwb.hh:171
jfile
Definition: cwb_job_obj.C:25
void PrintAnalysisInfo(CWB_STAGE stage, TString comment, TString info, bool out=true, bool log=true)
Definition: cwb.cc:2110
wavearray< double > hot[2]
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
i() int(T_cor *100))
wavearray< double > * getHoT()
param: no parameters
Definition: detector.hh:157
const int NIFO_MAX
Definition: wat.hh:4
double GetSegBegin()
Definition: cwb.hh:170
#define TWIN
bool IsSingleDetector()
Definition: cwb.hh:173
int k
std::vector< waveSegment > segList
Definition: network.hh:598
double dt
int index
Definition: network.hh:36
double Tb
char Name[16]
Definition: detector.hh:309
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
static vector< waveSegment > mergeSegLists(vector< waveSegment > &ilist1, vector< waveSegment > &ilist2)
Definition: Toolbox.cc:332
double stop
Definition: network.hh:38