]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDtrackletOflHelper.cxx
Removed useless function from AliMixInputEventHandler
[u/mrichter/AliRoot.git] / TRD / AliTRDtrackletOflHelper.cxx
CommitLineData
212b7ec8 1#include "TObjArray.h"
2#include "TMath.h"
3
4#include "AliLog.h"
5#include "AliMathBase.h"
6
7#include "AliTRDseedV1.h"
8#include "AliTRDcluster.h"
9#include "AliTRDpadPlane.h"
10#include "AliTRDtrackletOflHelper.h"
11
12ClassImp(AliTRDtrackletOflHelper)
13//___________________________________________________________________
14AliTRDtrackletOflHelper::AliTRDtrackletOflHelper()
15 :TObject()
16 ,fRow(-1)
17 ,fClusters(NULL)
18 ,fPadPlane(NULL)
19{
20// Default constructor
21 fCol[0]=144; fCol[1]=0;
22 fTBrange[0]=1; fTBrange[1]=21;
23}
24
25//___________________________________________________________________
26AliTRDtrackletOflHelper::AliTRDtrackletOflHelper(const AliTRDtrackletOflHelper &ref)
27 :TObject(ref)
28 ,fRow(-1)
29 ,fClusters(NULL)
30 ,fPadPlane(ref.fPadPlane)
31{
32// Copy constructor
33 fRow = ref.fRow; fCol[0]=ref.fCol[0]; fCol[1]=ref.fCol[1];
34 fTBrange[0] = ref.fTBrange[0];fTBrange[1] = ref.fTBrange[1];
35 Int_t n(0);
36 if(ref.fClusters){
37 if((n = ref.fClusters->GetEntriesFast())) {
38 fClusters = new TObjArray(n);
39 for(Int_t ic(n); ic--;) fClusters->AddAt(ref.fClusters->At(ic), ic);
40 }
41 }
42}
43
44//___________________________________________________________________
45AliTRDtrackletOflHelper& AliTRDtrackletOflHelper::operator=(const AliTRDtrackletOflHelper &rhs)
46{
47 if(this != &rhs){
48 TObject::operator=(rhs);
49 fPadPlane = rhs.fPadPlane;
50 fRow = rhs.fRow;
51 fCol[0] = rhs.fCol[0];
52 fCol[1] = rhs.fCol[1];
53 fTBrange[0]= rhs.fTBrange[0];
54 fTBrange[1]= rhs.fTBrange[1];
55 if(rhs.fClusters){
56 Int_t n(rhs.fClusters->GetEntriesFast());
57 if(!fClusters) fClusters = new TObjArray(n);
58 if(fClusters->GetEntriesFast() != n){
59 fClusters->Clear();
60 fClusters->Expand(n);
61 }
62 for(Int_t ic(n); ic--;) fClusters->AddAt(rhs.fClusters->At(ic), ic);
63 } else {
64 if(fClusters) delete fClusters;
65 }
66 }
67 return *this;
68}
69
70//___________________________________________________________________
71AliTRDtrackletOflHelper::~AliTRDtrackletOflHelper()
72{
73// Clean helper
74 if(fClusters) fClusters->Clear();
75 delete fClusters;
76}
77
78//___________________________________________________________________
79Int_t AliTRDtrackletOflHelper::Expand(TObjArray *cls, Int_t *mark, Int_t groupId)
80{
81// Allocate new clusters for this helper. If a subset of clusters is to be allocated
82// this can be specified via the identifier "groupId" which have to be compared
83// with individual elements in the array "mark".
84//
85// Return total no. of clusters in helper
86
87 if(!fPadPlane || !fClusters ){
88 AliError("Helper not initialized.");
89 return 0;
90 }
91 Int_t ncl(fClusters->GetEntriesFast());
92 Int_t mcl(cls->GetEntriesFast());
93 fClusters->Expand(mcl + ncl);
94 fCol[0]=144; fCol[1]=0;
95 AliTRDcluster *c(NULL);
96 for(Int_t icl(0), jcl(ncl); icl<mcl; icl++){
97 if(!(c=(AliTRDcluster*)cls->At(icl))) continue;
98 if(mark && mark[icl]!=groupId) continue;
99 fClusters->AddAt(c, jcl++);
100 Int_t col(c->GetPadCol());
101 Short_t *sig(c->GetSignals());
102 for(Int_t icol(0), jcol(col-3); icol<7; icol++, jcol++){
103 if(sig[icol]==0) continue;
104 if(jcol<fCol[0]) fCol[0]=jcol;
105 if(jcol>fCol[1]) fCol[1]=jcol;
106 }
107 }
108
109 AliDebug(1, Form("Segment[%d] Clusters[%2d] pad[%3d %3d].", groupId, fClusters->GetEntriesFast(), fCol[0], fCol[1]));
110 return fClusters->GetEntriesFast();
111}
112
113
114//___________________________________________________________________
115Int_t AliTRDtrackletOflHelper::Init(AliTRDpadPlane *p, TObjArray *cls, Int_t *mark, Int_t groupId)
116{
117// Allocate clusters for this helper. If a subset of clusters is to be allocated
118// this can be specified via the identifier "groupId" which have to be compared
119// with individual elements in the array "mark".
120//
121// Return no. of clusters allocated
122
123 if(!p){
124 AliError("PadPlane not initialized."); return 0;
125 }
126 if(!cls){
127 AliError("Cluster array not initialized."); return 0;
128 }
129 Int_t ncl(cls->GetEntriesFast());
130 if(ncl<2){
131 AliDebug(1, Form("Segment[%d] failed n[%d].", groupId, ncl)); return 0;
132 }
133 Int_t mcl(ncl);
134 if(mark){
135 mcl = 0;
136 for(Int_t icl(ncl); icl--;)
137 if(cls->At(icl) && mark[icl]==groupId) mcl++;
138 }
139 if(mcl<2){
140 AliDebug(1, Form("Segment[%d] failed n[%d] in group.", groupId, mcl)); return 0;
141 }
142 if(!fClusters) fClusters = new TObjArray(mcl);
143 else{
144 fClusters->Clear();
145 fClusters->Expand(mcl);
146 }
147 fCol[0]=144; fCol[1]=0; fRow=-1;
148 AliTRDcluster *c(NULL);
149 for(Int_t icl(0), jcl(0); icl<ncl; icl++){
150 if(!(c=(AliTRDcluster*)cls->At(icl))) continue;
151 if(mark && mark[icl]!=groupId) continue;
152 fClusters->AddAt(c, jcl++);
153 if(fRow<0) fRow = c->GetPadRow();
154 Int_t col(c->GetPadCol());
155 Short_t *sig(c->GetSignals());
156 for(Int_t icol(0), jcol(col-3); icol<7; icol++, jcol++){
157 if(sig[icol]==0) continue;
158 if(jcol<fCol[0]) fCol[0]=jcol;
159 if(jcol>fCol[1]) fCol[1]=jcol;
160 }
161 }
162 fPadPlane = p;
163
164 AliDebug(1, Form("Segment[%d] Clusters[%2d] pad[%3d %3d].", groupId, fClusters->GetEntriesFast(), fCol[0], fCol[1]));
165 return fClusters->GetEntriesFast();
166}
167
168//___________________________________________________________________
169Int_t AliTRDtrackletOflHelper::ClassifyTopology()
170{
171// Classify topology and return classification code
172// 0 - normal tracklet
173// 1 - delta ray candidate
174// 2 - secondary candidate
175// 3 - "elephant" candidate
176// 4 - unknown topology
177
178 if(!fClusters){
179 AliError("Helper not initialized. Missing clusters.");
180 return kUnknown;
181 }
182 Int_t ncl(fClusters->GetEntries());
183 if(!ncl){
184 AliError("Helper not initialized. No cluster allocated.");
185 return kUnknown;
186 }
187 const Int_t kRange(22); // DEFINE based on vd
188
189 // compute local occupancy
190 Int_t localOccupancy[AliTRDseedV1::kNtb]; memset(localOccupancy, 0, AliTRDseedV1::kNtb*sizeof(Int_t));
191 AliTRDcluster *c(NULL);
192 Double_t sy[kNcls]; Int_t mcl(0);
193 for(Int_t icl(ncl), time; icl--;){
194 c = (AliTRDcluster*)fClusters->At(icl);
195 time = c->GetPadTime();
196 if(time==0 || time>=kRange){
197 sy[icl] = -1; // mark clusters outsde drift volume
198 continue;
199 }
200 // protect against wrong error param.
201 if(c->GetSigmaY2() < 1.e-5) sy[icl] = 0.02;
202 else sy[icl] = TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.04);
203 localOccupancy[time]++;
204 mcl++;
205 }
206 if(!mcl){
207 AliWarning("No clusters in the active area.");
208 return kUnknown;
209 }
210 //compute number of 0 bins and high occupancy
211 Int_t goodOccupancy(0), highOccupancy(0), lowOccupancy(0);
212 for(Int_t itb(1); itb<kRange; itb++){
213 switch(localOccupancy[itb]){
214 case 0: lowOccupancy++; break;
215 case 1: goodOccupancy++; break;
216 default: highOccupancy++; break;
217 }
218 }
219 AliDebug(2, Form("H[%2d | %5.2f%%] L[%2d | %5.2f%%] N[%2d | %5.2f%% | %5.2f%%]",
220 highOccupancy, 100.*highOccupancy/kRange,
221 lowOccupancy, 100.*lowOccupancy/kRange,
222 goodOccupancy, 100.*goodOccupancy/kRange, 100.*goodOccupancy/mcl));
223
224 // filter
225 Double_t dy[kNcls];
226 if(goodOccupancy==mcl) return kNormal;
227 else if(Double_t(goodOccupancy)/kRange > 0.9){
228 if(highOccupancy == 0 ) {
229 if(!FitPSR(dy)) return kUnknown;
230 Int_t nin(0);
231 for(Int_t idy(ncl); idy--;){
232 if(sy[idy]<0.) continue;
233 if(dy[idy] > 3*sy[idy]) continue;
234 nin ++;
235 }
236 if(Double_t(nin)/mcl > 0.8) return kNormal;
237 else return kUnknown;
238 } else return kDeltaRay;
239 } else return kUnknown;
240}
241
242
243//___________________________________________________________________
244void AliTRDtrackletOflHelper::FindSolidCls(Bool_t *mark, Int_t *q)
245{
246// Find clusters produced by large fluctuations of energy deposits
247// Largest charge and well separation from neighbors
248
249 Int_t idx[AliTRDseedV1::kNtb];
250 TMath::Sort(Int_t(AliTRDseedV1::kNtb), q, idx, kTRUE);
251 Int_t qmax = Int_t(0.3*q[idx[0]]);
252 mark[0] = kFALSE;
253 for(Int_t icl(AliTRDseedV1::kNtb-5); icl<AliTRDseedV1::kNtb; icl++) mark[icl] = kFALSE;
254 for(Int_t icl(0); icl<AliTRDseedV1::kNtb; icl++){
255 Int_t jcl(idx[icl]);
256 if(!mark[jcl]) continue;
257 if(q[jcl-1]>q[jcl] || q[jcl+1]>q[jcl]){
258 mark[jcl] = kFALSE;
259 continue;
260 }
261 if(q[jcl] < qmax){
262 mark[jcl] = kFALSE;
263 continue;
264 }
265 for(Int_t kcl=TMath::Max(0, jcl-2); kcl<jcl+3; kcl++){
266 if(kcl==jcl) continue;
267 mark[kcl] = kFALSE;
268 }
269 }
270}
271
272//___________________________________________________________________
273Bool_t AliTRDtrackletOflHelper::FitPSR(Double_t dy[200], Bool_t useSolid)
274{
275// Fit tracklet in Pad System of Reference [PSR] to avoid uncertainty related to
276// Lorentz angle correction
277
278 Bool_t mark[200]; memset(mark, 0, 200*sizeof(Bool_t));
279 Int_t q[200]; memset(q, 0, 200*sizeof(Int_t));
280 Int_t ncl(fClusters->GetEntries());
281 AliTRDcluster *c(NULL);
282 for(Int_t icl(ncl); icl--;){
283 c = (AliTRDcluster*)fClusters->At(icl);
284 if(c->GetPadRow() != fRow) continue;
285 Int_t time(c->GetPadTime());
286 mark[time] = kTRUE; q[time] = Int_t(c->GetQ());
287 }
288 if(useSolid) FindSolidCls(mark, q);
289
290 Double_t
291 x[AliTRDseedV1::kNtb], y[AliTRDseedV1::kNtb], sy[AliTRDseedV1::kNtb],
292 xf[AliTRDseedV1::kNtb], yf[AliTRDseedV1::kNtb], syf[AliTRDseedV1::kNtb],
293 par[3];
294 Int_t jcl(0), kcl(0);
295 for(Int_t icl(0); icl<ncl; icl++){
296 c = (AliTRDcluster*)fClusters->At(icl);
297 if(c->GetPadRow() != fRow) continue;
298 Int_t col(c->GetPadCol());
299 Int_t time(c->GetPadTime());
300 Double_t center(c->GetCenter());
301 Double_t cw(fPadPlane->GetColSize(col));
302 //Double_t corr = AliTRDcluster::GetYcorr(AliTRDgeometry::GetLayer(det), center);
303 y[jcl] = fPadPlane->GetColPos(col) + (.5 + center)*cw /*+ corr*/;
304 x[jcl] = c->GetX();
305 sy[jcl]= TMath::Sqrt(c->GetSigmaY2());
306 if(mark[time]){
307 yf[kcl] = y[jcl];
308 xf[kcl] = x[jcl];
309 syf[kcl]= useSolid ? 0.5/time:sy[jcl];
310 kcl++;
311 }
312 jcl++;
313 }
314 Fit(kcl, xf, yf, syf, par);
315
316 for(Int_t icl(0); icl<jcl; icl++){
317 Double_t dx(x[icl] - par[2]);
318 dy[icl] = y[icl] - (par[0] + par[1]*dx);
319 }
320 return kTRUE;
321}
322
323//___________________________________________________________________
324Bool_t AliTRDtrackletOflHelper::Fit(Int_t n, Double_t *x, Double_t *y, Double_t *sy, Double_t *par, Double_t sCut, Double_t *cov)
325{
326// Iterative robust tracklet fit
327 if(n<3) return kFALSE;
328 //compute <x>
329 if(Int_t(par[2])==21122012){
330 par[2] = 0.;
331 } else {
332 par[2] = 0.; for(Int_t ic(n); ic--;) par[2] += x[ic]; par[2] /= n;
333 }
334 AliTRDtrackerV1::AliTRDLeastSquare &f=Fitter();
335 for(Int_t iter(0); iter<3; iter++){
336 f.Reset();
337 Int_t jp(0);
338 for(Int_t ip(0); ip<n; ip++){
339 Double_t dx(x[ip]-par[2]);
340 Double_t dy(y[ip] - (par[0] + par[1]*dx));
341 if(iter && TMath::Abs(dy)>sCut*sy[ip]) continue;
342 f.AddPoint(&dx, y[ip], sy[ip]);
343 jp++;
344 }
345 if(jp<3) continue;
346 if(!f.Eval()) continue;
347 par[0]=f.GetFunctionParameter(0);
348 par[1]=f.GetFunctionParameter(1);
349 if(cov) f.GetCovarianceMatrix(cov);
350 AliDebugGeneral("AliTRDtrackletOflHelper::Fit()", 2, Form("Iter[%d] Ncl[%2d/%2d] par[%f %f %f]", iter, jp, n, par[0], par[1], par[2]));
351 }
352 return kTRUE;
353}
354
355//___________________________________________________________________
356Bool_t AliTRDtrackletOflHelper::Fit(Double_t *par, Double_t sCut) const
357{
358// Wrapper for clusters attach to this for static Fit function
359//
360 Int_t n(fClusters->GetEntriesFast()), jc(0), dr(0);
361 Double_t corr = TMath::Tan(TMath::DegToRad()*fPadPlane->GetTiltingAngle())*
362 fPadPlane->GetLengthIPad();
363 Double_t x[kNcls], y[kNcls], sy[kNcls];
364 AliTRDcluster *c(NULL);
365 for(Int_t ic(0); ic<n; ic++){
366 c = (AliTRDcluster*)fClusters->At(ic);
367 if(!c->IsInChamber()) continue;
368 dr = c->GetPadRow() - fRow;
369 x[jc] = c->GetX();
370 y[jc] = c->GetY()+corr*dr;
371 sy[jc]= c->GetSigmaY2()>0?(TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.08)):0.08;
372 jc++;
373 }
374 return Fit(jc, x, y, sy, par, sCut);
375}
376
377//___________________________________________________________________
378void AliTRDtrackletOflHelper::GetColSignals(Int_t col, Int_t adc[32], Bool_t mainRow) const
379{
380 memset(adc, 0, 32*sizeof(Int_t));
381 if(col<fCol[0] || col>fCol[1]) return;
382
383 AliTRDcluster *cc(NULL);
384 for(Int_t ic(fClusters->GetEntriesFast()); ic--;){
385 cc = (AliTRDcluster*)fClusters->At(ic);
386 if((mainRow && cc->GetPadRow()!=fRow) ||
387 (!mainRow && cc->GetPadRow()==fRow)) continue;
388 Short_t *sig = cc->GetSignals();
389 Int_t padcol(cc->GetPadCol());
390 Int_t time(cc->GetPadTime());
391 for(Int_t icol(0), jcol(padcol-3); icol<7; icol++, jcol++){
392 if(jcol!=col) continue;
393 adc[time]+=sig[icol];
394 }
395 }
396}
397
398//___________________________________________________________________
399Int_t AliTRDtrackletOflHelper::GetRMS(Double_t &r, Double_t &m, Double_t &s, Double_t &xm) const
400{
401// Calculate Rotation[r], Mean y[m] (at mean radial position [xm]) and Sigma[s] (of a gaussian distribution in the tracklet [SR])
402// for clusters attach to this helper. The Rotation and Mean are calculated without tilt correction option.
403// It returns the number of clusters in 1 sigma cut.
404
405 Int_t n(fClusters->GetEntriesFast());
406 if(n==2){
407 return n;
408 }
409
410 Double_t par[3] = {0., 0., 0.};
411 Fit(par);
412 xm= par[2];
413 r = par[1];
414
415 Double_t corr = TMath::Tan(TMath::DegToRad()*fPadPlane->GetTiltingAngle())*
416 fPadPlane->GetLengthIPad();
417 Double_t y[kNcls];
418 AliTRDcluster *c(NULL);
419 for(Int_t ic(n); ic--;){
420 c = (AliTRDcluster*)fClusters->At(ic);
421 Double_t x(c->GetX() - xm);
422 Int_t dr(c->GetPadRow() - fRow);
423 y[ic] = c->GetY()+corr*dr - (par[0] + par[1]*x);
424 }
425 Double_t m1(0.);
426 AliMathBase::EvaluateUni(n, y, m1, s, 0);
427 m = par[0] + m1;
428 Int_t n0(0);
429 for(Int_t ic(n); ic--;){
430 c = (AliTRDcluster*)fClusters->At(ic);
431 Double_t sy = c->GetSigmaY2()>0?(TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.08)):0.08;
432 if(TMath::Abs(y[ic]-m1) <= sy) n0++;
433 }
434 return n0;
435}
436
437//___________________________________________________________________
438Double_t AliTRDtrackletOflHelper::GetSyMean() const
439{
440// Calculate mean uncertainty of clusters
441
442 Double_t sym(0.);
443 Int_t n(fClusters->GetEntriesFast());
444 AliTRDcluster *c(NULL);
445 for(Int_t ic(n); ic--;){
446 c = (AliTRDcluster*)fClusters->At(ic);
447 sym += c->GetSigmaY2()>0?(TMath::Min(TMath::Sqrt(c->GetSigmaY2()), 0.08)):0.08;
448 }
449 return sym/=n;
450}
451
452//___________________________________________________________________
453Int_t AliTRDtrackletOflHelper::Segmentation(Int_t n, Double_t *x, Double_t *y, Int_t *Index)
454{
455// Segmentation of clusters in the tracklet roads
456//
457// The user supply the coordinates of the clusters (x,y) and their number "n". Also
458// The array "Index" has to be allocated by the user with a size equal or larger than "n".
459// On return the function returns the number of segments found and the array "Index" i-th element
460// is filled with the index of the tracklet segment to which the i-th cluster was assigned too.
461//
462// Observation
463// The parameter which controls the segmentation is set inside the function "kGapSize" and it is used
464// for both "x" and "y" segmentations. An improvement can be a parametrization as function of angle of
465// incidence.
466//
467// author:
468// Alex Bercuci <a.bercuci@gsi.de>
469
470 if(!n || !x || !y || !Index){
471 AliErrorGeneral("AliTRDtrackletOflHelper::Segmentation()", "One of the input arrays non initialized.");
472 return 0;
473 }
474
475 const Double_t kGapSize(0.2); // cm
476 Int_t ng(0),
477 nc(0);
478 Double_t xx[100], dy;
479 Int_t idx[200], jdx[100], kdx[100];
480
481 TMath::Sort(n, y, idx);
482 for(Int_t iy(0); iy<n; iy++){
483 dy = iy>0?(TMath::Abs(y[idx[iy-1]]-y[idx[iy]])):0.;
484 if(dy>kGapSize){
485 TMath::Sort(nc, xx, jdx);
486 for(Int_t ic(0), jc0, jc1; ic<nc; ic++){
487 jc0 = ic>0?kdx[jdx[ic-1]]:0;
488 jc1 = kdx[jdx[ic]];
489 dy = TMath::Abs(y[jc0] - y[jc1]);
490 if(ic && dy>kGapSize) ng++;
491 Index[jc1] = ng;
492 AliDebugGeneral("AliTRDtrackletOflHelper::Segmentation()", 4, Form(" y[%2d]=%+f x[%+f] %2d -> %2d -> %2d ng[%2d]", jc1, y[jc1], x[jc1], ic, jdx[ic], jc1, ng));
493 }
494 ng++;
495 nc=0;
496 }
497 xx[nc] = x[idx[iy]];
498 kdx[nc]= idx[iy];
499 AliDebugGeneral("AliTRDtrackletOflHelper::Segmentation()", 4, Form("y[%2d]=%+f -> %2d", idx[iy], y[idx[iy]], nc));
500 nc++;
501 }
502 if(nc){
503 TMath::Sort(nc, xx, jdx);
504 for(Int_t ic(0), jc0, jc1; ic<nc; ic++){
505 jc0 = ic>0?kdx[jdx[ic-1]]:0;
506 jc1 = kdx[jdx[ic]];
507 dy = TMath::Abs(y[jc0] - y[jc1]);
508 if(ic && dy>kGapSize) ng++;
509 Index[jc1] = ng;
510 AliDebugGeneral("AliTRDtrackletOflHelper::Segmentation()", 4, Form(" y[%2d]=%+f x[%+f|%+f] %2d -> %2d -> %2d ng[%2d]\n", jc1, y[jc1], xx[jdx[ic]], x[jc1], ic, jdx[ic], jc1, ng));
511 }
512 ng++;
513 }
514 return ng;
515}
516
517//___________________________________________________________________
518void AliTRDtrackletOflHelper::SetTbRange(Float_t t0, Float_t vd)
519{
520// Set first time bin and total number of time bins corresponding to clusters
521// in chamber based on the calibrated info "t0" and drift velocity "vd"
522
523 // TO CHECK
524 fTBrange[0] = Int_t(t0*0.1);
525 fTBrange[1] = Int_t(vd*0.1);
526}
527
528#include "TH2.h"
529#include "TGraph.h"
530#include "TGraphErrors.h"
531#include "TCanvas.h"
532#include "TLegend.h"
533#include "TGaxis.h"
534//___________________________________________________________________
535void AliTRDtrackletOflHelper::View(TVirtualPad *vpad)
536{
537// Visualization support. Draw this tracklet segment
538
539
540 Int_t row(-1), det(-1);
541 Int_t n(fClusters->GetEntriesFast());
542 AliInfo(Form("Processing Clusters[%2d] Class[%d].", n ,ClassifyTopology()));
543
544 // prepare drawing objects
545 Int_t ncols(fPadPlane->GetNcols());
546 Double_t cw(fPadPlane->GetColSize(1));
547 TH2 *h2 = new TH2I("h2pm", ";y_{Local} [cm]; pad time;Charge", ncols, fPadPlane->GetColPos(1)-cw, fPadPlane->GetColPos(ncols-1)+cw, 31, -30.5, 0.5);
548 h2->SetMarkerColor(kWhite);h2->SetMarkerSize(2.);
549 h2->SetFillColor(9);
550
551 TGraph *gcls = new TGraph(n);
552 gcls->SetMarkerColor(kBlack);gcls->SetMarkerStyle(20);
553 TGraph *gclsLC = new TGraph(n);
554 gclsLC->SetMarkerColor(kBlack);gclsLC->SetMarkerStyle(28);
555 TGraph *gclsSC = new TGraph(n);
556 gclsSC->SetMarkerColor(kBlack);gclsSC->SetMarkerStyle(4);gclsSC->SetMarkerSize(1.5);
557 TGraph *gFP = new TGraph(n);
558 gFP->SetMarkerColor(kBlack);gFP->SetLineWidth(2);
559
560 // fill signal data
561 Bool_t map[200]; memset(map, 0, 200*sizeof(Bool_t));
562 Int_t qa[200]; memset(qa, 0, 200*sizeof(Int_t));
563 AliTRDcluster *cc(NULL); Int_t tm = 30; Double_t ym=0.;
564 for(Int_t ic(n); ic--;){
565 cc = (AliTRDcluster*)fClusters->At(ic);
566 if(cc->GetPadRow() != fRow) continue; //TODO extend for 2 pad rows
567 Short_t *sig = cc->GetSignals();
568 Int_t col(cc->GetPadCol());
569 det = cc->GetDetector();
570 row = cc->GetPadRow();
571 Int_t time(cc->GetPadTime());
572 map[time] = kTRUE; qa[time] = (Int_t)cc->GetQ();
573 for(Int_t ipad(0), jpad(col-2); ipad<7; ipad++, jpad++){
574 Int_t q = (Int_t)h2->GetBinContent(jpad, 31-time);
575 h2->SetBinContent(jpad, 31-time, q+sig[ipad]);
576 }
577 Double_t y0 = fPadPlane->GetColPos(col) + (.5 + cc->GetCenter())*cw;
578 //h2->GetXaxis()->GetBinCenter(pad)+cc->GetCenter();
579 gcls->SetPoint(ic, y0, -time);
580 if(time <= tm) {ym = cc->GetY() - y0; tm = time;}
581 }
582
583 // draw special clusters (solid and Lorentz corrected and fit)
584 FindSolidCls(map, qa);
585 Double_t ddy[200]; FitPSR(ddy, kTRUE);
586 Double_t dt, dy;
587 for(Int_t ic(0), jc(0); ic<n; ic++){
588 cc = (AliTRDcluster*)fClusters->At(ic);
589 if(cc->GetPadRow() != fRow) continue; //TODO extend for 2 pad rows
590 gcls->GetPoint(ic, dy, dt);
591 gclsLC->SetPoint(ic, cc->GetY()-ym, dt);
592 gFP->SetPoint(ic, dy-ddy[ic], dt);
593 if(map[cc->GetPadTime()]) gclsSC->SetPoint(jc++, dy, dt);
594 }
595
596
597 // prepare frame histo
598 h2->SetName(Form("h2s%03d%02d", det, row));
599 Int_t binSrt(fCol[0]+1), binSop(fCol[1]+1);
600 TH1 *h1(NULL);
601
602 // show everything
603 if(!vpad){
604 vpad = new TCanvas(Form("c%03d%02d", det, row), Form("D %03d [%02d_%d_%d] R[%02d]", det, AliTRDgeometry::GetSector(det), AliTRDgeometry::GetStack(det), AliTRDgeometry::GetLayer(det), row), 700, 500);
605 }
606 TVirtualPad *pp(NULL);
607 vpad->Divide(2,1,2.e-5,2.e-5);
608 pp = vpad->cd(1); pp->SetRightMargin(0.0001);pp->SetTopMargin(0.1);pp->SetBorderMode(0); pp->SetFillColor(kWhite);
609 h1 = h2->ProjectionY(); h1->GetYaxis()->SetTitle("Total Charge");
610 h1->SetTitle(Form("D%03d[%02d_%d_%d] R[%02d]",
611 det,AliTRDgeometry::GetSector(det), AliTRDgeometry::GetStack(det), AliTRDgeometry::GetLayer(det), row));
612 TAxis *ax(h1->GetXaxis()); for(Int_t ib(0), jb(1); ib<ax->GetNbins(); ib++, jb++) ax->SetBinLabel(jb, Form("%d", -Int_t(ax->GetBinCenter(jb))));
613 h1->Draw("hbar3");
614
615 pp = vpad->cd(2);
616 pp->SetRightMargin(0.15);pp->SetLeftMargin(0.0001);pp->SetTopMargin(0.1);
617 pp->SetBorderMode(0); pp->SetFillColor(kWhite);
618 ax=h2->GetXaxis();
619 ax->SetRange(TMath::Max(1, binSrt-1), TMath::Min(ncols, binSop+1));
620 h2->Draw("coltextz");
621 gPad->Update();
622 TGaxis *axis = new TGaxis(gPad->GetUxmin(),
623 gPad->GetUymax(),
624 gPad->GetUxmax(),
625 gPad->GetUymax(),
626 TMath::Max(0, binSrt-2)-0.5, TMath::Min(ncols-1, binSop)+0.5, 510,"-L");
627
628 axis->SetNdivisions(103+binSop-binSrt);
629 axis->SetTitle("pad col");
630 axis->Draw();
631
632 TLegend *leg = new TLegend(0.01, 0.1, 0.84, 0.21);
633 leg->SetBorderSize(1); leg->SetFillColor(kYellow-9);leg->SetTextSize(0.04);
634 gcls->Draw("p"); leg->AddEntry(gcls, "Cls. in Pad [SR]", "p");
635 gclsLC->Draw("p"); leg->AddEntry(gclsLC, "Lorentz Corr. Cls.", "p");
636 gclsSC->Draw("p"); leg->AddEntry(gclsSC, "Solid Cls.", "p");
637 gFP->Draw("l"); leg->AddEntry(gFP, "Fit in Pad [SR]", "l");
638 leg->Draw();
639}
640
641//___________________________________________________________________
642AliTRDtrackerV1::AliTRDLeastSquare& AliTRDtrackletOflHelper::Fitter()
643{
644 static AliTRDtrackerV1::AliTRDLeastSquare f;
645 return f;
646}