Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CWB_Plugin_DataConditioning.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 "regression.hh"
17 #include <vector>
18 
19 
20 #define REGRESSION_FILTER_LENGTH 8
21 #define REGRESSION_MATRIX_FRACTION 0.95
22 #define REGRESSION_SOLVE_EIGEN_THR 0.
23 #define REGRESSION_SOLVE_EIGEN_NUM 10
24 #define REGRESSION_SOLVE_REGULATOR 'h'
25 #define REGRESSION_APPLY_THR 0.8
26 
27 #define LPR_LEVELD 8
28 #define LPR_LEVELF 6
29 #define LPR_TLPR 120.
30 
31 #define USE_REGRESSION
32 //#define USE_LPR_FILTER
33 
34 void
36 //!DATA_CONDITIONING
37 // This plugin shows how to implement a custom Data Conditioning (only 2G)
38 
39  cout << endl;
40  cout << "-----> CWB_Plugin_DataConditioning.C" << endl;
41  cout << "ifo " << ifo.Data() << endl;
42  cout << "type " << type << endl;
43  cout << endl;
44 
45  if(type==CWB_PLUGIN_CONFIG) {
46  cfg->dcPlugin=true; // disable built in dc
47  }
48 
49  if(type==CWB_PLUGIN_DATA_CONDITIONING) {
50 
51  if(TString(cfg->analysis)!="2G") {
52  cout << "CWB_Plugin_DAtaConditioning.C -> implemented only for 2G" << endl;
53  gSystem->Exit(1);
54  }
55 
58  // import global variables
59  size_t gRATEANA; IMPORT(size_t,gRATEANA)
60  void* gMDC=NULL; IMPORT(void*,gMDC)
61  WSeries<double>* xM = (WSeries<double>*)gMDC;
62 #ifdef USE_LPR_FILTER // used in 1G
63  Meyer<double> S(1024,2); // set wavelet for production
64 #endif
65 
66  // get ifo id
67  int id=-1;
68  for(int n=0;n<cfg->nIFO;n++) if(ifo==net->ifoName[n]) {id=n;break;}
69  if(id<0) {cout << "Plugin : Error - bad ifo id" << endl; gSystem->Exit(1);}
70 
71  detector* pD = net->getifo(id);
72  WSeries<double>* pTF = pD->getTFmap();
74 
75  int layers = 0;
76  for(int level=cfg->l_high; level>=cfg->l_low; level--) {
77  for(int j=0;j<net->wdmMRA.nRes;j++)
78  if(layers<net->wdmMRA.layers[j]) layers=net->wdmMRA.layers[j];
79  }
80  WDM<double> WDMwhite(layers,layers,6,10); // set whitening WDM
81  layers = gRATEANA/8;
82  WDM<double> WDMlpr(layers,layers,6,10); // set LPE filter WDM
83 
84 #ifdef USE_LPR_FILTER // used in 1G
85  pTF->Forward(*hot,S,LPR_LEVELD);
86  pTF->lprFilter(2,0,LPR_TLPR,4.);
87  pTF->Inverse();
88  *hot = *pTF; // conditioned TS
89 #endif
90 
91 #ifdef USE_REGRESSION // used in 2G
92  // regression
93  pTF->Forward(*hot,WDMlpr);
94  regression rr(*pTF,const_cast<char*>("target"),1.,cfg->fHigh);
95  rr.add(*hot,const_cast<char*>("target"));
100  *hot = rr.getClean();
101 #endif
102 
103  // whitening
104  pTF->Forward(*hot,WDMwhite);
105  pTF->setlow(cfg->fLow);
106  pTF->sethigh(cfg->fHigh);
107  pD->white(cfg->whiteWindow,0,cfg->segEdge,cfg->whiteStride); // calculate noise rms
108  pD->nRMS.bandpass(16.,0.,1); // high pass filtering at 16Hz
109  pTF->white(pD->nRMS,1); // whiten 0 phase WSeries
110  pTF->white(pD->nRMS,-1); // whiten 90 phase WSeries
111  if(cfg->simulation) { // estimated whitened MDC parms
112  wM.Forward(*xM,WDMwhite);
113  wM.setlow(cfg->fLow);
114  wM.sethigh(cfg->fHigh);
115  // zero f<fLow to avoid whitening issues when psd noise is not well defined for f<fLow
116  int layers = wM.maxLayer();
117  for(int j=0;j<layers;j++) if(wM.frequency(j)<cfg->fLow) {
118  double layer = j+0.01; // -epsilon select 0 layer for 90 phase
119  wM.getLayer(y, layer);y=0;wM.putLayer(y, layer); // 0 phase
120  wM.getLayer(y,-layer);y=0;wM.putLayer(y,-layer); // 90 phase
121  }
122  // compute mdc snr and save whiten waveforms
123  pD->setsim(wM,net->getmdcTime(),cfg->iwindow/2.,cfg->segEdge,true);
124  }
125  WSeries<double> wtmp = *pTF;
126  pTF->Inverse();
127  wtmp.Inverse(-2);
128  *hot = *pTF;
129  *hot += wtmp;
130  *hot *= 0.5;
131 
132 #ifdef USE_LPR_FILTER // used in 1G
133  pTF->Forward(*hot,S,LPR_LEVELF);
134  pTF->lprFilter(2,0,LPR_TLPR,4.);
135  pTF->Inverse();
136  *hot = *pTF; // conditioned TS
137 #endif
138 
139  }
140  return;
141 }
monster wdmMRA
list of pixel pointers for MRA
Definition: network.hh:633
std::vector< char * > ifoName
Definition: network.hh:591
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
char analysis[8]
Definition: config.hh:99
double iwindow
Definition: config.hh:182
void white(double dT=0, int wtype=1, double offset=0., double stride=0.)
what it does: see algorithm description in wseries.hh
Definition: detector.hh:240
void setMatrix(double edge=0., double f=1.)
Definition: regression.cc:407
std::vector< double > * getmdcTime()
Definition: network.hh:404
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")
void bandpass(wavearray< DataType_t > &ts, double flow, double fhigh, int n=-1)
Definition: wseries.cc:295
#define REGRESSION_SOLVE_REGULATOR
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 whiteWindow
Definition: config.hh:171
double fLow
Definition: config.hh:122
bool dcPlugin
Definition: config.hh:348
#define LPR_TLPR
void setlow(double f)
Definition: wseries.hh:107
int layers
#define REGRESSION_APPLY_THR
size_t add(WSeries< double > &target, char *name, double fL=0., double fH=0.)
Definition: regression.cc:73
int j
Definition: cwb_net.C:10
char ifo[NIFO_MAX][8]
double Edge
Definition: network.hh:560
network ** net
NOISE_MDC_SIMULATION.
double segEdge
Definition: config.hh:146
jfile
Definition: cwb_job_obj.C:25
int simulation
Definition: config.hh:181
int l_high
Definition: config.hh:138
#define REGRESSION_SOLVE_EIGEN_THR
wavearray< double > hot[2]
#define IMPORT(TYPE, VAR)
Definition: cwb.hh:51
WSeries< double > nRMS
Definition: detector.hh:340
#define REGRESSION_FILTER_LENGTH
void apply(double threshold=0., char c='a')
Definition: regression.cc:691
WSeries< double > pTF[nRES]
Definition: revMonster.cc:8
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
CWB::mdc * gMDC
Definition: cwb_mkfad.C:44
int gRATEANA
void solve(double th, int nE=0, char c='s')
Definition: regression.cc:592
size_t setFilter(size_t)
Definition: regression.cc:258
wavearray< double > getClean()
Definition: regression.hh:117
WSeries< double > wM
Definition: cwb_job_obj.C:23
int l_low
Definition: config.hh:137
double whiteStride
Definition: config.hh:172
double fHigh
Definition: config.hh:123
WSeries< double > * getTFmap()
param: no parameters
Definition: detector.hh:161
Definition: Meyer.hh:18
int * layers
Definition: monster.hh:103
#define REGRESSION_SOLVE_EIGEN_NUM
virtual void lprFilter(double, int=0, double=0., double=0.)
Definition: wseries.cc:1108
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228
Meyer< double > S(1024, 2)
int nRes
Definition: monster.hh:99
int nIFO
Definition: config.hh:102
#define LPR_LEVELD
#define REGRESSION_MATRIX_FRACTION
#define LPR_LEVELF
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:273
wavearray< double > y
Definition: Test10.C:31
void putLayer(wavearray< DataType_t > &, double n)
param: n - layer number
Definition: wseries.cc:201
size_t setsim(WSeries< double > &, std::vector< double > *, double=5., double=8., bool saveWF=false)
Definition: detector.cc:1348
double frequency(int l)
Definition: wseries.cc:99
int maxLayer()
Definition: wseries.hh:121
detector ** pD
virtual WSeries< double > white(double, int, double=0., double=0.)
what it does: each wavelet layer is devided into k intervals.
Definition: wseries.cc:1128
void CWB_Plugin(TFile *jfile, CWB::config *cfg, network *net, WSeries< double > *x, TString ifo, int type)
COHERENCE.