CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
vDriftHistos.h
Go to the documentation of this file.
1 #ifndef vDriftHistos_H
2 #define vDriftHistos_H
3 
4 #include "TH1.h"
5 #include "TString.h"
6 #include "TFile.h"
7 #include "DTTMax.h"
8 #include <string>
9 
10 // A set of histograms on chamber angle and position
11 class h4DSegm{
12  public:
13  h4DSegm(std::string name_){
14  TString N = name_.c_str();
15  name=name_.c_str();
16  h4DSegmXPosInCham = new TH1F(N+"_h4DSegmXPosInCham",
17  "4D Segment x position (cm) in Chamber RF", 200, -200, 200);
18  h4DSegmYPosInCham = new TH1F(N+"_h4DSegmYPosInCham",
19  "4D Segment y position (cm) in Chamber RF", 200, -200, 200);
20  h4DSegmPhiAngleInCham = new TH1F(N+"_h4DSegmPhiAngleInCham",
21  "4D Segment phi angle (rad) in Chamber RF", 180, -180, 180);
22  h4DSegmThetaAngleInCham = new TH1F(N+"_h4DSegmThetaAngleInCham",
23  "4D Segment theta angle (rad) in Chamber RF", 180, -180, 180);
24  h4DSegmImpactAngleInCham = new TH1F(N+"_h4DSegmImpactAngleInCham",
25  "4D Segment impact angle (rad) in Chamber RF", 180, -180, 180);
26  }
27  h4DSegm(TString name_, TFile* file){
28  name=name_;
29  h4DSegmXPosInCham = (TH1F *) file->Get(name+"_h4DSegmXPosInCham");
30  h4DSegmYPosInCham = (TH1F *) file->Get(name+"_h4DSegmYPosInCham");
31  h4DSegmPhiAngleInCham = (TH1F *) file->Get(name+"_h4DSegmPhiAngleInCham");
32  h4DSegmThetaAngleInCham = (TH1F *) file->Get(name+"_h4DSegmThetaAngleInCham");
33  h4DSegmImpactAngleInCham = (TH1F *) file->Get(name+"_h4DSegmImpactAngleInCham");
34  }
36  delete h4DSegmXPosInCham;
37  delete h4DSegmYPosInCham;
38  delete h4DSegmPhiAngleInCham;
41  }
42 void Fill(float x, float y, float phi, float theta, float impact) {
43  h4DSegmXPosInCham->Fill(x);
44  h4DSegmYPosInCham->Fill(y);
45  h4DSegmPhiAngleInCham->Fill(phi);
46  h4DSegmThetaAngleInCham->Fill(theta);
47  h4DSegmImpactAngleInCham->Fill(impact);
48 }
49 void Fill(float x, float phi) {
50  h4DSegmXPosInCham->Fill(x);
51  h4DSegmPhiAngleInCham->Fill(phi);
52 }
53  void Write() {
54  h4DSegmXPosInCham->Write();
55  h4DSegmYPosInCham->Write();
56  h4DSegmPhiAngleInCham->Write();
57  h4DSegmThetaAngleInCham->Write();
58  h4DSegmImpactAngleInCham->Write();
59  }
60  public:
61 
67 
68  TString name;
69 };
70 
71 
72 // A set of histograms on SL angle and position
73 class h2DSegm{
74  public:
75  h2DSegm(std::string name_){
76  TString N = name_.c_str();
77  name=name_.c_str();
78  h2DSegmPosInCham = new TH1F(N+"_h2DSegmPosInCham",
79  "2D Segment position (cm) in Chamber RF", 200, -200, 200);
80  h2DSegmAngleInCham = new TH1F(N+"_h2DSegmAngleInCham",
81  "2D Segment angle (rad) in Chamber RF", 200, -2, 2);
82  h2DSegmCosAngleInCham = new TH1F(N+"_h2DSegmCosAngleInCham",
83  "2D Segment cos(angle) in Chamber RF", 200, -2, 2);
84  }
85  h2DSegm(TString name_, TFile* file){
86  name=name_;
87 
88  h2DSegmPosInCham = (TH1F *) file->Get(name+"_h2DSegmPosInCham");
89  h2DSegmAngleInCham = (TH1F *) file->Get(name+"_h2DSegmAngleInCham");
90  h2DSegmCosAngleInCham = (TH1F *) file->Get(name+"_h2DSegmCosAngleInCham");
91  }
93  delete h2DSegmPosInCham;
94  delete h2DSegmAngleInCham;
95  delete h2DSegmCosAngleInCham;
96  }
97  void Fill(float pos, float localAngle) {
98 
99  h2DSegmPosInCham->Fill(pos);
100  h2DSegmAngleInCham->Fill(atan(localAngle));
101  h2DSegmCosAngleInCham->Fill(cos(atan(localAngle)));
102  }
103  void Write() {
104  h2DSegmPosInCham->Write();
105  h2DSegmAngleInCham->Write();
106  h2DSegmCosAngleInCham->Write();
107  }
108  public:
109 
113 
114  TString name;
115 };
116 
117 // A set of histograms on SL Tmax
118 class hTMaxCell{
119  public:
120  hTMaxCell(TString name_){
121  name = name_;
122 
123  // book TMax histograms
124  hTmax123 = new TH1F (name+"_Tmax123", "Tmax123 value", 2000, -1000., 1000.);
125  hTmax124s72 = new TH1F (name+"_Tmax124_s72", "Tmax124 sigma=sqrt(7/2) value", 2000, -1000., 1000.);
126  hTmax124s78 = new TH1F (name+"_Tmax124_s78", "Tmax124 sigma=sqrt(7/8) value", 2000, -1000., 1000.);
127  hTmax134s72 = new TH1F (name+"_Tmax134_s72", "Tmax134 sigma=sqrt(7/2) value", 2000, -1000., 1000.);
128  hTmax134s78 = new TH1F (name+"_Tmax134_s78", "Tmax134 sigma=sqrt(7/8) value", 2000, -1000., 1000.);
129  hTmax234 = new TH1F (name+"_Tmax234", "Tmax234 value", 2000, -1000., 1000.);
130  hTmax_3t0 = new TH1F (name+"_3t0", "Tmax+3*Delta(t0)", 2000, -1000., 1000.);
131  hTmax_3t0_0 = new TH1F (name+"_3t0_0", "Tmax+3*Delta(t0); 3 hits", 2000, -1000., 1000.);
132  hTmax_3t0_1 = new TH1F (name+"_3t0_1", "Tmax+3*Delta(t0); one t<5ns", 2000, -1000., 1000.);
133  hTmax_3t0_2 = new TH1F (name+"_3t0_2", "Tmax+3*Delta(t0); one t<10ns", 2000, -1000., 1000.);
134  hTmax_3t0_3 = new TH1F (name+"_3t0_3", "Tmax+3*Delta(t0); one t<20ns", 2000, -1000., 1000.);
135  hTmax_3t0_4 = new TH1F (name+"_3t0_4", "Tmax+3*Delta(t0); one t<50ns", 2000, -1000., 1000.);
136  hTmax_3t0_5 = new TH1F (name+"_3t0_5", "Tmax+3*Delta(t0); all t>50ns", 2000, -1000., 1000.);
137  hTmax_2t0 = new TH1F (name+"_2t0", "Tmax+2*Delta(t0)", 2000, -1000., 1000.);
138  hTmax_2t0_0 = new TH1F (name+"_2t0_0", "Tmax+2*Delta(t0); 3 hits", 2000, -1000., 1000.);
139  hTmax_2t0_1 = new TH1F (name+"_2t0_1", "Tmax+2*Delta(t0); one t<5ns", 2000, -1000., 1000.);
140  hTmax_2t0_2 = new TH1F (name+"_2t0_2", "Tmax+2*Delta(t0); one t<10ns", 2000, -1000., 1000.);
141  hTmax_2t0_3 = new TH1F (name+"_2t0_3", "Tmax+2*Delta(t0); one t<20ns", 2000, -1000., 1000.);
142  hTmax_2t0_4 = new TH1F (name+"_2t0_4", "Tmax+2*Delta(t0); one t<50ns", 2000, -1000., 1000.);
143  hTmax_2t0_5 = new TH1F (name+"_2t0_5", "Tmax+2*Delta(t0); all t>50ns", 2000, -1000., 1000.);
144  hTmax_t0 = new TH1F (name+"_t0", "Tmax+Delta(t0)", 2000, -1000., 1000.);
145  hTmax_t0_0 = new TH1F (name+"_t0_0", "Tmax+Delta(t0); 3 hits", 2000, -1000., 1000.);
146  hTmax_t0_1 = new TH1F (name+"_t0_1", "Tmax+Delta(t0); one t<5ns", 2000, -1000., 1000.);
147  hTmax_t0_2 = new TH1F (name+"_t0_2", "Tmax+Delta(t0); one t<10ns", 2000, -1000., 1000.);
148  hTmax_t0_3 = new TH1F (name+"_t0_3", "Tmax+Delta(t0); one t<20ns", 2000, -1000., 1000.);
149  hTmax_t0_4 = new TH1F (name+"_t0_4", "Tmax+Delta(t0); one t<50ns", 2000, -1000., 1000.);
150  hTmax_t0_5 = new TH1F (name+"_t0_5", "Tmax+Delta(t0); all t>50ns", 2000, -1000., 1000.);
151  hTmax_0 = new TH1F (name+"_0", "Tmax", 2000, -1000., 1000.);
152  }
153 
154 
155  hTMaxCell (TString name_, TFile* file){
156  name=name_;
157  hTmax123 = (TH1F *) file->Get(name+"_Tmax123");
158  hTmax124s72 = (TH1F *) file->Get(name+"_Tmax124_s72");
159  hTmax124s78 = (TH1F *) file->Get(name+"_Tmax124_s78");
160  hTmax134s72 = (TH1F *) file->Get(name+"_Tmax134_s72");
161  hTmax134s78 = (TH1F *) file->Get(name+"_Tmax134_s78");
162  hTmax234 = (TH1F *) file->Get(name+"_Tmax234");
163  hTmax_3t0 = (TH1F *) file->Get(name+"_3t0");
164  hTmax_3t0_0 = (TH1F *) file->Get(name+"_3t0_0");
165  hTmax_3t0_1 = (TH1F *) file->Get(name+"_3t0_1");
166  hTmax_3t0_2 = (TH1F *) file->Get(name+"_3t0_2");
167  hTmax_3t0_3 = (TH1F *) file->Get(name+"_3t0_3");
168  hTmax_3t0_4 = (TH1F *) file->Get(name+"_3t0_4");
169  hTmax_3t0_5 = (TH1F *) file->Get(name+"_3t0_5");
170  hTmax_2t0 = (TH1F *) file->Get(name+"_2t0");
171  hTmax_2t0_0 = (TH1F *) file->Get(name+"_2t0_0");
172  hTmax_2t0_1 = (TH1F *) file->Get(name+"_2t0_1");
173  hTmax_2t0_2 = (TH1F *) file->Get(name+"_2t0_2");
174  hTmax_2t0_3 = (TH1F *) file->Get(name+"_2t0_3");
175  hTmax_2t0_4 = (TH1F *) file->Get(name+"_2t0_4");
176  hTmax_2t0_5 = (TH1F *) file->Get(name+"_2t0_5");
177  hTmax_t0 = (TH1F *) file->Get(name+"_t0");
178  hTmax_t0_1 = (TH1F *) file->Get(name+"_t0_1");
179  hTmax_t0_2 = (TH1F *) file->Get(name+"_t0_2");
180  hTmax_t0_3 = (TH1F *) file->Get(name+"_t0_3");
181  hTmax_t0_4 = (TH1F *) file->Get(name+"_t0_4");
182  hTmax_t0_5 = (TH1F *) file->Get(name+"_t0_5");
183  hTmax_0 = (TH1F *) file->Get(name+"_0");
184 
185  }
186 
187 
189  delete hTmax123;
190  delete hTmax124s72;
191  delete hTmax124s78;
192  delete hTmax134s72;
193  delete hTmax134s78;
194  delete hTmax234;
195  delete hTmax_3t0;
196  delete hTmax_3t0_0;
197  delete hTmax_3t0_1;
198  delete hTmax_3t0_2;
199  delete hTmax_3t0_3;
200  delete hTmax_3t0_4;
201  delete hTmax_3t0_5;
202  delete hTmax_2t0;
203  delete hTmax_2t0_0;
204  delete hTmax_2t0_1;
205  delete hTmax_2t0_2;
206  delete hTmax_2t0_3;
207  delete hTmax_2t0_4;
208  delete hTmax_2t0_5;
209  delete hTmax_t0;
210  delete hTmax_t0_0;
211  delete hTmax_t0_1;
212  delete hTmax_t0_2;
213  delete hTmax_t0_3;
214  delete hTmax_t0_4;
215  delete hTmax_t0_5;
216  delete hTmax_0;
217 
218  }
219 
220  void Fill(float tmax123, float tmax124, float tmax134, float tmax234,
221  dttmaxenums::SigmaFactor s124,// Give the factor relating the width of the TMax distribution
222  dttmaxenums::SigmaFactor s134,// and the cell resolution (for tmax123 and tmax234 is always= sqrt(3/2)).
223  unsigned t0_123, // Give the "quantity" of Delta(t0) included in the tmax
224  unsigned t0_124, // formula used for Tmax123 or Tma124 or Tma134 or Tma234
225  unsigned t0_134,
226  unsigned t0_234,
227  unsigned hSubGroup//different t0 hists(at least one hit within a given distance from the wire)
228  ) {
229 
230  if(tmax123 > 0.) {
231  hTmax123->Fill(tmax123);
232  if(t0_123==1) {
233  hTmax_t0->Fill(tmax123);
234  switch(hSubGroup) {
235  case 0: hTmax_t0_0->Fill(tmax123); break;
236  case 1: hTmax_t0_1->Fill(tmax123); break;
237  case 2: hTmax_t0_2->Fill(tmax123); break;
238  case 3: hTmax_t0_3->Fill(tmax123); break;
239  case 4: hTmax_t0_4->Fill(tmax123); break;
240  case 99: hTmax_t0_5->Fill(tmax123); break;
241  }
242  }
243  else {
244  hTmax_2t0->Fill(tmax123);
245  switch(hSubGroup) {
246  case 0: hTmax_2t0_0->Fill(tmax123); break;
247  case 1: hTmax_2t0_1->Fill(tmax123); break;
248  case 2: hTmax_2t0_2->Fill(tmax123); break;
249  case 3: hTmax_2t0_3->Fill(tmax123); break;
250  case 4: hTmax_2t0_4->Fill(tmax123); break;
251  case 99: hTmax_2t0_5->Fill(tmax123); break;
252  }
253  }
254  }
255  if(tmax124 > 0.) {
256  (s124==dttmaxenums::r72)? hTmax124s72->Fill(tmax124):hTmax124s78->Fill(tmax124);
257  if(t0_124==0)
258  hTmax_0->Fill(tmax124);
259  else if(t0_124==1) {
260  hTmax_t0->Fill(tmax124);
261  switch(hSubGroup) {
262  case 0: hTmax_t0_0->Fill(tmax124); break;
263  case 1: hTmax_t0_1->Fill(tmax124); break;
264  case 2: hTmax_t0_2->Fill(tmax124); break;
265  case 3: hTmax_t0_3->Fill(tmax124); break;
266  case 4: hTmax_t0_4->Fill(tmax124); break;
267  case 99: hTmax_t0_5->Fill(tmax124); break;
268  }
269  }
270  else if(t0_124== 2) {
271  hTmax_2t0->Fill(tmax124);
272  switch(hSubGroup) {
273  case 0: hTmax_2t0_0->Fill(tmax124); break;
274  case 1: hTmax_2t0_1->Fill(tmax124); break;
275  case 2: hTmax_2t0_2->Fill(tmax124); break;
276  case 3: hTmax_2t0_3->Fill(tmax124); break;
277  case 4: hTmax_2t0_4->Fill(tmax124); break;
278  case 99: hTmax_2t0_5->Fill(tmax124); break;
279  }
280  }
281  else if(t0_124==3) {
282  hTmax_3t0->Fill(tmax124);
283  switch(hSubGroup) {
284  case 0: hTmax_3t0_0->Fill(tmax124); break;
285  case 1: hTmax_3t0_1->Fill(tmax124); break;
286  case 2: hTmax_3t0_2->Fill(tmax124); break;
287  case 3: hTmax_3t0_3->Fill(tmax124); break;
288  case 4: hTmax_3t0_4->Fill(tmax124); break;
289  case 99: hTmax_3t0_5->Fill(tmax124);break;
290  }
291  }
292  }
293  if(tmax134 > 0.) {
294  (s134==dttmaxenums::r72)? hTmax134s72->Fill(tmax134):hTmax134s78->Fill(tmax134);
295  if(t0_134==0)
296  hTmax_0->Fill(tmax134);
297  else if(t0_134==1) {
298  hTmax_t0->Fill(tmax134);
299  switch(hSubGroup) {
300  case 0: hTmax_t0_0->Fill(tmax134); break;
301  case 1: hTmax_t0_1->Fill(tmax134); break;
302  case 2: hTmax_t0_2->Fill(tmax134); break;
303  case 3: hTmax_t0_3->Fill(tmax134); break;
304  case 4: hTmax_t0_4->Fill(tmax134); break;
305  case 99: hTmax_t0_5->Fill(tmax134); break;
306  }
307  }
308  else if(t0_134== 2) {
309  hTmax_2t0->Fill(tmax134);
310  switch(hSubGroup) {
311  case 0: hTmax_2t0_0->Fill(tmax134); break;
312  case 1: hTmax_2t0_1->Fill(tmax134); break;
313  case 2: hTmax_2t0_2->Fill(tmax134); break;
314  case 3: hTmax_2t0_3->Fill(tmax134); break;
315  case 4: hTmax_2t0_4->Fill(tmax134); break;
316  case 99: hTmax_2t0_5->Fill(tmax134); break;
317  }
318  }
319  else if(t0_134==3) {
320  hTmax_3t0->Fill(tmax134);
321  switch(hSubGroup) {
322  case 0: hTmax_3t0_0->Fill(tmax134); break;
323  case 1: hTmax_3t0_1->Fill(tmax134); break;
324  case 2: hTmax_3t0_2->Fill(tmax134); break;
325  case 3: hTmax_3t0_3->Fill(tmax134); break;
326  case 4: hTmax_3t0_4->Fill(tmax134); break;
327  case 99: hTmax_3t0_5->Fill(tmax134); break;
328  }
329  }
330  }
331  if(tmax234 > 0.) {
332  hTmax234->Fill(tmax234);
333  if(t0_234==1) {
334  hTmax_t0->Fill(tmax234);
335  switch(hSubGroup) {
336  case 0: hTmax_t0_0->Fill(tmax234); break;
337  case 1: hTmax_t0_1->Fill(tmax234); break;
338  case 2: hTmax_t0_2->Fill(tmax234); break;
339  case 3: hTmax_t0_3->Fill(tmax234); break;
340  case 4: hTmax_t0_4->Fill(tmax234); break;
341  case 99: hTmax_t0_5->Fill(tmax234);break;
342  }
343  }
344  else {
345  hTmax_2t0->Fill(tmax234);
346  switch(hSubGroup) {
347  case 0: hTmax_2t0_0->Fill(tmax234); break;
348  case 1: hTmax_2t0_1->Fill(tmax234); break;
349  case 2: hTmax_2t0_2->Fill(tmax234); break;
350  case 3: hTmax_2t0_3->Fill(tmax234); break;
351  case 4: hTmax_2t0_4->Fill(tmax234); break;
352  case 99: hTmax_2t0_5->Fill(tmax234); break;
353  }
354  }
355  }
356  }
357 
358  void Write() {
359  // write the Tmax histograms
360  hTmax123->Write();
361  hTmax124s72->Write();
362  hTmax124s78->Write();
363  hTmax134s72->Write();
364  hTmax134s78->Write();
365  hTmax234->Write();
366  hTmax_3t0->Write();
367  hTmax_3t0_0->Write();
368  hTmax_3t0_1->Write();
369  hTmax_3t0_2->Write();
370  hTmax_3t0_3->Write();
371  hTmax_3t0_4->Write();
372  hTmax_3t0_5->Write();
373  hTmax_2t0->Write();
374  hTmax_2t0_0->Write();
375  hTmax_2t0_1->Write();
376  hTmax_2t0_2->Write();
377  hTmax_2t0_3->Write();
378  hTmax_2t0_4->Write();
379  hTmax_2t0_5->Write();
380  hTmax_t0->Write();
381  hTmax_t0_0->Write();
382  hTmax_t0_1->Write();
383  hTmax_t0_2->Write();
384  hTmax_t0_3->Write();
385  hTmax_t0_4->Write();
386  hTmax_t0_5->Write();
387  hTmax_0->Write();
388 
389  }
390 
391  int GetT0Factor(TH1F* hist) {
392  unsigned t0 = 999;
393 
394  if (hist == hTmax_3t0)
395  t0 = 3;
396  else if(hist == hTmax_2t0)
397  t0 = 2;
398  else if(hist == hTmax_t0)
399  t0 = 1;
400  else if (hist == hTmax_0)
401  t0 = 0;
402 
403  return t0;
404  }
405 
406  TH1F* hTmax123;
407  TH1F* hTmax124s72;
408  TH1F* hTmax124s78;
409  TH1F* hTmax134s72;
410  TH1F* hTmax134s78;
411  TH1F* hTmax234;
412  TH1F* hTmax_3t0;
413  TH1F* hTmax_3t0_0;
414  TH1F* hTmax_3t0_1;
415  TH1F* hTmax_3t0_2;
416  TH1F* hTmax_3t0_3;
417  TH1F* hTmax_3t0_4;
418  TH1F* hTmax_3t0_5;
419  TH1F* hTmax_2t0;
420  TH1F* hTmax_2t0_0;
421  TH1F* hTmax_2t0_1;
422  TH1F* hTmax_2t0_2;
423  TH1F* hTmax_2t0_3;
424  TH1F* hTmax_2t0_4;
425  TH1F* hTmax_2t0_5;
426  TH1F* hTmax_t0;
427  TH1F* hTmax_t0_0;
428  TH1F* hTmax_t0_1;
429  TH1F* hTmax_t0_2;
430  TH1F* hTmax_t0_3;
431  TH1F* hTmax_t0_4;
432  TH1F* hTmax_t0_5;
433  TH1F* hTmax_0;
434 
435  TString name;
436 };
437 
438 #endif
TH1F * h4DSegmThetaAngleInCham
Definition: vDriftHistos.h:65
TH1F * hTmax134s78
Definition: vDriftHistos.h:410
TH1F * hTmax_2t0
Definition: vDriftHistos.h:419
TH1F * hTmax_t0_3
Definition: vDriftHistos.h:430
TH1F * h2DSegmPosInCham
Definition: vDriftHistos.h:110
TH1F * hTmax_3t0
Definition: vDriftHistos.h:412
hTMaxCell(TString name_, TFile *file)
Definition: vDriftHistos.h:155
int GetT0Factor(TH1F *hist)
Definition: vDriftHistos.h:391
void Fill(float x, float phi)
Definition: vDriftHistos.h:49
Geom::Theta< T > theta() const
TString name
Definition: vDriftHistos.h:114
void Fill(float x, float y, float phi, float theta, float impact)
Definition: vDriftHistos.h:42
TH1F * hTmax234
Definition: vDriftHistos.h:411
void Fill(float pos, float localAngle)
Definition: vDriftHistos.h:97
TH1F * hTmax_2t0_2
Definition: vDriftHistos.h:422
void Write()
Definition: vDriftHistos.h:103
TH1F * hTmax_3t0_4
Definition: vDriftHistos.h:417
TH1F * hTmax_t0
Definition: vDriftHistos.h:426
TH1F * h2DSegmAngleInCham
Definition: vDriftHistos.h:111
h4DSegm(TString name_, TFile *file)
Definition: vDriftHistos.h:27
TH1F * hTmax134s72
Definition: vDriftHistos.h:409
TH1F * hTmax_2t0_0
Definition: vDriftHistos.h:420
TH1F * hTmax124s72
Definition: vDriftHistos.h:407
TH1F * hTmax_2t0_4
Definition: vDriftHistos.h:424
TH1F * hTmax_0
Definition: vDriftHistos.h:433
TH1F * hTmax_2t0_3
Definition: vDriftHistos.h:423
TH1F * hTmax_t0_0
Definition: vDriftHistos.h:427
TH1F * h4DSegmImpactAngleInCham
Definition: vDriftHistos.h:66
TH1F * hTmax_t0_1
Definition: vDriftHistos.h:428
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
TH1F * h4DSegmYPosInCham
Definition: vDriftHistos.h:63
TH1F * hTmax124s78
Definition: vDriftHistos.h:408
TH1F * hTmax_2t0_1
Definition: vDriftHistos.h:421
void Write()
Definition: vDriftHistos.h:358
TString name
Definition: vDriftHistos.h:435
TH1F * hTmax_t0_5
Definition: vDriftHistos.h:432
TH1F * h4DSegmXPosInCham
Definition: vDriftHistos.h:62
void Fill(float tmax123, float tmax124, float tmax134, float tmax234, dttmaxenums::SigmaFactor s124, dttmaxenums::SigmaFactor s134, unsigned t0_123, unsigned t0_124, unsigned t0_134, unsigned t0_234, unsigned hSubGroup)
Definition: vDriftHistos.h:220
TH1F * hTmax_3t0_2
Definition: vDriftHistos.h:415
TH1F * hTmax_t0_2
Definition: vDriftHistos.h:429
#define N
Definition: blowfish.cc:9
TH1F * h2DSegmCosAngleInCham
Definition: vDriftHistos.h:112
TH1F * hTmax_2t0_5
Definition: vDriftHistos.h:425
h2DSegm(std::string name_)
Definition: vDriftHistos.h:75
TH1F * hTmax_3t0_3
Definition: vDriftHistos.h:416
TH1F * hTmax_3t0_1
Definition: vDriftHistos.h:414
TString name
Definition: vDriftHistos.h:68
TH1F * hTmax_3t0_0
Definition: vDriftHistos.h:413
h2DSegm(TString name_, TFile *file)
Definition: vDriftHistos.h:85
TH1F * h4DSegmPhiAngleInCham
Definition: vDriftHistos.h:64
TH1F * hTmax123
Definition: vDriftHistos.h:406
x
Definition: VDTMath.h:216
TH1F * hTmax_3t0_5
Definition: vDriftHistos.h:418
hTMaxCell(TString name_)
Definition: vDriftHistos.h:120
h4DSegm(std::string name_)
Definition: vDriftHistos.h:13
void Write()
Definition: vDriftHistos.h:53
Definition: DDAxes.h:10
TH1F * hTmax_t0_4
Definition: vDriftHistos.h:431