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