CbmRoot
Simple/LxTrackAnaSegments.cxx
Go to the documentation of this file.
1 #include "LxTrackAnaSegments.h"
2 #include "LxTrackAna.h"
3 #include "TH1.h"
4 #include "TH2.h"
5 #include <cmath>
6 #include <dirent.h>
7 #include <sys/stat.h>
8 #include <sys/types.h>
9 
10 using namespace std;
11 
16 
19 
22 
27 
30 
33 
36 
39 
42 
43 static TH2F* muchXTxCovHisto[LXSTATIONS - 1];
44 static TH2F* muchYTyCovHisto[LXSTATIONS - 1];
45 
50 
53 
54 static TString particleType("jpsi");
55 
57  particleType = v;
58 
59  if (v == "omega") stationsInAlgo = 5;
60 }
61 
63  char name[64];
64  char title[256];
65 
66  for (Int_t i = 0; i < LXSTATIONS; ++i) {
67  sprintf(name, "muchInStationXDispLeft_%d", i);
68  sprintf(
69  title, "X dispersion from central to left layer inside station: %d", i);
70  muchInStationXDispLeft[i] = new TH1F(name, title, 100, -3.0, 3.0);
71  muchInStationXDispLeft[i]->StatOverflows();
72 
73  sprintf(name, "muchInStationXDispRight_%d", i);
74  sprintf(
75  title, "X dispersion from central to right layer inside station: %d", i);
76  muchInStationXDispRight[i] = new TH1F(name, title, 100, -3.0, 3.0);
77  muchInStationXDispRight[i]->StatOverflows();
78 
79  sprintf(name, "muchInStationYDispLeft_%d", i);
80  sprintf(
81  title, "Y dispersion from central to left layer inside station: %d", i);
82  muchInStationYDispLeft[i] = new TH1F(name, title, 100, -3.0, 3.0);
83  muchInStationYDispLeft[i]->StatOverflows();
84 
85  sprintf(name, "muchInStationYDispRight_%d", i);
86  sprintf(
87  title, "Y dispersion from central to right layer inside station: %d", i);
88  muchInStationYDispRight[i] = new TH1F(name, title, 100, -3.0, 3.0);
89  muchInStationYDispRight[i]->StatOverflows();
90 
91  sprintf(name, "muchInStationXDispRL_%d", i);
92  sprintf(
93  title, "X dispersion on left layer predicted by right station: %d", i);
94  muchInStationXDispRL[i] = new TH1F(name, title, 100, -0.1, 0.1);
95  muchInStationXDispRL[i]->StatOverflows();
96 
97  sprintf(name, "muchInStationYDispRL_%d", i);
98  sprintf(
99  title, "Y dispersion on left layer predicted by right station: %d", i);
100  muchInStationYDispRL[i] = new TH1F(name, title, 100, -0.1, 0.1);
101  muchInStationYDispRL[i]->StatOverflows();
102 
103  sprintf(name, "muchInStationTxBreak_%d", i);
104  sprintf(title, "Tx break inside station: %d", i);
105  muchInStationTxBreak[i] = new TH1F(name, title, 100, -0.02, 0.02);
106  muchInStationTxBreak[i]->StatOverflows();
107 
108  sprintf(name, "muchInStationTyBreak_%d", i);
109  sprintf(title, "Ty break inside station: %d", i);
110  muchInStationTyBreak[i] = new TH1F(name, title, 100, -0.02, 0.02);
111  muchInStationTyBreak[i]->StatOverflows();
112 
113  if (i > 0) {
114  sprintf(name, "muchLongSegmentTxHisto_%d", i);
115  sprintf(
116  title,
117  "Tx tangents distribution for segments between stations: %d and %d",
118  i - 1,
119  i);
120  muchLongSegmentTxHisto[i - 1] = new TH1F(name, title, 100, -.15, .15);
121  muchLongSegmentTxHisto[i - 1]->StatOverflows();
122 
123  sprintf(name, "muchLongSegmentTyHisto_%d", i);
124  sprintf(
125  title,
126  "Ty tangents distribution for segments between stations: %d and %d",
127  i - 1,
128  i);
129  muchLongSegmentTyHisto[i - 1] = new TH1F(name, title, 100, -.15, .15);
130  muchLongSegmentTyHisto[i - 1]->StatOverflows();
131 
132  sprintf(name, "muchClusterXDispHisto_%d", i);
133  sprintf(title,
134  "X coordinate dispersion for cluster segments between stations: "
135  "%d and %d",
136  i - 1,
137  i);
138  muchClusterXDispHisto[i - 1] = new TH1F(name, title, 100, .0, 3.0);
139  muchClusterXDispHisto[i - 1]->StatOverflows();
140 
141  sprintf(name, "muchClusterYDispHisto_%d", i);
142  sprintf(title,
143  "Y coordinate dispersion for cluster segments between stations: "
144  "%d and %d",
145  i - 1,
146  i);
147  muchClusterYDispHisto[i - 1] = new TH1F(name, title, 100, .0, 3.0);
148  muchClusterYDispHisto[i - 1]->StatOverflows();
149 
150  sprintf(name, "muchClusterTxDispHisto_%d", i);
151  sprintf(title,
152  "Tx tangent dispersion for cluster segments between stations: %d "
153  "and %d",
154  i - 1,
155  i);
156  muchClusterTxDispHisto[i - 1] = new TH1F(name, title, 100, .0, .05);
157  muchClusterTxDispHisto[i - 1]->StatOverflows();
158 
159  sprintf(name, "muchClusterTyDispHisto_%d", i);
160  sprintf(title,
161  "Ty tangent dispersion for cluster segments between stations: %d "
162  "and %d",
163  i - 1,
164  i);
165  muchClusterTyDispHisto[i - 1] = new TH1F(name, title, 100, .0, .05);
166  muchClusterTyDispHisto[i - 1]->StatOverflows();
167 
168  sprintf(name, "muchOutStationTxBreakLeft_%d", i);
169  sprintf(title,
170  "Tx break between right segment of station and left tip of the "
171  "interstation segment: %d",
172  i);
174  new TH1F(name, title, 100, -0.15, 0.15);
175  muchOutStationTxBreakLeft[i - 1]->StatOverflows();
176 
177  sprintf(name, "muchOutStationTxBreakRight_%d", i);
178  sprintf(title,
179  "Tx break between left segment of station and right tip of the "
180  "interstation segment: %d",
181  i);
183  new TH1F(name, title, 100, -0.15, 0.15);
184  muchOutStationTxBreakRight[i - 1]->StatOverflows();
185 
186  sprintf(name, "muchOutStationTyBreakLeft_%d", i);
187  sprintf(title,
188  "Ty break between right segment of station and left tip of the "
189  "interstation segment: %d",
190  i);
192  new TH1F(name, title, 100, -0.15, 0.15);
193  muchOutStationTyBreakLeft[i - 1]->StatOverflows();
194 
195  sprintf(name, "muchOutStationTyBreakRight_%d", i);
196  sprintf(title,
197  "Ty break between left segment of station and right tip of the "
198  "interstation segment: %d",
199  i);
201  new TH1F(name, title, 100, -0.15, 0.15);
202  muchOutStationTyBreakRight[i - 1]->StatOverflows();
203 
204  sprintf(name, "muchOutStationXDispByTriplet_%d", i);
205  sprintf(title,
206  "X dispersion of prediction by triplet angle for station: %d",
207  i);
209  new TH1F(name, title, 100, -10.0, 10.0);
210  muchOutStationXDispByTriplet[i - 1]->StatOverflows();
211 
212  sprintf(name, "muchOutStationYDispByTriplet_%d", i);
213  sprintf(title,
214  "Y dispersion of prediction by triplet angle for station: %d",
215  i);
217  new TH1F(name, title, 100, -10.0, 10.0);
218  muchOutStationYDispByTriplet[i - 1]->StatOverflows();
219 
220  sprintf(name, "muchOutStationXDispByVertex_%d", i);
221  sprintf(
222  title,
223  "X dispersion of prediction by an angle to vertex for station: %d",
224  i);
226  new TH1F(name, title, 100, -10.0, 10.0);
227  muchOutStationXDispByVertex[i - 1]->StatOverflows();
228 
229  sprintf(name, "muchOutStationYDispByVertex_%d", i);
230  sprintf(
231  title,
232  "Y dispersion of prediction by an angle to vertex for station: %d",
233  i);
235  new TH1F(name, title, 100, -10.0, 10.0);
236  muchOutStationYDispByVertex[i - 1]->StatOverflows();
237 
238  if (i < LXSTATIONS - 1) {
239  sprintf(name, "muchSegmentTxBreakHisto_%d", i);
240  sprintf(title,
241  "Tx tangents breaks distribution for adjacent segments on "
242  "station: %d",
243  i);
244  muchSegmentTxBreakHisto[i - 1] = new TH1F(name, title, 100, -.15, .15);
245  muchSegmentTxBreakHisto[i - 1]->StatOverflows();
246 
247  sprintf(name, "muchSegmentTyBreakHisto_%d", i);
248  sprintf(title,
249  "Ty tangents breaks distribution for adjacent segments on "
250  "station: %d",
251  i);
252  muchSegmentTyBreakHisto[i - 1] = new TH1F(name, title, 100, -.15, .15);
253  muchSegmentTyBreakHisto[i - 1]->StatOverflows();
254  }
255 
256  sprintf(name, "muchStationTxDispHisto_%d", i);
257  sprintf(title,
258  "Tx tangents dispersion for segments between stations: %d and %d",
259  i - 1,
260  i);
261  muchStationTxDispHisto[i - 1] = new TH1F(name, title, 100, -.05, .05);
262  muchStationTxDispHisto[i - 1]->StatOverflows();
263 
264  sprintf(name, "muchStationTyDispHisto_%d", i);
265  sprintf(title,
266  "Ty tangents dispersion for segments between stations: %d and %d",
267  i - 1,
268  i);
269  muchStationTyDispHisto[i - 1] = new TH1F(name, title, 100, -.05, .05);
270  muchStationTyDispHisto[i - 1]->StatOverflows();
271  }
272 
273  if (i < LXSTATIONS - 1) {
274  sprintf(name, "muchXTxCovHisto_%d", i);
275  sprintf(title, "muchXTxCovHisto on %d", i);
276  muchXTxCovHisto[i] =
277  new TH2F(name, title, 100, -5.0, 5.0, 100, -0.15, 0.15);
278  muchXTxCovHisto[i]->StatOverflows();
279 
280  sprintf(name, "muchYTyCovHisto_%d", i);
281  sprintf(title, "muchYTyCovHisto on %d", i);
282  muchYTyCovHisto[i] =
283  new TH2F(name, title, 100, -5.0, 5.0, 100, -0.15, 0.15);
284  muchYTyCovHisto[i]->StatOverflows();
285  }
286  }
287 }
288 
289 static void SaveHisto(TH1* histo) {
290  char dir_name[256];
291  sprintf(dir_name, "configuration.%s", particleType.Data());
292  DIR* dir = opendir(dir_name);
293 
294  if (dir)
295  closedir(dir);
296  else
297  mkdir(dir_name, 0700);
298 
299  char name[256];
300  sprintf(name, "%s/%s.root", dir_name, histo->GetName());
301  TFile fh(name, "RECREATE");
302  histo->Write();
303  fh.Close();
304  delete histo;
305 }
306 
308  for (Int_t i = 0; i < LXSTATIONS; ++i) {
317 
318  if (i > 0) {
321 
326 
331 
336 
337  if (i < LXSTATIONS - 1) {
340  }
341 
344  }
345 
346  if (i < LXSTATIONS - 1) {
349  }
350  }
351 }
352 
354  vector<LxSimpleTrack*>& tracks = owner.allTracks;
355 
356  for (vector<LxSimpleTrack*>::iterator i = tracks.begin(); i != tracks.end();
357  ++i) {
358  LxSimpleTrack* track = *i;
359 
360  if (0 > track->motherId && (13 == track->pdgCode || -13 == track->pdgCode))
361  StatForTrack(track);
362  }
363 }
364 
370 
371  LxSimpleSegment() : tx(0), ty(0) {}
373  : source(s)
374  , end(e)
375  , tx((e.x - s.x) / (e.z - s.z))
376  , ty((e.y - s.y) / (e.z - s.z)) {}
377 };
378 
380  for (Int_t i = 0; i < stationsInAlgo; ++i) {
381  for (Int_t j = 0; j < LXLAYERS; ++j) {
382  if (track->muchPoints[i][j].empty()) return;
383  }
384  }
385 
386  LxSimplePoint p1;
387  LxSimplePoint p2;
388  LxSimplePoint p3;
389  scaltype deltaZ;
390  scaltype deltaZ2;
391  scaltype tx;
392  scaltype tx2;
393  scaltype ty;
394  scaltype ty2;
395  scaltype stTx;
396  scaltype stTy;
397  scaltype stTxP;
398  scaltype stTyP;
399 
400  for (Int_t i = 0; i < LXSTATIONS; ++i) {
401  if (track->muchPoints[i][0].empty() || track->muchPoints[i][1].empty()
402  || track->muchPoints[i][2].empty())
403  continue;
404 
405  p1 = track->muchPoints[i][1].front();
406  scaltype txEst = p1.x / p1.z;
407  scaltype tyEst = p1.y / p1.z;
408 
409  p2 = track->muchPoints[i][0].front();
410  deltaZ = p2.z - p1.z;
411  scaltype xEst = p1.x + txEst * deltaZ;
412  scaltype yEst = p1.y + tyEst * deltaZ;
413  muchInStationXDispLeft[i]->Fill(p2.x - xEst);
414  muchInStationYDispLeft[i]->Fill(p2.y - yEst);
415  tx = (p2.x - p1.x) / deltaZ;
416  ty = (p2.y - p1.y) / deltaZ;
417 
418  p2 = track->muchPoints[i][2].front();
419  deltaZ = p2.z - p1.z;
420  xEst = p1.x + txEst * deltaZ;
421  yEst = p1.y + tyEst * deltaZ;
422  muchInStationXDispRight[i]->Fill(p2.x - xEst);
423  muchInStationYDispRight[i]->Fill(p2.y - yEst);
424 
425  tx2 = (p2.x - p1.x) / deltaZ;
426  ty2 = (p2.y - p1.y) / deltaZ;
427 
428  muchInStationTxBreak[i]->Fill(tx2 - tx);
429  muchInStationTyBreak[i]->Fill(ty2 - ty);
430 
431  stTxP = stTx;
432  stTyP = stTy;
433  p1 = track->muchPoints[i][0].front();
434  deltaZ = p2.z - p1.z;
435  stTx = (p2.x - p1.x) / deltaZ;
436  stTy = (p2.y - p1.y) / deltaZ;
437  muchInStationXDispRL[i]->Fill(p1.x - p2.x + tx2 * deltaZ);
438  muchInStationYDispRL[i]->Fill(p1.y - p2.y + ty2 * deltaZ);
439 
440  if (i > 0) {
441  p1 = track->muchPoints[i - 1][LXMIDDLE].front();
442  p2 = track->muchPoints[i][LXMIDDLE].front();
443  deltaZ = p2.z - p1.z;
444  tx = (p2.x - p1.x) / deltaZ;
445  muchLongSegmentTxHisto[i - 1]->Fill(tx - p2.x / p2.z);
446  ty = (p2.y - p1.y) / deltaZ;
447  muchLongSegmentTyHisto[i - 1]->Fill(ty - p2.y / p2.z);
448 
449  muchOutStationTxBreakRight[i - 1]->Fill(tx - stTx);
450  muchOutStationTyBreakRight[i - 1]->Fill(ty - stTy);
451 
452  muchOutStationTxBreakLeft[i - 1]->Fill(tx - stTxP);
453  muchOutStationTyBreakLeft[i - 1]->Fill(ty - stTyP);
454 
455  muchOutStationXDispByTriplet[i - 1]->Fill(p1.x - p2.x + stTx * deltaZ);
456  muchOutStationYDispByTriplet[i - 1]->Fill(p1.y - p2.y + stTy * deltaZ);
457 
458  muchOutStationXDispByVertex[i - 1]->Fill(p1.x - p2.x
459  + (p2.x / p2.z) * deltaZ);
460  muchOutStationYDispByVertex[i - 1]->Fill(p1.y - p2.y
461  + (p2.y / p2.z) * deltaZ);
462 
463  // Rather complex part for implementation: calculate the dispersion characteristics for segment clusters.
464  scaltype maxXdisp = 0;
465  scaltype maxYdisp = 0;
466  scaltype maxTxdisp = 0;
467  scaltype maxTydisp = 0;
468 
469  for (list<LxSimplePoint>::iterator l0 =
470  track->muchPoints[i - 1][0].begin();
471  l0 != track->muchPoints[i - 1][0].end();
472  ++l0) {
473  for (list<LxSimplePoint>::iterator l1 =
474  track->muchPoints[i - 1][1].begin();
475  l1 != track->muchPoints[i - 1][1].end();
476  ++l1) {
477  for (list<LxSimplePoint>::iterator l2 =
478  track->muchPoints[i - 1][2].begin();
479  l2 != track->muchPoints[i - 1][2].end();
480  ++l2) {
481  for (list<LxSimplePoint>::iterator r0 =
482  track->muchPoints[i][0].begin();
483  r0 != track->muchPoints[i][0].end();
484  ++r0) {
485  for (list<LxSimplePoint>::iterator r1 =
486  track->muchPoints[i][1].begin();
487  r1 != track->muchPoints[i][1].end();
488  ++r1) {
489  for (list<LxSimplePoint>::iterator r2 =
490  track->muchPoints[i][2].begin();
491  r2 != track->muchPoints[i][2].end();
492  ++r2) {
493  LxSimplePoint lPoints[LXLAYERS] = {*l0, *l1, *l2};
494  LxSimplePoint rPoints[LXLAYERS] = {*r0, *r1, *r2};
495  LxSimpleSegment segments[LXLAYERS * LXLAYERS] = {};
496 
497  for (Int_t j = 0; j < LXLAYERS; ++j) {
498  for (Int_t k = 0; k < LXLAYERS; ++k)
499  segments[j * LXLAYERS + k] =
500  LxSimpleSegment(rPoints[j], lPoints[k]);
501  }
502 
503  for (Int_t j = 0; j < LXLAYERS * LXLAYERS - 1; ++j) {
504  LxSimpleSegment s1 = segments[j];
505 
506  for (Int_t k = j + 1; k < LXLAYERS * LXLAYERS; ++k) {
507  LxSimpleSegment s2 = segments[k];
508  scaltype diffZ = s1.source.z - s2.source.z;
509  scaltype dtx = abs(s2.tx - s1.tx);
510  scaltype dty = abs(s2.ty - s1.ty);
511  scaltype dx =
512  abs(s2.source.x + s2.tx * diffZ - s1.source.x);
513  scaltype dy =
514  abs(s2.source.y + s2.ty * diffZ - s1.source.y);
515 
516  if (maxXdisp < dx) maxXdisp = dx;
517 
518  if (maxYdisp < dy) maxYdisp = dy;
519 
520  if (maxTxdisp < dtx) maxTxdisp = dtx;
521 
522  if (maxTydisp < dty) maxTydisp = dty;
523  }
524  }
525  }
526  }
527  }
528  }
529  }
530  } // for (list<LxSimplePoint>::iterator l0 = track->muchPoints[i - 1][0].begin(); l0 != track->muchPoints[i - 1][0].end(); ++l0)
531 
532  muchClusterXDispHisto[i - 1]->Fill(maxXdisp);
533  muchClusterYDispHisto[i - 1]->Fill(maxYdisp);
534  muchClusterTxDispHisto[i - 1]->Fill(maxTxdisp);
535  muchClusterTyDispHisto[i - 1]->Fill(maxTydisp);
536 
537  if (i < LXSTATIONS - 1) {
538  p1 = track->muchPoints[i - 1][LXMIDDLE].front();
539  p2 = track->muchPoints[i][LXMIDDLE].front();
540  deltaZ = p2.z - p1.z;
541  tx = (p2.x - p1.x) / deltaZ;
542  ty = (p2.y - p1.y) / deltaZ;
543  p3 = track->muchPoints[i + 1][LXMIDDLE].front();
544  deltaZ2 = p3.z - p2.z;
545  tx2 = (p3.x - p2.x) / deltaZ2;
546  muchSegmentTxBreakHisto[i - 1]->Fill(tx2 - tx);
547  ty2 = (p3.y - p2.y) / deltaZ2;
548  muchSegmentTyBreakHisto[i - 1]->Fill(ty2 - ty);
549  }
550 
551  for (Int_t j = 0; j < LXLAYERS * LXLAYERS; ++j) {
552  p1 = track->muchPoints[i - 1][j % LXLAYERS].front();
553  p2 = track->muchPoints[i][j / LXLAYERS].front();
554  deltaZ = p2.z - p1.z;
555  tx = (p2.x - p1.x) / deltaZ;
556  ty = (p2.y - p1.y) / deltaZ;
557 
558  for (Int_t k = j + 1; k < LXLAYERS * LXLAYERS; ++k) {
559  p1 = track->muchPoints[i - 1][k % LXLAYERS].front();
560  p2 = track->muchPoints[i][k / LXLAYERS].front();
561  deltaZ = p2.z - p1.z;
562  tx2 = (p2.x - p1.x) / deltaZ;
563  muchStationTxDispHisto[i - 1]->Fill(tx2 - tx);
564  ty2 = (p2.y - p1.y) / deltaZ;
565  muchStationTyDispHisto[i - 1]->Fill(ty2 - ty);
566  }
567  }
568  }
569 
570  if (i < LXSTATIONS - 1) {
571  p1 = track->muchPoints[i + 1][LXMIDDLE].front();
572  p2 = track->muchPoints[i][LXMIDDLE].front();
573  deltaZ = p2.z - p1.z;
574  scaltype deltaX = p2.x - p1.x - p1.tx * deltaZ;
575  scaltype deltaY = p2.y - p1.y - p1.ty * deltaZ;
576  scaltype deltaTx = p2.tx - p1.tx;
577  scaltype deltaTy = p2.ty - p1.ty;
578  muchXTxCovHisto[i]->Fill(deltaX, deltaTx);
579  muchYTyCovHisto[i]->Fill(deltaY, deltaTy);
580  }
581  }
582 }
LxSimpleSegment::ty
scaltype ty
Definition: Simple/LxTrackAnaSegments.cxx:369
muchYTyCovHisto
static TH2F * muchYTyCovHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:44
LxSimpleTrack::motherId
Int_t motherId
Definition: Simple/LxTrackAna.h:27
muchStationTxDispHisto
static TH1F * muchStationTxDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:40
muchOutStationYDispByTriplet
static TH1F * muchOutStationYDispByTriplet[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:29
muchLongSegmentTxHisto
static TH1F * muchLongSegmentTxHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:34
SaveHisto
static void SaveHisto(TH1 *histo)
Definition: Simple/LxTrackAnaSegments.cxx:289
scaltype
#define scaltype
Definition: CbmGlobalTrackingDefs.h:17
muchSegmentTxBreakHisto
static TH1F * muchSegmentTxBreakHisto[LXSTATIONS - 2]
Definition: Simple/LxTrackAnaSegments.cxx:37
muchInStationXDispRight
static TH1F * muchInStationXDispRight[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:13
LxTrackAnaSegments::owner
LxTrackAna & owner
Definition: Simple/LxTrackAnaSegments.h:21
muchOutStationTxBreakRight
static TH1F * muchOutStationTxBreakRight[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:24
muchClusterYDispHisto
static TH1F * muchClusterYDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:47
LxTrackAnaSegments::SetParticleType
void SetParticleType(TString v)
Definition: Simple/LxTrackAnaSegments.cxx:56
LxSimpleSegment::source
LxSimplePoint source
Definition: Simple/LxTrackAnaSegments.cxx:366
muchInStationXDispLeft
static TH1F * muchInStationXDispLeft[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:12
muchStationTyDispHisto
static TH1F * muchStationTyDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:41
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
muchInStationXDispRL
static TH1F * muchInStationXDispRL[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:17
LXLAYERS
#define LXLAYERS
Definition: Simple/LxSettings.h:8
LxSimplePoint::y
scaltype y
Definition: Simple/LxTrackAna.h:16
LXMIDDLE
#define LXMIDDLE
Definition: Simple/LxSettings.h:10
LxSimpleTrack
Definition: Simple/LxTrackAna.h:25
muchOutStationXDispByVertex
static TH1F * muchOutStationXDispByVertex[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:31
muchOutStationTyBreakRight
static TH1F * muchOutStationTyBreakRight[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:26
muchInStationYDispLeft
static TH1F * muchInStationYDispLeft[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:14
muchOutStationXDispByTriplet
static TH1F * muchOutStationXDispByTriplet[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:28
LxTrackAna.h
particleType
static TString particleType("jpsi")
LxSimpleSegment::tx
scaltype tx
Definition: Simple/LxTrackAnaSegments.cxx:368
LxSimplePoint::z
scaltype z
Definition: Simple/LxTrackAna.h:17
muchOutStationTxBreakLeft
static TH1F * muchOutStationTxBreakLeft[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:23
tracks
TClonesArray * tracks
Definition: Analyze_matching.h:17
LxTrackAnaSegments::BuildStatistics
void BuildStatistics()
Definition: Simple/LxTrackAnaSegments.cxx:353
muchInStationTxBreak
static TH1F * muchInStationTxBreak[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:20
LxSimpleSegment
Definition: Simple/LxTrackAnaSegments.cxx:365
LxTrackAna::allTracks
std::vector< LxSimpleTrack * > allTracks
Definition: Simple/LxTrackAna.h:118
muchInStationYDispRL
static TH1F * muchInStationYDispRL[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:18
muchOutStationYDispByVertex
static TH1F * muchOutStationYDispByVertex[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:32
LxTrackAnaSegments::StatForTrack
void StatForTrack(LxSimpleTrack *track)
Definition: Simple/LxTrackAnaSegments.cxx:379
muchClusterXDispHisto
static TH1F * muchClusterXDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:46
muchSegmentTyBreakHisto
static TH1F * muchSegmentTyBreakHisto[LXSTATIONS - 2]
Definition: Simple/LxTrackAnaSegments.cxx:38
muchClusterTxDispHisto
static TH1F * muchClusterTxDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:48
LxSimplePoint::x
scaltype x
Definition: Simple/LxTrackAna.h:15
LxTrackAnaSegments::stationsInAlgo
Int_t stationsInAlgo
Definition: Simple/LxTrackAnaSegments.h:22
muchInStationYDispRight
static TH1F * muchInStationYDispRight[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:15
LxTrackAnaSegments::Init
void Init()
Definition: Simple/LxTrackAnaSegments.cxx:62
LxTrackAnaSegments.h
LxTrackAna
Definition: Simple/LxTrackAna.h:66
muchXTxCovHisto
static TH2F * muchXTxCovHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:43
v
__m128 v
Definition: L1/vectors/P4_F32vec4.h:1
LxSimplePoint::ty
scaltype ty
Definition: Simple/LxTrackAna.h:19
x
Double_t x
Definition: CbmMvdSensorDigiToHitTask.cxx:68
muchInStationTyBreak
static TH1F * muchInStationTyBreak[LXSTATIONS]
Definition: Simple/LxTrackAnaSegments.cxx:21
y
Double_t y
Definition: CbmMvdSensorDigiToHitTask.cxx:68
LxSimpleSegment::end
LxSimplePoint end
Definition: Simple/LxTrackAnaSegments.cxx:367
LxSimpleTrack::muchPoints
std::list< LxSimplePoint > muchPoints[LXSTATIONS][LXLAYERS]
Definition: Simple/LxTrackAna.h:55
LxTrackAnaSegments::LxTrackAnaSegments
LxTrackAnaSegments(LxTrackAna &o)
Definition: Simple/LxTrackAnaSegments.cxx:51
LxSimplePoint::tx
scaltype tx
Definition: Simple/LxTrackAna.h:18
LxSimpleSegment::LxSimpleSegment
LxSimpleSegment()
Definition: Simple/LxTrackAnaSegments.cxx:371
LxTrackAnaSegments::Finish
void Finish()
Definition: Simple/LxTrackAnaSegments.cxx:307
LXSTATIONS
#define LXSTATIONS
Definition: Simple/LxSettings.h:9
muchLongSegmentTyHisto
static TH1F * muchLongSegmentTyHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:35
LxSimpleTrack::pdgCode
Int_t pdgCode
Definition: Simple/LxTrackAna.h:26
muchOutStationTyBreakLeft
static TH1F * muchOutStationTyBreakLeft[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:25
LxSimpleSegment::LxSimpleSegment
LxSimpleSegment(LxSimplePoint s, LxSimplePoint e)
Definition: Simple/LxTrackAnaSegments.cxx:372
muchClusterTyDispHisto
static TH1F * muchClusterTyDispHisto[LXSTATIONS - 1]
Definition: Simple/LxTrackAnaSegments.cxx:49
LxSimplePoint
Definition: Simple/LxTrackAna.h:14