Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_Coherence.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 "config.hh"
7 #include "network.hh"
8 #include "wavearray.hh"
9 #include "TString.h"
10 #include "TObjArray.h"
11 #include "TObjString.h"
12 #include "TRandom.h"
13 #include "TComplex.h"
14 #include "TMath.h"
15 #include "mdc.hh"
16 #include <vector>
17 
18 //!COHERENCE
19 
20 void
22 
23 // This plugin implements the built-in coherence stage (only 2G), user can provide its own coherence stage implementation
24 
25  cout << endl;
26  cout << "-----> plugins/CWB_Plugin_Coherence.C" << endl;
27  cout << "ifo " << ifo.Data() << endl;
28  cout << "type " << type << endl;
29  cout << endl;
30 
31  if(type==CWB_PLUGIN_CONFIG) {
32  cfg->cohPlugin=true; // disable built-in coherence stage
33  }
34 
35  if(type==CWB_PLUGIN_ICOHERENCE) {
36 
37  cout << "type==CWB_PLUGIN_ICOHERENCE" << endl;
38 
39  // import ifactor
40  int gIFACTOR=-1; IMPORT(int,gIFACTOR)
41  cout << "-----> CWB_Plugin_Sim4.C -> " << " gIFACTOR : " << gIFACTOR << endl;
42 
43  int rateANA=cfg->inRate>>cfg->levelR;
44  int nIFO = net->ifoListSize();
45  detector* pD[NIFO_MAX]; // pointers to detectors
46  for(int n=0;n<nIFO;n++) pD[n] = net->getifo(n);
47  double mTau=net->getDelay(const_cast<char*>("MAX")); // maximum time delay
48 
49  int upN = rateANA/1024; if(upN<1) upN=1; // calculate upsample factor
50  int nRES = net->wdmListSize(); // number of resolution levels
51 
52  double Eo;
53  netcluster* pwc;
54 
55  for(int i=0; i<nRES; i++) { // loop over TF resolutions
56  // print level infos
57  int level=cfg->l_high-i;
58  int layers = level>0 ? 1<<level : 0;
59  int rate = rateANA>>level;
60  cout << "level : " << level << "\t rate(hz) : " << rate
61  << "\t layers : " << layers << "\t df(hz) : " << rateANA/2./double(1<<level)
62  << "\t dt(ms) : " << 1000./rate << endl;
63 
64  // produce TF maps with max over the sky energy
65  for(int n=0; n<nIFO; n++) {
66  WDM<double>* pwdm = net->wdmList[i];
67  wavearray<double>* hot = net->getifo(n)->getHoT();
68  net->getifo(n)->getTFmap()->maxEnergy(*hot,*pwdm,mTau,upN);
69  //if(singleDetector) {
70  // *(net->getifo(1)->getTFmap()) = *(net->getifo(0)->getTFmap());
71  // break;
72  //}
73  // restore the frequency boundaries changed by the maxEnergy call
74  net->getifo(n)->getTFmap()->setlow(cfg->fLow);
75  net->getifo(n)->getTFmap()->sethigh(cfg->fHigh);
76  }
77 
78  Eo = net->THRESHOLD(cfg->bpp); // threshold on pixel energy
79  cout<<"thresholds in units of noise variance: Eo="<<Eo<<" Emax="<<Eo*2<<endl;
80 
81  double TL = net->setVeto(cfg->iwindow);
82  cout<<"live time in zero lag: "<<TL<<endl; // set veto array
83  if(TL <= 0.) {cout<<"livetime is zero : exit"<<endl;exit(1);} // exit if live time is zero
84 
85  // select pixels
86  if(cfg->simulation) {cout<<"ifactor|clusters|pixels ";cout.flush();}
87  else {cout<<"lag|clusters|pixels "; cout.flush();}
88  int csize_tot=0;int psize_tot=0;
89  for(int j=0; j<(int)net->nLag; j++) {
90 
91  net->getNetworkPixels(j,Eo,Eo*2);
92  net->cluster(1,1);
93  pwc = net->getwc(j);
94 
95  // store cluster into temporary job file
96  int cycle = cfg->simulation ? gIFACTOR : Long_t(pwc->shift);
97  pwc->write(jfile,"coherence","clusters",0,cycle);
98  pwc->write(jfile,"coherence","clusters",-1,cycle);
99  cout<<cycle<<"|"<<pwc->csize()<<"|"<<pwc->size()<<" ";cout.flush();
100  csize_tot+=pwc->csize(); psize_tot+=pwc->size();
101 
102  pwc->clear();
103  }
104  }
105  }
106 
107  return;
108 }
void sethigh(double f)
Definition: wseries.hh:114
detector * getifo(size_t n)
param: detector index
Definition: network.hh:418
CWB::config * cfg
Definition: TestCWB_Plugin.C:5
double iwindow
Definition: config.hh:182
size_t write(const char *file, int app=0)
Definition: netcluster.cc:2989
size_t nLag
Definition: network.hh:555
long getNetworkPixels(int LAG, double Eo, double DD=1., TH1F *hist=NULL)
Definition: network.cc:60
int n
Definition: cwb_net.C:10
cout<< "skymap size : "<< L<< endl;for(int l=0;l< L;l++) sm.set(l, l);sm > const_cast< char * >("skymap.dat")
TString("c")
bool cohPlugin
Definition: config.hh:349
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 fLow
Definition: config.hh:122
netcluster * pwc
Definition: cwb_job_obj.C:20
void setlow(double f)
Definition: wseries.hh:107
int layers
int j
Definition: cwb_net.C:10
i drho i
size_t cluster(int kt, int kf)
param: time gap in pixels return: number of reconstructed clusters
Definition: network.hh:303
size_t wdmListSize()
Definition: network.hh:428
char ifo[NIFO_MAX][8]
network ** net
NOISE_MDC_SIMULATION.
size_t ifoListSize()
Definition: network.hh:413
#define nIFO
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *net, WSeries< double > *x, TString ifo, int type)
COHERENCE.
jfile
Definition: cwb_job_obj.C:25
double setVeto(double=5.)
param: time window around injections
Definition: network.cc:3456
double shift
Definition: netcluster.hh:364
int simulation
Definition: config.hh:181
int l_high
Definition: config.hh:138
wavearray< double > hot[2]
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
double maxEnergy(wavearray< DataType_t > &ts, Wavelet &w, double=0, int=1, int=0, TH1F *=NULL)
Definition: wseries.cc:486
i() int(T_cor *100))
wavearray< double > * getHoT()
param: no parameters
Definition: detector.hh:157
const int NIFO_MAX
Definition: wat.hh:4
double getDelay(const char *c="")
Definition: network.cc:2787
std::vector< WDM< double > * > wdmList
Definition: network.hh:599
int gIFACTOR
double THRESHOLD(double bpp)
param: selected fraction of LTF pixels assuming Gaussian noise
Definition: network.cc:2584
size_t rateANA
Definition: test_config1.C:21
size_t size()
Definition: netcluster.hh:129
double fHigh
Definition: config.hh:123
WSeries< double > * getTFmap()
param: no parameters
Definition: detector.hh:161
size_t csize()
Definition: netcluster.hh:133
double mTau
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:421
const int nRES
Definition: revMonster.cc:7
double bpp
Definition: config.hh:115
size_t inRate
Definition: config.hh:114
int levelR
Definition: config.hh:134
detector ** pD
exit(0)
void clear()
Definition: netcluster.hh:409