CbmRoot
CbmTrdModuleRecT.cxx
Go to the documentation of this file.
1 #include "CbmTrdModuleRecT.h"
2 #include "CbmTrdCluster.h"
3 #include "CbmTrdDigi.h"
4 #include "CbmTrdFASP.h"
5 #include "CbmTrdHit.h"
6 #include "CbmTrdParModDigi.h"
7 
8 #include <FairLogger.h>
9 
10 #include <TGeoPhysicalNode.h>
11 
12 #include <TClonesArray.h>
13 #include <TF1.h>
14 #include <TGraphErrors.h>
15 
16 #include <iostream>
17 #include <vector>
18 
19 using std::cout;
20 using std::endl;
21 using std::vector;
22 
23 Double_t CbmTrdModuleRecT::fgDT[] = {4.181e-6, 1586, 24};
24 TGraphErrors* CbmTrdModuleRecT::fgEdep = NULL;
25 TGraphErrors* CbmTrdModuleRecT::fgT = NULL;
26 TF1* CbmTrdModuleRecT::fgPRF = NULL;
27 //_______________________________________________________________________________
29  : CbmTrdModuleRec()
30  , fConfigMap(0)
31  , fT0(0)
32  , fBuffer()
33  , vs()
34  , vse()
35  , vt()
36  , vx()
37  , vxe() {}
38 
39 //_______________________________________________________________________________
40 CbmTrdModuleRecT::CbmTrdModuleRecT(Int_t mod, Int_t ly, Int_t rot)
41  : CbmTrdModuleRec(mod, ly, rot)
42  , fConfigMap(0)
43  , fT0(0)
44  , fBuffer()
45  , vs()
46  , vse()
47  , vt()
48  , vx()
49  , vxe() {
50  //printf("AddModuleT @ %d\n", mod); Config(1,0);
51  SetNameTitle(Form("TrdRecT%d", mod),
52  "Reconstructor for triangular read-out.");
53 }
54 
55 //_______________________________________________________________________________
57 
58 //_______________________________________________________________________________
59 Bool_t CbmTrdModuleRecT::AddDigi(const CbmTrdDigi* d, Int_t id) {
60  /* Add digi to cluster fragments. At first clusters are ordered on pad rows and time.
61  * No channel ordering is assumed. The time condition for a digi to enter a cluster
62  * chunk is to have abs(dt)<5 wrt cluster t0
63  */
64 
65  if (CWRITE()) cout << "add @" << id << " " << d->ToString();
66 
67  Int_t ch = d->GetAddressChannel(), col, row = GetPadRowCol(ch, col), dtime;
68  Double_t t, r = d->GetCharge(t, dtime);
69  Int_t tm = d->GetTimeDAQ() - fT0, terminator(0);
70  if (dtime < 0)
71  tm += dtime; // correct for the difference between tilt and rect
72  if (r < 1)
73  terminator = 1;
74  else if (t < 1)
75  terminator = -1;
76 
77  if (CWRITE())
78  printf(
79  "row[%2d] col[%2d] tm[%2d] terminator[%d]\n", row, col, tm, terminator);
80  CbmTrdCluster* cl(NULL);
81 
82  // get the link to saved clusters
83  std::map<Int_t, std::list<CbmTrdCluster*>>::iterator it = fBuffer.find(row);
84  // check for saved
85  if (it != fBuffer.end()) {
86  Bool_t kINSERT(kFALSE);
87  std::list<CbmTrdCluster*>::iterator itc = fBuffer[row].begin();
88  for (; itc != fBuffer[row].end(); itc++) {
89  if (CWRITE()) cout << (*itc)->ToString();
90 
91  UShort_t tc = (*itc)->GetStartTime();
92  Int_t dt = tc - tm;
93 
94  if (dt < -5)
95  continue;
96  else if (dt < 5) {
97  if (CWRITE()) printf("match time dt=%d\n", dt);
98  if ((*itc)->IsChannelInRange(ch) == 0) {
99  if (!(*itc)->AddDigi(id, ch, terminator, dt)) break;
100  kINSERT = kTRUE;
101  if (CWRITE()) cout << " => Cl (Ad) " << (*itc)->ToString();
102  break;
103  }
104  } else {
105  if (CWRITE()) printf("break for time dt=%d\n", dt);
106  break;
107  }
108  }
109 
110  if (!kINSERT) {
111  if (itc != fBuffer[row].end() && itc != fBuffer[row].begin()) {
112  itc--;
113  fBuffer[row].insert(
114  itc, cl = new CbmTrdCluster(fModAddress, id, ch, row, tm));
115  if (CWRITE()) cout << " => Cl (I) " << cl->ToString();
116  } else {
117  fBuffer[row].push_back(
118  cl = new CbmTrdCluster(fModAddress, id, ch, row, tm));
119  if (CWRITE()) cout << " => Cl (Pb) " << cl->ToString();
120  }
121  cl->SetTrianglePads();
122  if (terminator < 0)
123  cl->SetProfileStart();
124  else if (terminator > 0)
125  cl->SetProfileStop();
126  }
127  } else {
128  fBuffer[row].push_back(cl =
129  new CbmTrdCluster(fModAddress, id, ch, row, tm));
130  cl->SetTrianglePads();
131  if (terminator < 0)
132  cl->SetProfileStart();
133  else if (terminator > 0)
134  cl->SetProfileStop();
135  if (CWRITE()) cout << " => Cl (Nw) " << cl->ToString();
136  }
137 
138  return kTRUE;
139 }
140 
141 //_______________________________________________________________________________
143  Int_t nch(0);
144  for (std::map<Int_t, std::list<CbmTrdCluster*>>::const_iterator ir =
145  fBuffer.cbegin();
146  ir != fBuffer.cend();
147  ir++) {
148  for (std::list<CbmTrdCluster*>::const_iterator itc = (*ir).second.cbegin();
149  itc != (*ir).second.cend();
150  itc++)
151  nch += (*itc)->GetNCols();
152  }
153  return nch;
154 }
155 
156 //_______________________________________________________________________________
158  CbmTrdCluster* cl(NULL);
159 
160  // get the link to saved clusters
161  Int_t ncl(0);
162  std::list<CbmTrdCluster*>::iterator itc0, itc1;
163  for (std::map<Int_t, std::list<CbmTrdCluster*>>::iterator ir =
164  fBuffer.begin();
165  ir != fBuffer.end();
166  ir++) {
167  itc0 = (*ir).second.begin();
168  while ((*ir).second.size() > 1
169  && itc0 != (*ir).second.end()) { // try merge clusters
170  itc1 = itc0;
171  itc1++;
172 
173  Bool_t kMERGE = kFALSE;
174  while (itc1 != (*ir).second.end()) {
175  if (CWRITE())
176  cout << " base cl[0] : " << (*itc0)->ToString()
177  << " + cl[1] : " << (*itc1)->ToString();
178  if ((*itc0)->Merge((*itc1))) {
179  if (CWRITE()) cout << " SUM : " << (*itc0)->ToString();
180  kMERGE = kTRUE;
181  delete (*itc1);
182  itc1 = (*ir).second.erase(itc1);
183  if (itc1 == (*ir).second.end()) break;
184  } else
185  itc1++;
186  }
187  if (!kMERGE) itc0++;
188  }
189 
190  for (itc0 = (*ir).second.begin(); itc0 != (*ir).second.end(); itc0++) {
191  cl = (*itc0);
192  cl = new ((*fClusters)[ncl++]) CbmTrdCluster(*cl);
193  cl->SetTrianglePads();
194  delete (*itc0);
195  }
196  }
197  fBuffer.clear();
198 
199  //printf("fClusters[%p] nCl[%d]\n", (void*)fClusters, fClusters->GetEntriesFast());
200  //LOG(info) << GetName() << "::FindClusters : " << ncl;
201  return ncl;
202 }
203 
204 //_______________________________________________________________________________
205 Bool_t CbmTrdModuleRecT::MakeHits() { return kTRUE; }
206 
207 //_______________________________________________________________________________
209  /* Steering routine for classifying hits and apply further analysis
210  * -> hit deconvolution (see Deconvolute)
211  * -> hit merging for row-cross (see RowCross)
212  */
213 
214  // Int_t nhits=fHits->GetEntriesFast();
215  // //if(CWRITE())
216  // LOG(info) << "\n"<<GetName()<<"::Finalize("<<nhits<<")";
217  // CbmTrdHit *hit(NULL), *hitp(NULL);
218  // for(Int_t ih(0); ih<nhits; ih++){
219  // hit = (CbmTrdHit*)((*fHits)[ih]);
220  // for(Int_t jh(ih+1); jh<nhits; jh++){
221  // hitp = (CbmTrdHit*)((*fHits)[jh]);
222  // //if(CWRITE())
223  // cout<<ih<<" "<<hit->ToString();
224  // cout<<"-> "<<jh<<" "<<hitp->ToString();
225  // }
226  // }
227  return kTRUE;
228 }
229 
230 #include <TCanvas.h>
231 #include <TH1.h>
232 #define NANODE 9
233 //_______________________________________________________________________________
235  const CbmTrdCluster* cl,
236  std::vector<const CbmTrdDigi*>* digis) {
237  if (!fgEdep) { // first use initialization of PRF helppers
238  LOG(info) << GetName() << "::MakeHit: Init static helpers. ";
239  fgEdep = new TGraphErrors;
240  fgEdep->SetLineColor(kBlue);
241  fgEdep->SetLineWidth(2);
242  fgT = new TGraphErrors;
243  fgT->SetLineColor(kBlue);
244  fgT->SetLineWidth(2);
245  fgPRF = new TF1("prf", "[0]*exp(-0.5*((x-[1])/[2])**2)", -10, 10);
246  fgPRF->SetLineColor(kRed);
247  fgPRF->SetParNames("E", "x", "prf");
248  }
249  TH1* hf(NULL);
250  TCanvas* cvs(NULL);
251  if (CDRAW()) {
252  cvs = new TCanvas("c", "TRD Anode Hypothesis", 10, 600, 1000, 500);
253  cvs->Divide(2, 1, 1.e-5, 1.e-5);
254  TVirtualPad* vp = cvs->cd(1);
255  vp->SetLeftMargin(0.06904908);
256  vp->SetRightMargin(0.009969325);
257  vp->SetTopMargin(0.01402806);
258  vp = cvs->cd(2);
259  vp->SetLeftMargin(0.06904908);
260  vp->SetRightMargin(0.009969325);
261  vp->SetTopMargin(0.01402806);
262  hf = new TH1I("hf", ";x [pw];s [ADC chs]", 100, -3, 3);
263  hf->GetYaxis()->SetRangeUser(-50, 4500);
264  hf->Draw("p");
265  }
266 
267  if (CWRITE()) cout << cl->ToString();
268 
269  // check cluster integrity and do digi merging if needed
270  vector<Bool_t> vmask(digis->size(), 0); // masks in case of merged Digis
271  vector<CbmTrdDigi*> vdgM;
272  if (cl->GetNCols() != digis->size() && !MergeDigis(digis, &vdgM, &vmask)) {
273  cout << cl->ToString();
274  for (vector<const CbmTrdDigi*>::iterator i = digis->begin();
275  i != digis->end();
276  i++)
277  cout << (*i)->ToString();
278  return NULL;
279  }
280 
281  // Read in all digis information()
282  ULong64_t t0 = fT0 + cl->GetStartTime(); // absolute hit time (prompt signal)
283  Int_t n0(0), ovf(0), cM;
284  if (!(n0 = LoadDigis(digis, &vdgM, &vmask, t0, cM))) {
285  cout << cl->ToString();
286  for (vector<const CbmTrdDigi*>::iterator i = digis->begin();
287  i != digis->end();
288  i++)
289  cout << (*i)->ToString();
290  return NULL;
291  }
292  if (n0 < 0) {
293  ovf = 1;
294  n0 *= -1;
295  }
296 
297  // analyse digis topology; no of signal types, maximum, etc
298  Bool_t tM(kTRUE); // maximum type tilt=1; rect=0
299  Int_t nL(0); // signal index for the max signal
300  Int_t col, row = GetPadRowCol(cl->GetStartCh(), col);
301  Double_t max(0.), // maximum signal
302  LS(0.), // left side sum of signals
303  S(0.); // sum of signals
304  Int_t nr(0), nt(0);
305  for (Int_t is(1); is <= n0; is++) {
306  if (vxe[is] > 0) { // select tilted coupling
307  nt++;
308  S += vs[is];
309  if (vs[is] > max) {
310  max = vs[is];
311  tM = kTRUE;
312  nL = is;
313  LS += vs[is];
314  }
315  } else { // select rectangular coupling
316  nr++;
317  S += vs[is];
318  if (vs[is] > max) {
319  max = vs[is];
320  tM = kFALSE;
321  nL = is;
322  LS += vs[is];
323  }
324  }
325  }
326  col += cM;
327  S -= LS;
328  LS -= max;
329  // evaluate asymmetry
330  Int_t lr(0), // max signal left-right asymmetry wrt central pad
331  tr(0), // left-right innequality for symmetric clusters
332  nR(n0 + 1 - nL); // no of signals to the right of maximum
333  if (nL < nR)
334  lr = 1;
335  else if (nL > nR)
336  lr = -1;
337  if (!lr && (n0 % 2)) tr = (LS < S ? -1 : 1);
338  // compute x and y offset from center pad
339  Double_t dx(0.), dy(0.), edx(0.21650635),
340  edy(0.77942286); // fixed error parametrization
341  switch (n0) {
342  case 1:
343  if (nt) {
344  dx = -0.5;
345  dy = 0;
346  } // T
347  else {
348  dx = 0.5;
349  dy = 0;
350  } // R
351  break;
352  case 2:
353  if (cl->HasOpenStart() && cl->HasOpenStop()) {
354  dx = cM ? -1. : 0.;
355  dy = -0.5;
356  } // RT
357  else {
358  dx = 0.;
359  dy = 0.5;
360  } // TR
361  break;
362  case 3:
363  if (tM && lr) {
364  dx = cM ? -1. : 0.;
365  dy = GetYoffset(n0);
366  } // TRT asymm
367  else if (!tM && !lr) {
368  dx = 0.;
369  dy = GetYoffset(n0);
370  } // TRT symm
371  else if (tM && !lr) {
372  dx = GetXoffset(n0);
373  dy = 0.;
374  } // RTR symm
375  else if (!tM && lr) {
376  dx = GetXoffset(n0);
377  dy = cM ? 0.5 : -0.5;
378  } // RTR asymm
379  break;
380  default:
381  dx = GetXoffset(n0);
382  dy = GetYoffset(n0);
383  break;
384  }
385  if (dx < -0.5
386  && cM > 0) { // shift graph representation to fit dx[pw] in [-0.5, 0.5]
387  dx += 1.;
388  col -= 1;
389  for (UInt_t idx(0); idx < vx.size(); idx++)
390  vx[idx] += 1;
391  }
392  if (dx > 0.5) { // dirty fix for compound clusters TODO
393  Int_t ishift = Int_t(dx - 0.5) + 1;
394  dx -= ishift;
395  col += ishift;
396  for (UInt_t idx(0); idx < vx.size(); idx++)
397  vx[idx] -= ishift;
398  }
399  dy = dx - dy; // only on natural scalling !
400  // go to cm scale
401  dx *= fDigiPar->GetPadSizeX(0);
402  dy *= fDigiPar->GetPadSizeY(0);
403 
404  // apply systematic correction for x (MC derived)
405  Int_t typ = 0; // [0] center hit type
406  // [1] side hit type
407  if ((n0 == 5 && ((tM && !lr) || (!tM && lr))) || // TRTRT symm/asymm
408  n0 == 4 || (n0 == 3 && ((tM && !lr) || (!tM && lr != 0))))
409  typ = 1; // RTR symm/asymm
410  Double_t xcorr(0.);
411  Int_t nbins((NBINSCORRX - 1) >> 1), ii = nbins + TMath::Nint(dx / fgCorrXdx),
412  i0, i1;
413  if (ii < 0 || ii >= NBINSCORRX)
414  LOG(warn) << GetName() << "::MakeHit : Idx " << ii
415  << " outside range for displacement " << dx << ".";
416  else {
417  if (dx < fgCorrXdx * ii) {
418  i0 = TMath::Max(0, ii - 1);
419  i1 = ii;
420  } else {
421  i0 = ii;
422  i1 = TMath::Min(NBINSCORRX - 1, ii + 1);
423  }
424  Double_t DDx = (fgCorrXval[typ][i1] - fgCorrXval[typ][i0]),
425  a = DDx / fgCorrXdx, b = fgCorrXval[typ][i0] - DDx * (i0 - nbins);
426  xcorr = 0.1 * (b + a * dx);
427  }
428  dx += xcorr;
429  dy += xcorr;
430  if (dx > 0.5 * fDigiPar->GetPadSizeX(0))
431  dx = 0.5 * fDigiPar->GetPadSizeX(0);
432  else if (dx < -0.5 * fDigiPar->GetPadSizeX(0))
433  dx = -0.5 * fDigiPar->GetPadSizeX(0);
434 
435  if (dy > 0.5 * fDigiPar->GetPadSizeY(0)) { //
436  //printf("BEFORE dy[%+6.4f] dx[%+6.4f] {n[%d] max[%c] lr[%+d]}\n", dy, dx, n0, tM?'T':'R', lr);
437  dy -= fDigiPar->GetPadSizeY(0);
438  }
439  if (dy < -0.5 * fDigiPar->GetPadSizeY(0)) { //
440  //printf("(BEFORE) dy[%+6.4f] dx[%+6.4f] {n[%d] max[%c] lr[%+d]}\n", dy, dx, n0, tM?'T':'R', lr);
441  dy += fDigiPar->GetPadSizeY(0);
442  }
443 
444  // process y offset
445  // apply systematic correction for y (MC derived)
446  Float_t fdy(1.), yoff(0.);
447  switch (n0) {
448  case 3:
449  fdy = fgCorrYval[0][0];
450  yoff = fgCorrYval[0][1];
451  if (tM && !lr)
452  dy -= tr * 0.5 * fDigiPar->GetPadSizeY(0);
453  else if (lr)
454  dy -= 0.5 * fDigiPar->GetPadSizeY(0);
455  ;
456  break;
457  case 4:
458  fdy = fgCorrYval[1][0];
459  yoff = fgCorrYval[1][1];
460  if ((!tM && lr == 1) || (tM && lr == -1)) yoff *= -1;
461  break;
462  case 5:
463  case 7:
464  case 9:
465  case 11:
466  fdy = fgCorrYval[2][0];
467  yoff = fgCorrYval[2][1];
468  break;
469  case 6:
470  case 8:
471  case 10:
472  fdy = fgCorrYval[3][0];
473  yoff = fgCorrYval[3][1];
474  if ((!tM && lr == 1) || (tM && lr == -1)) yoff *= -1;
475  break;
476  }
477  dy *= fdy;
478  dy += yoff;
479  if (dy > 0.5 * fDigiPar->GetPadSizeY(0))
480  dy = 0.5 * fDigiPar->GetPadSizeY(0);
481  else if (dy < -0.5 * fDigiPar->GetPadSizeY(0))
482  dy = -0.5 * fDigiPar->GetPadSizeY(0);
483 
484  // get anode wire offset
485  Int_t ia(0);
486  Float_t ya; // anode position in local pad coordinates
487  for (; ia <= NANODE; ia++) {
488  ya = -1.35 + ia * 0.3;
489  if (dy > ya + 0.15) continue;
490  break;
491  }
492 
493  // Error parametrization X : parabolic model on cl size
494  Double_t parX[] = {0.713724, -0.318667, 0.0366036};
495  Int_t nex = TMath::Min(n0, 7);
496  edx = parX[0] + parX[1] * nex + parX[2] * nex * nex;
497  Double_t parY[] = {0.0886413, 0., 0.0435141};
498  edy = parY[0] + parY[2] * dy * dy;
499 
500  if (CWRITE()) {
501  printf("row[%2d] col[%2d] sz[%d%c] M[%d%c] dx[mm]=%6.3f dy[mm]=%6.3f "
502  "t0[clk]=%llu OVF[%c]\n",
503  row,
504  col,
505  n0,
506  (lr ? (lr < 0 ? 'L' : 'R') : 'C'),
507  cM,
508  (tM ? 'T' : 'R'),
509  10 * dx,
510  10 * dy,
511  t0,
512  (ovf ? 'y' : 'n'));
513  for (UInt_t idx(0); idx < vs.size(); idx++) {
514  printf("%2d%cdt[%2d] s[ADC]=%6.1f+-%6.1f x[pw]=%5.2f+-%5.2f\n",
515  idx,
516  (UInt_t(nL) == idx ? '*' : ' '),
517  vt[idx],
518  vs[idx],
519  vse[idx],
520  vx[idx],
521  vxe[idx]);
522  }
523  }
524 
525  // compute energy
526  for (UInt_t idx(0); idx < vs.size(); idx++) {
527  if (vxe[idx] > 0) {
528  fgEdep->SetPoint(idx, vx[idx] + dy / fDigiPar->GetPadSizeY(0), vs[idx]);
529  fgEdep->SetPointError(idx, vxe[idx], vse[idx]);
530  } else {
531  fgEdep->SetPoint(idx, vx[idx], vs[idx]);
532  fgEdep->SetPointError(idx, vxe[idx], vse[idx]);
533  }
534  }
535  Double_t xc = vx[n0 + 2];
536  for (Int_t ip(vs.size()); ip < fgEdep->GetN(); ip++) {
537  //fgEdep->RemovePoint(ip);
538  xc += 0.5;
539  fgEdep->SetPoint(ip, xc, 0);
540  fgEdep->SetPointError(ip, 0., 300);
541  }
542  if (CWRITE()) fgEdep->Print();
543 
544  Double_t e(0.), chi(100), xlo(*vx.begin()), xhi(*vx.rbegin());
545  fgPRF->SetParameter(0, max);
546  fgPRF->FixParameter(1, dx / fDigiPar->GetPadSizeX(0));
547  fgPRF->SetParameter(2, 0.65);
548  fgPRF->SetParLimits(2, 0.45, 10.5);
549  fgEdep->Fit(fgPRF, "QBN", "goff", xlo - 0.5, xhi + 0.5);
550  if (!fgPRF->GetNDF()) return NULL;
551  chi = fgPRF->GetChisquare() / fgPRF->GetNDF();
552  e = fgPRF->Integral(xlo - 0.5, xhi + 0.5);
553 
554  // apply MC correction
555  Float_t gain0 = 210.21387; //(XeCO2 @ 1900V)
556  // Float_t grel[3] = {1., 0.98547803, 0.93164071},
557  // goff[3][3] = {
558  // {0.05714, -0.09, -0.09},
559  // {0., -0.14742, -0.14742},
560  // {0., -0.29, -0.393}
561  // };
562  // Int_t ian=0;
563  // if(TMath::Abs(dy)<=0.3) ian=0;
564  // else if(TMath::Abs(dy)<=0.6) ian=1;
565  // else if(TMath::Abs(dy)<=0.9) ian=2;
566  // Int_t isize=0;
567  // if(n0<=3) isize=0;
568  // else if(n0<=4) isize=1;
569  // else isize=2;
570  Float_t gain = gain0; //*grel[ian];
571  e /= gain; // apply effective gain
572  //e+=goff[ian][isize]; // apply missing energy offset
573 
574  TVector3 local_pad, local_pad_err;
575  Int_t srow, sector = fDigiPar->GetSectorRow(row, srow);
576  fDigiPar->GetPadPosition(sector, col, srow, local_pad, local_pad_err);
577  //printf("r[%2d] c[%2d] max[%d] lr[%d] n0[%d] cM[%d] nM[%d] dx[%7.4f] dy[%7.4f] loc[%6.2f %6.2f %6.2f] err[%6.2f %6.2f %6.2f] e[%f] chi[%f]\n", row, col, mtyp, lr, n0, cM, nM, dx, dy, local_pad[0], local_pad[1], local_pad[2], local_pad_err[0], local_pad_err[1], local_pad_err[2], e, chi);
578  Double_t local[3] = {local_pad[0] + dx, local_pad[1] + dy, local_pad[2]},
579  global[3], globalErr[3] = {edx, edy, 0.};
580  LocalToMaster(local, global);
581 
582  // process time profile
583  for (Int_t idx(1); idx <= n0; idx++) {
584  Double_t dtFEE = fgDT[0] * (vs[idx] - fgDT[1]) * (vs[idx] - fgDT[1])
586  if (vxe[idx] > 0) vx[idx] += dy / fDigiPar->GetPadSizeY(0);
587  fgT->SetPoint(idx - 1, vx[idx], vt[idx] - dtFEE);
588  }
589  xc = vx[n0 + 2];
590  for (Int_t ip(n0); ip < fgT->GetN(); ip++) {
591  fgT->SetPoint(ip, xc, 0);
592  xc += 0.5;
593  }
594 
595  Double_t time(-21.),
596  edt(26.33), // should be parametrized as function of da TODO
597  tdrift(100); // should depend on Ua
598  if (n0 > 1 && (fgT->Fit("pol1", "QC", "goff") == 0)) {
599  TF1* f = fgT->GetFunction("pol1");
600  time = f->GetParameter(0) - fgDT[2];
601  if (TMath::IsNaN(time)) time = -21;
602  //dtime += TMath::Abs(f->GetParameter(1)*(vx[n0+1] - vx[1]));
603  }
604 
605  if (CDRAW()) {
606  Double_t rangex(vx[0] - 0.25), rangeX(vx[n0 + 2] + 0.25);
607  cvs->cd(1);
608  hf->Draw("p");
609  hf->GetXaxis()->SetRangeUser(rangex, rangeX);
610  hf->SetTitle(Form(
611  "%d[%d] row[%d] col[%2d] an[%+d] m[%+4.2f] s[%4.2f] E[%7.2f] chi2[%7.2f]",
612  ic,
613  Int_t(vs.size()),
614  row,
615  col,
616  ia,
617  fgPRF->GetParameter(1),
618  fgPRF->GetParameter(2),
619  e,
620  chi));
621  hf->GetXaxis()->SetRangeUser(rangex, rangeX);
622  hf->GetYaxis()->SetRangeUser(-100., 4500);
623  fgEdep->Draw("pl");
624  fgPRF->Draw("same");
625 
626  cvs->cd(2);
627  hf = (TH1*) hf->DrawClone("p");
628  hf->GetXaxis()->SetRangeUser(rangex, rangeX);
629  hf->GetYaxis()->SetRangeUser(-10, 20);
630  //hf->SetTitle(Form("%d row[%d] col[%2d] an[%+d] m[%+4.2f] s[%4.2f] E[%7.2f] chi2[%7.2f]", ic, row, col, ia, fgPRF->GetParameter(1), fgPRF->GetParameter(2), fgPRF->Integral(xlo, xhi), fgPRF->GetChisquare()/fgPRF->GetNDF()));
631  // hf->GetXaxis()->SetRangeUser(xlo-0.25, xhi+0.25);
632  // //hf->GetYaxis()->SetRangeUser(0., 4500);
633  fgT->Draw("pl");
634  cvs->Modified();
635  cvs->Update();
636  cvs->SaveAs(Form("cl_%02d_A%d.gif", ic, ia));
637  }
638 
639  Int_t nofHits = fHits->GetEntriesFast();
640  CbmTrdHit* hit = new ((*fHits)[nofHits])
642  global,
643  globalErr,
644  0., // sxy chi,
645  ic,
646  e, // energy
647  CbmTrdDigi::Clk(CbmTrdDigi::kFASP) * (t0 + time) - tdrift + 30.29,
648  edt);
649  hit->SetClassType();
650  hit->SetMaxType(tM);
651  if (ovf) hit->SetOverFlow();
652  if (CWRITE()) cout << hit->ToString();
653  return hit;
654 }
655 
656 //_______________________________________________________________________________
657 Double_t CbmTrdModuleRecT::GetXoffset(Int_t n0) const {
658  Double_t dx(0.), R(0.);
659  for (Int_t ir(1); ir <= n0; ir++) {
660  if (vxe[ir] > 0) continue; // select rectangular coupling
661  R += vs[ir];
662  dx += vs[ir] * vx[ir];
663  }
664  if (TMath::Abs(R) > 0) return dx / R;
665  LOG(warn) << GetName() << "::GetDx : Unexpected null sum.";
666  return 0.;
667 }
668 
669 //_______________________________________________________________________________
670 Double_t CbmTrdModuleRecT::GetYoffset(Int_t n0) const {
671  Double_t dy(0.), T(0.);
672  for (Int_t it(1); it <= n0; it++) {
673  if (vxe[it] > 0) { // select tilted coupling
674  T += vs[it];
675  dy += vs[it] * vx[it];
676  }
677  }
678  if (TMath::Abs(T) > 0) return dy / T;
679  LOG(warn) << GetName() << "::GetDy : Unexpected null sum.";
680  return 0.;
681 }
682 
683 //_______________________________________________________________________________
684 Int_t CbmTrdModuleRecT::LoadDigis(vector<const CbmTrdDigi*>* digis,
685  vector<CbmTrdDigi*>* vdgM,
686  vector<Bool_t>* vmask,
687  ULong64_t& t0,
688  Int_t& cM) {
689  /* Load digis information in working vectors.
690  * The digis as provided by the digitizer are replaced by the merged one
691  * according to the vmask map. Digis are represented in the normal coordinate system of
692  * (pad width [pw], DAQ time [clk], signal [ADC chs]) with rectangular signals at integer
693  * positions.
694  */
695  vs.clear();
696  vse.clear();
697  vx.clear();
698  vxe.clear();
699  vt.clear();
700 
701  cM = 0; // relative position of maximum signal
702  Int_t n0(0), // number of measured signals
703  ovf(1), // over-flow flag for at least one of the digis
704  dt;
705  Char_t ddt, // signal time offset wrt prompt
706  dt0(0); // cluster time offset wrt arbitrary t0
707  Double_t r, re(100.), // rect signal
708  t, te(100.), // tilt signal
709  err, // final error parametrization for signal
710  xc(-0.5), // running signal-pad position
711  max(0.); // max signal
712  Int_t j(0), col0(-1), col1(0);
713  const CbmTrdDigi* dg(NULL);
714  vector<CbmTrdDigi*>::iterator idgM = vdgM->begin();
715  for (vector<const CbmTrdDigi*>::iterator i = digis->begin();
716  i != digis->end();
717  i++, j++) {
718  dg = (*i);
719  if ((*vmask)[j]) {
720  dg = (*idgM);
721  idgM++;
722  }
723  if (CWRITE()) cout << dg->ToString();
724  r = dg->GetCharge(t, dt);
725  if (t0 == 0)
726  t0 = dg->GetTimeDAQ(); // set arbitrary t0 to avoid double digis loop
727  if (col0 < 0) GetPadRowCol(dg->GetAddressChannel(), col0); // initialilze
728  ddt = dg->GetTimeDAQ() - t0;
729  if (ddt < dt0) dt0 = ddt;
730 
731  // check column wise organization
732  GetPadRowCol(dg->GetAddressChannel(), col1);
733  if (col0 + j != col1) {
734  LOG(error) << GetName()
735  << "::LoadDigis : digis in cluster not in increasing order !";
736  return 0;
737  }
738 
739  // process tilt signal
740  if (t > 0) {
741  err = te;
742  if (t > 4094) {
743  ovf = -1;
744  err = 150.;
745  }
747  n0++;
748  if (t > max) {
749  max = t;
750  cM = j;
751  }
752  } else
753  err = 300.;
754  vt.push_back(ddt);
755  vs.push_back(t);
756  vse.push_back(err);
757  vx.push_back(xc);
758  vxe.push_back(0.035);
759  xc += 0.5;
760 
761  // process rect signal
762  ddt += dt;
763  if (ddt < dt0) dt0 = ddt;
764  if (r > 0) {
765  err = re;
766  if (r > 4094) {
767  ovf = -1;
768  err = 150.;
769  }
771  n0++;
772  if (r > max) {
773  max = r;
774  cM = j;
775  }
776  } else
777  err = 300.;
778  vt.push_back(ddt);
779  vs.push_back(r);
780  vse.push_back(err);
781  vx.push_back(xc);
782  vxe.push_back(0.);
783  xc += 0.5;
784  }
785 
786  // remove merged digis if they were created
787  if (idgM != vdgM->end())
788  LOG(warn) << GetName()
789  << "::LoadDigis : not all merged digis have been consumed !";
790  for (idgM = vdgM->begin(); idgM != vdgM->end(); idgM++)
791  delete (*idgM);
792 
793  // add front and back anchor points if needed
794  if (TMath::Abs(vs[0]) > 1.e-3) {
795  xc = vx[0];
796  ddt = vt[0];
797  vs.insert(vs.begin(), 0);
798  vse.insert(vse.begin(), 300);
799  vt.insert(vt.begin(), ddt);
800  vx.insert(vx.begin(), xc - 0.5);
801  vxe.insert(vxe.begin(), 0);
802  }
803  Int_t n(vs.size() - 1);
804  if (TMath::Abs(vs[n]) > 1.e-3) {
805  xc = vx[n] + 0.5;
806  ddt = vt[n];
807  vs.push_back(0);
808  vse.push_back(300);
809  vt.push_back(ddt);
810  vx.push_back(xc);
811  vxe.push_back(0.035);
812  }
813  // recenter time and space profile
814  t0 += dt0;
815  for (UInt_t idx(0); idx < vx.size(); idx++) {
816  vt[idx] -= dt0;
817  vx[idx] -= cM;
818  }
819  return ovf * n0;
820 }
821 
822 //_______________________________________________________________________________
823 Bool_t CbmTrdModuleRecT::MergeDigis(vector<const CbmTrdDigi*>* digis,
824  vector<CbmTrdDigi*>* vdgM,
825  vector<Bool_t>* vmask) {
826  /* Merge digis in the cluster if their topology within it allows it although cuts in the
827  * digi merger procedure (CbmTrdFASP::WriteDigi()) were not fulfilled.
828  * Normally this are boundary signals with large time delays wrt neighbors
829  */
830 
831  CbmTrdDigi* dgM(NULL);
832  if (digis->size() < 2) { // sanity check ... just in case
833  LOG(warn) << GetName() << "::MergeDigis : Bad digi config for cluster :";
834  return kFALSE;
835  }
836 
837  Double_t r, t;
838  Int_t colR, colT, dt, contor(0);
839  Bool_t kFOUND(0);
840  for (vector<const CbmTrdDigi*>::iterator idig = digis->begin(),
841  jdig = idig + 1;
842  jdig != digis->end();
843  idig++, jdig++, contor++) {
844  const CbmTrdDigi* dgT((*idig)); // tilt signal digi
845  const CbmTrdDigi* dgR((*jdig)); // rect signal digi
846  GetPadRowCol(dgR->GetAddressChannel(), colR);
847  GetPadRowCol(dgT->GetAddressChannel(), colT);
848 
849  if (colR != colT) continue;
850 
851  dgM = new CbmTrdDigi(*dgT);
852  r = dgR->GetCharge(t, dt);
853  dgT->GetCharge(t, dt);
854  dt = dgR->GetTimeDAQ() - dgT->GetTimeDAQ();
855  dgM->SetCharge(t, r, dt);
856  Int_t rtrg(dgR->GetTriggerType() & 2), ttrg(dgT->GetTriggerType() & 1);
857  dgM->SetTriggerType(rtrg | ttrg); //merge the triggers
858  if (CWRITE()) {
859  cout << "MERGE" << endl;
860  cout << dgT->ToString();
861  cout << dgR->ToString();
862  cout << "TO" << endl;
863  cout << dgM->ToString();
864  cout << "..." << endl;
865  }
866  kFOUND = 1;
867 
868  vdgM->push_back(dgM);
869  (*vmask)[contor] = 1;
870  jdig = digis->erase(jdig);
871  if (jdig == digis->end()) break;
872  }
873 
874  if (!kFOUND) {
875  LOG(warn) << GetName()
876  << "::MergeDigis : Digi to pads matching failed for cluster :";
877  return kFALSE;
878  }
879  return kTRUE;
880 }
881 
882 Float_t CbmTrdModuleRecT::fgCorrXdx = 0.005;
884  {0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
885  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
886  0.000, 0.000, 0.000, -0.144, -0.091, -0.134, -0.185, -0.120, -0.115,
887  -0.125, -0.125, -0.124, -0.124, -0.122, -0.120, -0.119, -0.116, -0.114,
888  -0.113, -0.111, -0.109, -0.107, -0.105, -0.102, -0.100, -0.098, -0.097,
889  -0.093, -0.091, -0.089, -0.087, -0.084, -0.082, -0.079, -0.077, -0.074,
890  -0.072, -0.068, -0.065, -0.062, -0.059, -0.056, -0.053, -0.049, -0.046,
891  -0.043, -0.039, -0.036, -0.032, -0.029, -0.025, -0.022, -0.018, -0.015,
892  -0.011, -0.007, -0.003, 0.000, 0.003, 0.007, 0.011, 0.014, 0.018,
893  0.022, 0.025, 0.029, 0.032, 0.036, 0.039, 0.043, 0.046, 0.049,
894  0.053, 0.056, 0.059, 0.062, 0.065, 0.068, 0.071, 0.074, 0.077,
895  0.080, 0.082, 0.084, 0.087, 0.090, 0.091, 0.094, 0.096, 0.098,
896  0.100, 0.103, 0.104, 0.107, 0.108, 0.110, 0.113, 0.115, 0.116,
897  0.120, 0.121, 0.121, 0.123, 0.125, 0.124, 0.127, 0.140, 0.119,
898  0.114, 0.028, 0.049, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
899  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
900  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000},
901  {0.003, 0.013, 0.026, 0.039, 0.052, 0.065, 0.078, 0.091, 0.104,
902  0.118, 0.132, 0.145, 0.160, 0.174, 0.189, 0.203, 0.219, 0.234,
903  0.250, 0.267, 0.283, 0.301, 0.319, 0.338, 0.357, 0.375, 0.398,
904  0.419, 0.440, 0.464, 0.491, 0.514, 0.541, 0.569, 0.598, 0.623,
905  0.660, 0.696, 0.728, 0.763, 0.804, 0.847, 0.888, 0.930, 0.988,
906  1.015, 1.076, 1.128, 1.167, 1.228, 1.297, 1.374, 1.443, 1.511,
907  1.564, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
908  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
909  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
910  0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
911  0.000, -1.992, -1.884, -1.765, -1.703, -1.609, -1.572, -1.493, -1.426,
912  -1.356, -1.309, -1.243, -1.202, -1.109, -1.069, -1.026, -0.970, -0.933,
913  -0.881, -0.844, -0.803, -0.767, -0.721, -0.691, -0.659, -0.629, -0.596,
914  -0.569, -0.541, -0.514, -0.490, -0.466, -0.441, -0.419, -0.397, -0.377,
915  -0.357, -0.337, -0.319, -0.301, -0.283, -0.267, -0.250, -0.234, -0.218,
916  -0.203, -0.189, -0.174, -0.160, -0.145, -0.131, -0.119, -0.104, -0.091,
917  -0.078, -0.064, -0.052, -0.039, -0.026, -0.013, -0.002}};
918 Float_t CbmTrdModuleRecT::fgCorrYval[NBINSCORRY][2] = {{2.421729, 0.},
919  {1.537359, 0.483472},
920  {1.1752, 0.},
921  {1.154183, -0.090229}};
922 
CbmTrdModuleRecT::vs
std::vector< Double_t > vs
Definition: CbmTrdModuleRecT.h:87
CbmTrdParModDigi::GetPadSizeY
Double_t GetPadSizeY(Int_t i) const
Definition: CbmTrdParModDigi.h:37
CbmTrdModuleRecT::fgDT
static Double_t fgDT[3]
discretized correction params
Definition: CbmTrdModuleRecT.h:97
CbmTrdModuleRecT::fgPRF
static TF1 * fgPRF
data handler for cluster PRF
Definition: CbmTrdModuleRecT.h:99
CbmTrdDigi::GetTimeDAQ
ULong64_t GetTimeDAQ() const
Getter for global DAQ time [clk]. Differs for each ASIC. In FASP case DAQ time is already stored in f...
Definition: CbmTrdDigi.h:131
CbmTrdModuleRecT::MakeHits
virtual Bool_t MakeHits()
Steering routine for building hits.
Definition: CbmTrdModuleRecT.cxx:205
f
float f
Definition: L1/vectors/P4_F32vec4.h:24
NBINSCORRX
#define NBINSCORRX
Definition: CbmTrdModuleRecT.h:8
NANODE
#define NANODE
Definition: CbmTrdModuleRecT.cxx:232
CbmTrdModuleAbstract::fModAddress
UShort_t fModAddress
unique identifier for current module
Definition: CbmTrdModuleAbstract.h:78
CbmTrdModuleRecT::~CbmTrdModuleRecT
virtual ~CbmTrdModuleRecT()
Definition: CbmTrdModuleRecT.cxx:56
CbmTrdFASP.h
CbmTrdModuleRecT::vxe
std::vector< Double_t > vxe
working copy of signal relative positions
Definition: CbmTrdModuleRecT.h:92
CbmTrdModuleRecT::Finalize
virtual Bool_t Finalize()
Finalize hits (merge RC hits, etc)
Definition: CbmTrdModuleRecT.cxx:208
CbmTrdModuleRecT::AddDigi
virtual Bool_t AddDigi(const CbmTrdDigi *d, Int_t id)
Add digi to local module.
Definition: CbmTrdModuleRecT.cxx:59
CbmTrdDigi::SetCharge
void SetCharge(Float_t c)
Charge setter for SPADIC ASIC.
Definition: CbmTrdDigi.cxx:206
CbmTrdModuleRecT::FindClusters
virtual Int_t FindClusters()
Finalize clusters.
Definition: CbmTrdModuleRecT.cxx:157
CbmTrdModuleRecT::fgCorrXval
static Float_t fgCorrXval[2][NBINSCORRX]
step of the discretized correction LUT
Definition: CbmTrdModuleRecT.h:95
CbmTrdDigi::GetAddressChannel
Int_t GetAddressChannel() const
Getter read-out id.
Definition: CbmTrdDigi.cxx:119
CbmTrdModuleRecT.h
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmTrdHit
data class for a reconstructed Energy-4D measurement in the TRD
Definition: CbmTrdHit.h:35
CbmTrdModuleRecT::CWRITE
Bool_t CWRITE() const
Definition: CbmTrdModuleRecT.h:57
CbmTrdHit::SetOverFlow
void SetOverFlow(Bool_t set=kTRUE)
Mark overflow in one or more digits which define the hit.
Definition: CbmTrdHit.h:88
CbmTrdModuleRecT::GetYoffset
Double_t GetYoffset(Int_t n0) const
Definition: CbmTrdModuleRecT.cxx:670
CbmTrdModuleRec
Abstract class for module wise cluster finding and hit reconstruction.
Definition: CbmTrdModuleRec.h:16
CbmTrdDigi::kFASP
@ kFASP
Definition: CbmTrdDigi.h:16
CbmTrdModuleRecT::MergeDigis
Bool_t MergeDigis(std::vector< const CbmTrdDigi * > *digis, std::vector< CbmTrdDigi * > *vdgM, std::vector< Bool_t > *vmask)
Merge R/T signals to digis if topological conditions in cluster are fulfilled.
Definition: CbmTrdModuleRecT.cxx:823
CbmTrdCluster::ToString
virtual std::string ToString() const
Extended functionality.
Definition: CbmTrdCluster.cxx:201
CbmTrdModuleRecT::vt
std::vector< Char_t > vt
working copy of signal errors from cluster
Definition: CbmTrdModuleRecT.h:89
CbmTrdModuleRecT::fgCorrYval
static Float_t fgCorrYval[NBINSCORRY][2]
discretized correction LUT
Definition: CbmTrdModuleRecT.h:96
CbmTrdCluster
Data Container for TRD clusters.
Definition: CbmTrdCluster.h:23
CbmTrdParModDigi::GetPadSizeX
Double_t GetPadSizeX(Int_t i) const
Definition: CbmTrdParModDigi.h:36
CbmTrdDigi.h
CbmTrdDigi::Clk
static Float_t Clk(CbmTrdAsicType ty)
DAQ clock accessor for each ASIC.
Definition: CbmTrdDigi.h:87
CbmTrdModuleRecT::fgCorrXdx
static Float_t fgCorrXdx
working copy of signal relative position errors
Definition: CbmTrdModuleRecT.h:94
CbmTrdCluster::SetTrianglePads
void SetTrianglePads(Bool_t set=kTRUE)
Definition: CbmTrdCluster.h:96
d
double d
Definition: P4_F64vec2.h:24
CbmTrdCluster::GetStartCh
UShort_t GetStartCh() const
Definition: CbmTrdCluster.h:68
CbmTrdModuleRecT::vx
std::vector< Double_t > vx
working copy of signal relative timing
Definition: CbmTrdModuleRecT.h:90
CbmTrdParModDigi.h
CbmTrdParModDigi::GetPadPosition
void GetPadPosition(const Int_t sector, const Int_t col, const Int_t row, TVector3 &padPos, TVector3 &padPosErr) const
Definition: CbmTrdParModDigi.cxx:718
CbmTrdModuleRecT::vse
std::vector< Double_t > vse
working copy of signals from cluster
Definition: CbmTrdModuleRecT.h:88
CbmTrdHit.h
Class for hits in TRD detector.
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmTrdCluster::SetProfileStart
void SetProfileStart(Bool_t set=kTRUE)
Definition: CbmTrdCluster.h:99
CbmTrdModuleRecT::CDRAW
Bool_t CDRAW() const
Definition: CbmTrdModuleRecT.h:56
CbmTrdModuleRecT::GetOverThreshold
virtual Int_t GetOverThreshold() const
Count RO channels (R or T) with data.
Definition: CbmTrdModuleRecT.cxx:142
CbmTrdModuleAbstract::fDigiPar
const CbmTrdParModDigi * fDigiPar
read-out description of module
Definition: CbmTrdModuleAbstract.h:83
CbmTrdCluster::SetProfileStop
void SetProfileStop(Bool_t set=kTRUE)
Definition: CbmTrdCluster.h:102
CbmTrdModuleRecT::MakeHit
virtual CbmTrdHit * MakeHit(Int_t cId, const CbmTrdCluster *c, std::vector< const CbmTrdDigi * > *digis)
Steering routine for converting cluster to hit.
Definition: CbmTrdModuleRecT.cxx:234
CbmTrdModuleRecT
Triangular pad module; Cluster finding and hit reconstruction algorithms.
Definition: CbmTrdModuleRecT.h:16
CbmTrdCluster::GetNCols
UShort_t GetNCols() const
Definition: CbmTrdCluster.h:64
CbmTrdModuleAbstract::GetPadRowCol
virtual Int_t GetPadRowCol(Int_t address, Int_t &c) const
Addressing read-out pads based on module address.
Definition: CbmTrdModuleAbstract.h:100
CbmTrdHit::SetMaxType
void SetMaxType(Bool_t set=kTRUE)
Extra bool definition for the hit (e.g. the type of maximum for triangular pads).
Definition: CbmTrdHit.h:100
CbmTrdDigi::SetTriggerType
void SetTriggerType(const Int_t ttype)
Set digi trigger type (SPADIC only)
Definition: CbmTrdDigi.cxx:237
CbmTrdCluster::HasOpenStop
Bool_t HasOpenStop() const
Definition: CbmTrdCluster.h:72
CbmTrdCluster::GetStartTime
UShort_t GetStartTime() const
Definition: CbmTrdCluster.h:69
CbmTrdHit::ToString
virtual std::string ToString() const
Inherited from CbmBaseHit.
Definition: CbmTrdHit.cxx:39
CbmTrdModuleRecT::LoadDigis
Int_t LoadDigis(std::vector< const CbmTrdDigi * > *digis, std::vector< CbmTrdDigi * > *vdgM, std::vector< Bool_t > *vmask, ULong64_t &t0, Int_t &cM)
Load digis info into local data structures.
Definition: CbmTrdModuleRecT.cxx:684
CbmTrdDigi::GetTriggerType
Int_t GetTriggerType() const
Channel trigger type. SPADIC specific see CbmTrdTriggerType.
Definition: CbmTrdDigi.h:135
CbmTrdCluster::HasOpenStart
Bool_t HasOpenStart() const
Definition: CbmTrdCluster.h:71
CbmTrdModuleRecT::fgEdep
static TGraphErrors * fgEdep
FASP delay wrt signal.
Definition: CbmTrdModuleRecT.h:98
CbmTrdModuleAbstract::LocalToMaster
virtual void LocalToMaster(Double_t in[3], Double_t out[3])
Definition: CbmTrdModuleAbstract.cxx:33
CbmTrdHit::SetClassType
void SetClassType(Bool_t set=kTRUE)
Type of pad layout used in reconstruction R[0], T[1].
Definition: CbmTrdHit.h:96
CbmTrdDigi
Definition: CbmTrdDigi.h:14
CbmTrdCluster.h
Data Container for TRD clusters.
NBINSCORRY
#define NBINSCORRY
Definition: CbmTrdModuleRecT.h:9
CbmTrdFASP::GetBaselineCorr
static Float_t GetBaselineCorr()
Return the baseline value in ADC ch.
Definition: CbmTrdFASP.h:38
max
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: L1/vectors/P4_F32vec4.h:36
CbmTrdModuleRecT::fgT
static TGraphErrors * fgT
fitter for cluster PRF
Definition: CbmTrdModuleRecT.h:100
CbmTrdDigi::GetCharge
Double_t GetCharge() const
Charge getter for SPADIC.
Definition: CbmTrdDigi.cxx:133
CbmTrdModuleRec::fHits
TClonesArray * fHits
module wise storage of reconstructed hits
Definition: CbmTrdModuleRec.h:67
CbmTrdModuleRecT::GetXoffset
Double_t GetXoffset(Int_t n0) const
Definition: CbmTrdModuleRecT.cxx:657
CbmTrdParModDigi::GetSectorRow
Int_t GetSectorRow(Int_t growId, Int_t &srowId) const
Find the sector wise row given the module row. Inverse of GetModuleRow()
Definition: CbmTrdParModDigi.cxx:411
CbmTrdModuleRecT::fT0
ULong64_t fT0
Definition: CbmTrdModuleRecT.h:84
CbmTrdModuleRecT::CbmTrdModuleRecT
CbmTrdModuleRecT()
Default constructor.
Definition: CbmTrdModuleRecT.cxx:28
CbmTrdModuleRecT::fBuffer
std::map< Int_t, std::list< CbmTrdCluster * > > fBuffer
start time of event/time slice [clk]
Definition: CbmTrdModuleRecT.h:86
CbmTrdDigi::ToString
std::string ToString() const
String representation of a TRD digi. Account for digi type and specific information.
Definition: CbmTrdDigi.cxx:243