Logo coherent WaveBurst  
Library Reference Guide
Logo
wavenoise.cc
Go to the documentation of this file.
1 /*
2 # Copyright (C) 2019 Gabriele Vedovato, Sergey Klimenko, Sergey Klimenko
3 #
4 # This program is free software: you can redistribute it and/or modify
5 # it under the terms of the GNU General Public License as published by
6 # the Free Software Foundation, either version 3 of the License, or
7 # (at your option) any later version.
8 #
9 # This program is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 # GNU General Public License for more details.
13 #
14 # You should have received a copy of the GNU General Public License
15 # along with this program. If not, see <https://www.gnu.org/licenses/>.
16 */
17 
18 
19 
20 #include "wavenoise.hh"
21 #include "TH2.h"
22 #include "TStyle.h"
23 #include "TCanvas.h"
24 
25 ClassImp(wavenoise) // used by THtml doc
26 
28 {
29  fChain= NULL;
30  fCurrent= a.fCurrent;
31  layer= a.layer;
32  ifo= a.ifo;
33  rms= a.rms;
34  frequency= a.frequency;
35  time= a.time;
36  gps= a.gps;
37  return *this;
38 }
39 
40 // Set branch addresses
41 void wavenoise::Init(TTree *tree)
42 {
43  if (tree == 0) return;
44  fChain = tree;
45  fCurrent = -1;
46  fChain->SetMakeClass(1);
47 
48  fChain->SetBranchAddress("layer",&layer);
49  fChain->SetBranchAddress("ifo",&ifo);
50  fChain->SetBranchAddress("rms",&rms);
51  fChain->SetBranchAddress("frequency",&frequency);
52  fChain->SetBranchAddress("time",&time);
53  fChain->SetBranchAddress("gps",&gps);
54 
55  Notify();
56 }
57 
58 
59 //++++++++++++++++++++++++++++++++++++++++++++++
60 // set noise wavenoise tree
61 //++++++++++++++++++++++++++++++++++++++++++++++
63 {
64  TTree* waveTree = new TTree("noise","noise");
65 
66  //==================================
67  // Define trigger tree
68  //==================================
69 
70  waveTree->Branch("layer", &layer, "layer/I");
71  waveTree->Branch("ifo", &ifo, "ifo/I");
72  waveTree->Branch("rms", &rms, "rms/D");
73  waveTree->Branch("frequency", &frequency, "frequency/F");
74  waveTree->Branch("time", &time, "time/D");
75  waveTree->Branch("gps", &gps, "gps/D");
76 
77  return waveTree;
78 }
79 
80 
81 //++++++++++++++++++++++++++++++++++++++++++++++++++++
82 // dump noise wavenoise into tree
83 //++++++++++++++++++++++++++++++++++++++++++++++++++++
84 void wavenoise::output(TTree* waveTree, WSeries<double>* p, int ifo_rms, double band_rms)
85 {
86  size_t i,j;
87  size_t n = p->maxLayer()+1; // number of layers
88  double rate_rms = p->wrate();
90 
91 //Fill tree
92 
93  this->gps = p->start();
94 
95  for(i=0; i<n; i++){
96  p->getLayer(a,i);
97  for(j=0; j<a.size(); j++){
98  this->layer= i;
99  this->ifo= ifo_rms;
100  this->rms= a.data[j];
101  this->time= this->gps+j/rate_rms;
102  this->frequency= (i+0.5)*band_rms/n;
103  waveTree->Fill();
104  }
105  }
106 }
107 
108 
109 
111 {
112  // Called when loading a new file.
113  // Get branch pointers.
114  b_layer = fChain->GetBranch("layer");
115  b_ifo = fChain->GetBranch("ifo");
116  b_rms = fChain->GetBranch("rms");
117  b_frequency = fChain->GetBranch("frequency");
118  b_time = fChain->GetBranch("time");
119  b_gps = fChain->GetBranch("gps");
120  return kTRUE;
121 }
122 
124 {
125  if (!fChain) return 0;
126  return fChain->GetEntry(entry);
127 };
128 
130 {
131 // Print contents of entry.
132 // If entry is not specified, print current entry
133  if (!fChain) return;
134  fChain->Show(entry);
135 }
136 
137 /*
138 Int_t wavenoise::LoadTree(Int_t entry)
139 {
140 // Set the environment to read one entry
141  if (!fChain) return -5;
142  Int_t centry = fChain->LoadTree(entry);
143  if (centry < 0) return centry;
144  if (fChain->IsA() != TChain::Class()) return centry;
145  TChain *chain = (TChain*)fChain;
146  if (chain->GetTreeNumber() != fCurrent) {
147  fCurrent = chain->GetTreeNumber();
148  Notify();
149  }
150  return centry;
151 }
152 
153 Int_t wavenoise::Cut(Int_t entry)
154 {
155 // This function may be called from Loop.
156 // returns 1 if entry is accepted.
157 // returns -1 otherwise.
158  return 1;
159 }
160 
161 void wavenoise::Loop()
162 {
163 // In a ROOT session, you can do:
164 // Root > .L wavenoise.C
165 // Root > wavenoise t
166 // Root > t.GetEntry(12); // Fill t data members with entry number 12
167 // Root > t.Show(); // Show rmss of entry 12
168 // Root > t.Show(16); // Read and show values of entry 16
169 // Root > t.Loop(); // Loop on all entries
170 //
171 
172 // This is the loop skeleton where:
173 // jentry is the global entry number in the chain
174 // ientry is the entry number in the current Tree
175 // Note that the argument to GetEntry must be:
176 // jentry for TChain::GetEntry
177 // ientry for TTree::GetEntry and TBranch::GetEntry
178 //
179 // To read only selected branches, Insert statements like:
180 // METHOD1:
181 // fChain->SetBranchStatus("*",0); // disable all branches
182 // fChain->SetBranchStatus("branchname",1); // activate branchname
183 // METHOD2: replace line
184 // fChain->GetEntry(jentry); //read all branches
185 //by b_branchname->GetEntry(ientry); //read only this branch
186  if (fChain == 0) return;
187 
188  Int_t nentries = Int_t(fChain->GetEntriesFast());
189 
190  Int_t nbytes = 0, nb = 0;
191  for (Int_t jentry=0; jentry<nentries;jentry++) {
192  Int_t ientry = LoadTree(jentry);
193  if (ientry < 0) break;
194  nb = fChain->GetEntry(jentry); nbytes += nb;
195  // if (Cut(ientry) < 0) continue;
196  }
197 }
198 */
199 
200 
201 
202 
203 
204 
205 
206 
TTree * tree
Definition: TimeSortTree.C:20
Double_t time
Definition: wavenoise.hh:48
wavearray< double > a(hp.size())
int n
Definition: cwb_net.C:28
double frequency
TBranch * b_layer
Definition: wavenoise.hh:53
TBranch * b_gps
Definition: wavenoise.hh:58
Float_t frequency
Definition: wavenoise.hh:47
virtual void start(double s)
Definition: wavearray.hh:137
int j
Definition: cwb_net.C:28
i drho i
Int_t ifo
Definition: wavenoise.hh:45
char ifo[NIFO_MAX][8]
void wrate(double r)
Definition: wseries.hh:120
virtual size_t size() const
Definition: wavearray.hh:145
TTree * setTree()
Definition: wavenoise.cc:62
Int_t fCurrent
pointer to the analyzed TTree or TChain
Definition: wavenoise.hh:39
TBranch * b_ifo
Definition: wavenoise.hh:54
Int_t GetEntry(Int_t)
Definition: wavenoise.cc:123
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:193
TBranch * b_time
Definition: wavenoise.hh:57
TTree * fChain
Definition: wavenoise.hh:38
void Init(TTree *)
Definition: wavenoise.cc:41
TBranch * b_rms
Definition: wavenoise.hh:55
void Show(Int_t entry=-1)
Definition: wavenoise.cc:129
double * entry
Definition: cwb_setcuts.C:224
Int_t layer
current Tree number in a TChain
Definition: wavenoise.hh:44
void output(TTree *, WSeries< double > *, int=0, double=1.)
Definition: wavenoise.cc:84
double gps
TBranch * b_frequency
Definition: wavenoise.hh:56
DataType_t * data
Definition: wavearray.hh:319
Double_t gps
Definition: wavenoise.hh:49
Double_t rms
Definition: wavenoise.hh:46
Bool_t Notify()
Definition: wavenoise.cc:110
int maxLayer()
Definition: wseries.hh:139