]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSVertexer3D.cxx
new 3D vertexer based on SPD recppoints
[u/mrichter/AliRoot.git] / ITS / AliITSVertexer3D.cxx
CommitLineData
70c95f95 1/**************************************************************************
2 * Copyright(c) 2006-2008, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15#include <AliESDVertex.h>
16#include <AliITSVertexer3D.h>
17#include <AliStrLine.h>
18#include <AliVertexerTracks.h>
19#include <Riostream.h>
20#include <TH3F.h>
21#include <TTree.h>
22#include<TClonesArray.h>
23#include "AliRunLoader.h"
24#include "AliITSLoader.h"
25#include "AliITSgeom.h"
26#include "AliITSDetTypeRec.h"
27#include "AliITSRecPoint.h"
28/////////////////////////////////////////////////////////////////
29// this class implements a method to determine
30// the 3 coordinates of the primary vertex
31// for p-p collisions
32// It can be used successfully with Pb-Pb collisions
33////////////////////////////////////////////////////////////////
34
35ClassImp(AliITSVertexer3D)
36
37//______________________________________________________________________
38AliITSVertexer3D::AliITSVertexer3D():AliITSVertexer(),
39fLines(),
40fVert3D(),
41fCoarseDiffPhiCut(0.),
42fMaxRCut(0.),
43fZCutDiamond(0.),
44fDCAcut(0.),
45fDiffPhiMax(0.) {
46 // Default constructor
47 SetCoarseDiffPhiCut();
48 SetMaxRCut();
49 SetZCutDiamond();
50 SetDCAcut();
51 SetDiffPhiMax();
52
53}
54
55//______________________________________________________________________
56AliITSVertexer3D::AliITSVertexer3D(TString fn): AliITSVertexer(fn),
57fLines(),
58fVert3D(),
59fCoarseDiffPhiCut(0.),
60fMaxRCut(0.),
61fZCutDiamond(0.),
62fDCAcut(0.),
63fDiffPhiMax(0.) {
64 // Standard constructor
65 fLines = new TClonesArray("AliStrLine",1000);
66 SetCoarseDiffPhiCut();
67 SetMaxRCut();
68 SetZCutDiamond();
69 SetDCAcut();
70 SetDiffPhiMax();
71}
72
73//______________________________________________________________________
74AliITSVertexer3D::~AliITSVertexer3D() {
75 // Destructor
76 fVert3D = 0; // this object is not owned by this class
77 if(fLines){
78 fLines->Delete();
79 delete fLines;
80 }
81}
82
83//______________________________________________________________________
84AliESDVertex* AliITSVertexer3D::FindVertexForCurrentEvent(Int_t evnumber){
85 // Defines the AliESDVertex for the current event
86 if (fDebug) cout<<"\n \n FindVertexForCurrentEvent - 3D - PROCESSING EVENT "<<evnumber<<endl;
87 fVert3D = 0;
88 if(fLines)fLines->Clear();
89
90 Int_t nolines = FindTracklets(evnumber,0);
91 fCurrentVertex = 0;
92 if(nolines<2)return fCurrentVertex;
93 Int_t rc=Prepare3DVertex();
94 if(rc==0)Find3DVertex();
95 if(fVert3D){
96 if(fLines) fLines->Delete();
97 nolines = FindTracklets(evnumber,1);
98 if(nolines<2)return fCurrentVertex;
99 rc=Prepare3DVertex();
100 if(rc==0)Find3DVertex();
101 }
102
103 if(fVert3D){
104 fCurrentVertex = new AliESDVertex();
105 fCurrentVertex->SetTitle("vertexer: 3D");
106 fCurrentVertex->SetName("Vertex");
107 fCurrentVertex->SetXv(fVert3D->GetXv());
108 fCurrentVertex->SetYv(fVert3D->GetYv());
109 fCurrentVertex->SetZv(fVert3D->GetZv());
110 fCurrentVertex->SetDispersion(fVert3D->GetDispersion());
111 fCurrentVertex->SetNContributors(fVert3D->GetNContributors());
112 }
113 FindMultiplicity(evnumber);
114 return fCurrentVertex;
115}
116
117//______________________________________________________________________
118Int_t AliITSVertexer3D::FindTracklets(Int_t evnumber, Int_t optCuts){
119 // All the possible combinations between recpoints on layer 1and 2 are
120 // considered. Straight lines (=tracklets)are formed.
121 // The tracklets are processed in Prepare3DVertex
122 AliRunLoader *rl =AliRunLoader::GetRunLoader();
123 AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
124 AliITSgeom* geom = itsLoader->GetITSgeom();
125 itsLoader->LoadRecPoints();
126 rl->GetEvent(evnumber);
127
128 AliITSDetTypeRec detTypeRec;
129
130 TTree *tR = itsLoader->TreeR();
131 detTypeRec.SetTreeAddressR(tR);
132 TClonesArray *itsRec = 0;
133 // lc and gc are local and global coordinates for layer 1
134 Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
135 Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
136 // lc2 and gc2 are local and global coordinates for layer 2
137 Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
138 Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
139
140 itsRec = detTypeRec.RecPoints();
141 TBranch *branch;
142 branch = tR->GetBranch("ITSRecPoints");
143
144 // Set values for cuts
145 Float_t xbeam=0., ybeam=0.;
146 Float_t deltaPhi=fCoarseDiffPhiCut;
147 if(optCuts){
148 xbeam=fVert3D->GetXv();
149 ybeam=fVert3D->GetYv();
150 deltaPhi = fDiffPhiMax;
151 }
152
153 Int_t nrpL1 = 0; // number of rec points on layer 1
154 Int_t nrpL2 = 0; // number of rec points on layer 2
155
156 // By default irstL1=0 and lastL1=79
157 Int_t irstL1 = geom->GetStartDet(0);
158 Int_t lastL1 = geom->GetModuleIndex(2,1,1)-1;
159 for(Int_t module= irstL1; module<=lastL1;module++){ // count number of recopints on layer 1
160 branch->GetEvent(module);
161 nrpL1+= itsRec->GetEntries();
162 detTypeRec.ResetRecPoints();
163 }
164 //By default irstL2=80 and lastL2=239
165 Int_t irstL2 = geom->GetModuleIndex(2,1,1);
166 Int_t lastL2 = geom->GetLastDet(0);
167 for(Int_t module= irstL2; module<=lastL2;module++){ // count number of recopints on layer 2
168 branch->GetEvent(module);
169 nrpL2+= itsRec->GetEntries();
170 detTypeRec.ResetRecPoints();
171 }
172 if(nrpL1 == 0 || nrpL2 == 0){
173 itsLoader->UnloadRecPoints();
174 return -1;
175 }
176 if(fDebug) printf("RecPoints on Layer 1,2 = %d, %d\n",nrpL1,nrpL2);
177
178 Double_t a[3]={xbeam,ybeam,0.};
179 Double_t b[3]={xbeam,ybeam,10.};
180 AliStrLine zeta(a,b,kTRUE);
181
182 Int_t nolines = 0;
183 // Loop on modules of layer 1
184 for(Int_t modul1= irstL1; modul1<=lastL1;modul1++){ // Loop on modules of layer 1
185 branch->GetEvent(modul1);
186 Int_t nrecp1 = itsRec->GetEntries();
187 TClonesArray *prpl1 = new TClonesArray("AliITSRecPoint",nrecp1);
188 prpl1->SetOwner();
189 TClonesArray &rpl1 = *prpl1;
190 for(Int_t j=0;j<nrecp1;j++){
191 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
192 new(rpl1[j])AliITSRecPoint(*recp);
193 }
194 detTypeRec.ResetRecPoints();
195 for(Int_t j=0;j<nrecp1;j++){
196 AliITSRecPoint *recp = (AliITSRecPoint*)prpl1->At(j);
197 // Local coordinates of this recpoint
198 lc[0]=recp->GetDetLocalX();
199 lc[2]=recp->GetDetLocalZ();
200 geom->LtoG(modul1,lc,gc); // global coordinates
201 Double_t phi1 = TMath::ATan2(gc[1]-ybeam,gc[0]-xbeam);
202 if(phi1<0)phi1=2*TMath::Pi()+phi1;
203 for(Int_t modul2= irstL2; modul2<=lastL2;modul2++){
204 branch->GetEvent(modul2);
205 Int_t nrecp2 = itsRec->GetEntries();
206 for(Int_t j2=0;j2<nrecp2;j2++){
207 recp = (AliITSRecPoint*)itsRec->At(j2);
208 lc2[0]=recp->GetDetLocalX();
209 lc2[2]=recp->GetDetLocalZ();
210 geom->LtoG(modul2,lc2,gc2);
211 Double_t phi2 = TMath::ATan2(gc2[1]-ybeam,gc2[0]-xbeam);
212 if(phi2<0)phi2=2*TMath::Pi()+phi2;
213 Double_t diff = TMath::Abs(phi2-phi1);
214 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
215 if(diff>deltaPhi)continue;
216 AliStrLine line(gc,gc2,kTRUE);
217 Double_t cp[3];
218 Int_t retcode = line.Cross(&zeta,cp);
219 if(retcode<0)continue;
220 Double_t dca = line.GetDCA(&zeta);
221 if(dca<0.) continue;
222 if(dca>fMaxRCut)continue;
223 if(cp[2]>fZCutDiamond || cp[2]<-fZCutDiamond)continue;
224 MakeTracklet(gc,gc2,nolines);
225 }
226 detTypeRec.ResetRecPoints();
227 }
228 }
229 delete prpl1;
230 }
231 if(nolines == 0)return -2;
232 return nolines;
233}
234
235//______________________________________________________________________
236void AliITSVertexer3D::FindVertices(){
237 // computes the vertices of the events in the range FirstEvent - LastEvent
238 AliRunLoader *rl = AliRunLoader::GetRunLoader();
239 AliITSLoader* itsLoader = (AliITSLoader*) rl->GetLoader("ITSLoader");
240 itsLoader->ReloadRecPoints();
241 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
242 rl->GetEvent(i);
243 FindVertexForCurrentEvent(i);
244 if(fCurrentVertex){
245 WriteCurrentVertex();
246 }
247 }
248}
249
250
251//______________________________________________________________________
252void AliITSVertexer3D::Find3DVertex(){
253 // Determines the vertex position
254 // same algorithm as in AliVertexerTracks::StrLinVertexFinderMinDist(0)
255 // adapted to pure AliStrLine objects
256
257 Int_t knacc = fLines->GetEntriesFast();
258 Double_t (*vectP0)[3]=new Double_t [knacc][3];
259 Double_t (*vectP1)[3]=new Double_t [knacc][3];
260
261 Double_t sum[3][3];
262 Double_t dsum[3]={0,0,0};
263 for(Int_t i=0;i<3;i++)
264 for(Int_t j=0;j<3;j++)sum[i][j]=0;
265 for(Int_t i=0; i<knacc; i++){
266 AliStrLine *line1 = (AliStrLine*)fLines->At(i);
267
268 Double_t p0[3],cd[3];
269 line1->GetP0(p0);
270 line1->GetCd(cd);
271 Double_t p1[3]={p0[0]+cd[0],p0[1]+cd[1],p0[2]+cd[2]};
272 vectP0[i][0]=p0[0];
273 vectP0[i][1]=p0[1];
274 vectP0[i][2]=p0[2];
275 vectP1[i][0]=p1[0];
276 vectP1[i][1]=p1[1];
277 vectP1[i][2]=p1[2];
278
279 Double_t matr[3][3];
280 Double_t dknow[3];
281 AliVertexerTracks::GetStrLinDerivMatrix(p0,p1,matr,dknow);
282
283 for(Int_t iii=0;iii<3;iii++){
284 dsum[iii]+=dknow[iii];
285 for(Int_t lj=0;lj<3;lj++) sum[iii][lj]+=matr[iii][lj];
286 }
287 }
288
289 Double_t vett[3][3];
290 Double_t det=AliVertexerTracks::GetDeterminant3X3(sum);
291 Double_t initPos[3];
292 Double_t sigma = 0.;
293 for(Int_t i=0;i<3;i++)initPos[i]=0.;
294
295 Int_t goodvert=0;
296
297 if(det!=0){
298 for(Int_t zz=0;zz<3;zz++){
299 for(Int_t ww=0;ww<3;ww++){
300 for(Int_t kk=0;kk<3;kk++) vett[ww][kk]=sum[ww][kk];
301 }
302 for(Int_t kk=0;kk<3;kk++) vett[kk][zz]=dsum[kk];
303 initPos[zz]=AliVertexerTracks::GetDeterminant3X3(vett)/det;
304 }
305
306 Double_t rvert=TMath::Sqrt(initPos[0]*initPos[0]+initPos[1]*initPos[1]);
307 Double_t zvert=initPos[2];
308 if(rvert<fMaxRCut && TMath::Abs(zvert)<fZCutDiamond){
309 for(Int_t i=0; i<knacc; i++){
310 Double_t p0[3]={0,0,0},p1[3]={0,0,0};
311 for(Int_t ii=0;ii<3;ii++){
312 p0[ii]=vectP0[i][ii];
313 p1[ii]=vectP1[i][ii];
314 }
315 sigma+=AliVertexerTracks::GetStrLinMinDist(p0,p1,initPos);
316 }
317 sigma=TMath::Sqrt(sigma);
318 goodvert=1;
319 }
320 }
321 delete vectP0;
322 delete vectP1;
323 if(!goodvert){
324 Warning("Find3DVertex","Finder did not succed");
325 for(Int_t i=0;i<3;i++)initPos[i]=0.;
326 sigma=999;
327 knacc=-1;
328 }
329 if(fVert3D){
330 fVert3D-> SetXYZ(initPos);
331 fVert3D->SetDispersion(sigma);
332 fVert3D->SetNContributors(knacc);
333 }else{
334 fVert3D = new AliVertex(initPos,sigma,knacc);
335 }
336 if(fDebug) fVert3D->Print();
337}
338
339
340
341//______________________________________________________________________
342Int_t AliITSVertexer3D::Prepare3DVertex(){
343 // Finds the 3D vertex information using tracklets
344 Int_t retcode = -1;
345
346 Int_t nbr=50;
347 Float_t rl=-fMaxRCut;
348 Float_t rh=fMaxRCut;
349 Int_t nbz=100;
350 Float_t zl=-fZCutDiamond;
351 Float_t zh=fZCutDiamond;
352 Float_t binsizer=(rh-rl)/nbr;
353 Float_t binsizez=(zh-zl)/nbz;
354 TH3F *h3d = new TH3F("h3d","xyz distribution",nbr,rl,rh,nbr,rl,rh,nbz,zl,zh);
355
356
357 // cleanup of the TCLonesArray of tracklets (i.e. fakes are removed)
358 Int_t *validate = new Int_t [fLines->GetEntriesFast()];
359 for(Int_t i=0; i<fLines->GetEntriesFast();i++)validate[i]=0;
360 // Int_t itot=fLines->GetEntriesFast()-1-1;
361 for(Int_t i=0; i<fLines->GetEntriesFast()-1;i++){
362 //if(validate[i]==1)continue;
363 AliStrLine *l1 = (AliStrLine*)fLines->At(i);
364 for(Int_t j=i+1;j<fLines->GetEntriesFast();j++){
365 AliStrLine *l2 = (AliStrLine*)fLines->At(j);
366 if(l1->GetDCA(l2) > fDCAcut) continue;
367 Double_t point[3];
368 Int_t retc = l1->Cross(l2,point);
369 if(retc<0)continue;
370 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
371 if(TMath::Abs(point[2])>fZCutDiamond)continue;
372 if(rad>fMaxRCut)continue;
373 validate[i]++;
374 validate[j]++;
375 h3d->Fill(point[0],point[1],point[2]);
376 }
377 }
378
379
380
381 for(Int_t i=0; i<fLines->GetEntriesFast()-1;i++){ // remove tracklets sharing one recpoint
382 if(validate[i]==0)continue;
383 AliStrLine *l1 = (AliStrLine*)fLines->At(i);
384 for(Int_t j=i+1;j<fLines->GetEntriesFast();j++){
385 if(validate[j]==0)continue;
386 AliStrLine *l2 = (AliStrLine*)fLines->At(j);
387 if(l1->GetDCA(l2) > 0.00001) continue;
388 Double_t point[3];
389 Int_t retc = l1->Cross(l2,point);
390 if(retc<0)continue;
391 Double_t rad=TMath::Sqrt(point[0]*point[0]+point[1]*point[1]);
392 if(rad<fMaxRCut)continue;
393 validate[i]--;
394 validate[j]--;
395 }
396 }
397
398 Int_t numbtracklets=0;
399 for(Int_t i=0; i<fLines->GetEntriesFast();i++)if(validate[i]>=1)numbtracklets++;
400 if(numbtracklets<2){delete [] validate; delete h3d; return retcode; }
401
402 for(Int_t i=0; i<fLines->GetEntriesFast();i++){
403 if(validate[i]<1)fLines->RemoveAt(i);
404 }
405 fLines->Compress();
406 if (fDebug) cout<<"Number of tracklets (after compress) "<<fLines->GetEntriesFast()<<endl;
407 delete [] validate;
408
409
410 // finds region of origin"
411 Int_t nobi = 0;
412 Float_t cont = 0.;
413 Float_t cont2 = 0.;
414 TAxis *xax = h3d->GetXaxis();
415 TAxis *yax = h3d->GetYaxis();
416 TAxis *zax = h3d->GetZaxis();
417 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
418 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
419 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
420 nobi++;
421 cont+=h3d->GetBinContent(i,j,k);
422 cont2+=h3d->GetBinContent(i,j,k)*h3d->GetBinContent(i,j,k);
423 }
424 }
425 }
426 cont/=nobi;
427 cont2/=nobi;
428 cont2 = TMath::Sqrt((cont2 - cont*cont)/(nobi-1));
429 Double_t peak[3]={0.,0.,0.};
430 Float_t contref = 0.;
431 for(Int_t i=xax->GetFirst();i<=xax->GetLast();i++){
432 Float_t xval = xax->GetBinCenter(i);
433 for(Int_t j=yax->GetFirst();j<=yax->GetLast();j++){
434 Float_t yval = yax->GetBinCenter(j);
435 for(Int_t k=zax->GetFirst();k<=zax->GetLast();k++){
436 Float_t bc = h3d->GetBinContent(i,j,k);
437 Float_t zval = zax->GetBinCenter(k);
438 if(bc>contref){
439 contref = bc;
440 peak[2] = zval;
441 peak[1] = yval;
442 peak[0] = xval;
443 }
444 }
445 }
446 }
447 delete h3d;
448
449 // Second selection loop
450 Float_t bs=(binsizer+binsizez)/2.;
451 for(Int_t i=0; i<fLines->GetEntriesFast()-1;i++){
452 AliStrLine *l1 = (AliStrLine*)fLines->At(i);
453 if(l1->GetDistFromPoint(peak)>2.5*bs)fLines->RemoveAt(i);
454 }
455 fLines->Compress();
456 if (fDebug) cout<<"Number of tracklets (after 2nd compression) "<<fLines->GetEntriesFast()<<endl;
457
458
459 if(fLines->GetEntriesFast()>1){
460 Find3DVertex(); // find a first candidate for the primary vertex
461 // make a further selection on tracklets based on this first candidate
462 fVert3D->GetXYZ(peak);
463 if (fDebug) cout<<"FIRST V candidate: "<<peak[0]<<"; "<<peak[1]<<"; "<<peak[2]<<endl;
464 for(Int_t i=0; i<fLines->GetEntriesFast()-1;i++){
465 AliStrLine *l1 = (AliStrLine*)fLines->At(i);
466 if(l1->GetDistFromPoint(peak)> fDCAcut)fLines->RemoveAt(i);
467 }
468 fLines->Compress();
469 if (fDebug) cout<<"Number of tracklets (after 3rd compression) "<<fLines->GetEntriesFast()<<endl;
470 if(fLines->GetEntriesFast()>1){
471 retcode = 0; // this last tracklet selection will be used
472 delete fVert3D;
473 fVert3D = 0;
474 }
475 else {
476 retcode =1; // the previous tracklet selection will be used
477 }
478 }
479 else {
480 retcode = 0;
481 }
482 return retcode;
483}
484
485 //______________________________________________________________________
486void AliITSVertexer3D::MakeTracklet(Double_t *pA, Double_t *pB, Int_t &nolines) {
487 // makes a tracklet
488 TClonesArray &lines = *fLines;
489 if(nolines == 0){
490 if(fLines->GetEntriesFast()>0)fLines->Clear();
491 }
492 if(fLines->GetEntriesFast()==fLines->GetSize()){
493 Int_t newsize=(Int_t) 1.5*fLines->GetEntriesFast();
494 fLines->Expand(newsize);
495 }
496
497 new(lines[nolines++])AliStrLine(pA,pB,kTRUE);
498}
499
500 //______________________________________________________________________
501void AliITSVertexer3D::MakeTracklet(Float_t *pA, Float_t *pB, Int_t &nolines) {// Makes a tracklet
502 //
503 Double_t a[3],b[3];
504 for(Int_t i=0;i<3;i++){
505 a[i] = pA[i];
506 b[i] = pB[i];
507 }
508 MakeTracklet(a,b,nolines);
509}
510
511//________________________________________________________
512void AliITSVertexer3D::PrintStatus() const {
513 // Print current status
514 cout <<"=======================================================\n";
515 cout << "Loose cut on Delta Phi "<<fCoarseDiffPhiCut<<endl;
516 cout << "Cut on tracklet DCA to Z axis "<<fMaxRCut<<endl;
517 cout << "Cut on diamond (Z) "<<fZCutDiamond<<endl;
518 cout << "Cut on DCA - tracklet to tracklet and to vertex "<<fDCAcut<<endl;
519 cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
520
521
522 cout <<"=======================================================\n";
523}