Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
wavenoise.cc
Go to the documentation of this file.
1 
2 #include "wavenoise.hh"
3 #include "TH2.h"
4 #include "TStyle.h"
5 #include "TCanvas.h"
6 
7 ClassImp(wavenoise) // used by THtml doc
8 
9 wavenoise& wavenoise::operator<<(wavenoise& a)
10 {
11  fChain= NULL;
12  fCurrent= a.fCurrent;
13  layer= a.layer;
14  ifo= a.ifo;
15  rms= a.rms;
16  frequency= a.frequency;
17  time= a.time;
18  gps= a.gps;
19  return *this;
20 }
21 
22 // Set branch addresses
23 void wavenoise::Init(TTree *tree)
24 {
25  if (tree == 0) return;
26  fChain = tree;
27  fCurrent = -1;
28  fChain->SetMakeClass(1);
29 
30  fChain->SetBranchAddress("layer",&layer);
31  fChain->SetBranchAddress("ifo",&ifo);
32  fChain->SetBranchAddress("rms",&rms);
33  fChain->SetBranchAddress("frequency",&frequency);
34  fChain->SetBranchAddress("time",&time);
35  fChain->SetBranchAddress("gps",&gps);
36 
37  Notify();
38 }
39 
40 
41 //++++++++++++++++++++++++++++++++++++++++++++++
42 // set noise wavenoise tree
43 //++++++++++++++++++++++++++++++++++++++++++++++
45 {
46  TTree* waveTree = new TTree("noise","noise");
47 
48  //==================================
49  // Define trigger tree
50  //==================================
51 
52  waveTree->Branch("layer", &layer, "layer/I");
53  waveTree->Branch("ifo", &ifo, "ifo/I");
54  waveTree->Branch("rms", &rms, "rms/D");
55  waveTree->Branch("frequency", &frequency, "frequency/F");
56  waveTree->Branch("time", &time, "time/D");
57  waveTree->Branch("gps", &gps, "gps/D");
58 
59  return waveTree;
60 }
61 
62 
63 //++++++++++++++++++++++++++++++++++++++++++++++++++++
64 // dump noise wavenoise into tree
65 //++++++++++++++++++++++++++++++++++++++++++++++++++++
66 void wavenoise::output(TTree* waveTree, WSeries<double>* p, int ifo_rms, double band_rms)
67 {
68  size_t i,j;
69  size_t n = p->maxLayer()+1; // number of layers
70  double rate_rms = p->wrate();
72 
73 //Fill tree
74 
75  this->gps = p->start();
76 
77  for(i=0; i<n; i++){
78  p->getLayer(a,i);
79  for(j=0; j<a.size(); j++){
80  this->layer= i;
81  this->ifo= ifo_rms;
82  this->rms= a.data[j];
83  this->time= this->gps+j/rate_rms;
84  this->frequency= (i+0.5)*band_rms/n;
85  waveTree->Fill();
86  }
87  }
88 }
89 
90 
91 
93 {
94  // Called when loading a new file.
95  // Get branch pointers.
96  b_layer = fChain->GetBranch("layer");
97  b_ifo = fChain->GetBranch("ifo");
98  b_rms = fChain->GetBranch("rms");
99  b_frequency = fChain->GetBranch("frequency");
100  b_time = fChain->GetBranch("time");
101  b_gps = fChain->GetBranch("gps");
102  return kTRUE;
103 }
104 
106 {
107  if (!fChain) return 0;
108  return fChain->GetEntry(entry);
109 };
110 
112 {
113 // Print contents of entry.
114 // If entry is not specified, print current entry
115  if (!fChain) return;
116  fChain->Show(entry);
117 }
118 
119 /*
120 Int_t wavenoise::LoadTree(Int_t entry)
121 {
122 // Set the environment to read one entry
123  if (!fChain) return -5;
124  Int_t centry = fChain->LoadTree(entry);
125  if (centry < 0) return centry;
126  if (fChain->IsA() != TChain::Class()) return centry;
127  TChain *chain = (TChain*)fChain;
128  if (chain->GetTreeNumber() != fCurrent) {
129  fCurrent = chain->GetTreeNumber();
130  Notify();
131  }
132  return centry;
133 }
134 
135 Int_t wavenoise::Cut(Int_t entry)
136 {
137 // This function may be called from Loop.
138 // returns 1 if entry is accepted.
139 // returns -1 otherwise.
140  return 1;
141 }
142 
143 void wavenoise::Loop()
144 {
145 // In a ROOT session, you can do:
146 // Root > .L wavenoise.C
147 // Root > wavenoise t
148 // Root > t.GetEntry(12); // Fill t data members with entry number 12
149 // Root > t.Show(); // Show rmss of entry 12
150 // Root > t.Show(16); // Read and show values of entry 16
151 // Root > t.Loop(); // Loop on all entries
152 //
153 
154 // This is the loop skeleton where:
155 // jentry is the global entry number in the chain
156 // ientry is the entry number in the current Tree
157 // Note that the argument to GetEntry must be:
158 // jentry for TChain::GetEntry
159 // ientry for TTree::GetEntry and TBranch::GetEntry
160 //
161 // To read only selected branches, Insert statements like:
162 // METHOD1:
163 // fChain->SetBranchStatus("*",0); // disable all branches
164 // fChain->SetBranchStatus("branchname",1); // activate branchname
165 // METHOD2: replace line
166 // fChain->GetEntry(jentry); //read all branches
167 //by b_branchname->GetEntry(ientry); //read only this branch
168  if (fChain == 0) return;
169 
170  Int_t nentries = Int_t(fChain->GetEntriesFast());
171 
172  Int_t nbytes = 0, nb = 0;
173  for (Int_t jentry=0; jentry<nentries;jentry++) {
174  Int_t ientry = LoadTree(jentry);
175  if (ientry < 0) break;
176  nb = fChain->GetEntry(jentry); nbytes += nb;
177  // if (Cut(ientry) < 0) continue;
178  }
179 }
180 */
181 
182 
183 
184 
185 
186 
187 
188 
TTree * tree
Definition: TimeSortTree.C:20
Double_t time
Definition: wavenoise.hh:30
virtual size_t size() const
Definition: wavearray.hh:127
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:10
double frequency
TBranch * b_layer
Definition: wavenoise.hh:35
TBranch * b_gps
Definition: wavenoise.hh:40
Float_t frequency
Definition: wavenoise.hh:29
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
i drho i
Int_t ifo
Definition: wavenoise.hh:27
char ifo[NIFO_MAX][8]
void wrate(double r)
Definition: wseries.hh:102
TTree * setTree()
Definition: wavenoise.cc:44
double time[6]
Definition: cbc_plots.C:435
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: wavenoise.hh:21
TBranch * b_ifo
Definition: wavenoise.hh:36
Int_t GetEntry(Int_t)
Definition: wavenoise.cc:105
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
TBranch * b_time
Definition: wavenoise.hh:39
TTree * fChain
Definition: wavenoise.hh:20
void Init(TTree *)
Definition: wavenoise.cc:23
TBranch * b_rms
Definition: wavenoise.hh:37
void Show(Int_t entry=-1)
Definition: wavenoise.cc:111
double * entry
Definition: cwb_setcuts.C:206
Int_t layer
current Tree number in a TChain
Definition: wavenoise.hh:26
void output(TTree *, WSeries< double > *, int=0, double=1.)
Definition: wavenoise.cc:66
double gps
TBranch * b_frequency
Definition: wavenoise.hh:38
DataType_t * data
Definition: wavearray.hh:301
Double_t gps
Definition: wavenoise.hh:31
Double_t rms
Definition: wavenoise.hh:28
Bool_t Notify()
Definition: wavenoise.cc:92
int maxLayer()
Definition: wseries.hh:121