]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerCosmics.cxx
New Class AliITSTrackleterSPDEff:
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerCosmics.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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
16 #include <TClonesArray.h>
17 #include "AliLog.h"
18 #include "AliESDVertex.h"
19 #include "AliRunLoader.h"
20 #include "AliITSLoader.h"
21 #include "AliITSgeomTGeo.h"
22 #include "AliITSRecPoint.h"
23 #include "AliITSReconstructor.h"
24 #include "AliITSVertexerCosmics.h"
25 #include "AliStrLine.h"
26
27 //------------------------------------------------------------------------
28 // This class implements a method to construct a "fake" primary
29 // vertex for cosmic events in which the muon crosses one of 5 inner
30 // ITS layers. A fake primary vertex is needed for the reconstruction,
31 // with e.g. AliITStrackerSA, of the two tracks produced by the muon 
32 // in the ITS.
33 //   We build pairs of clusters on a given layer and define the fake vertex as
34 // the mid-point of the straight line joining the two clusters.
35 //   We use the innermost layer that has at least two clusters.
36 //   We reject the background by requiring at least one cluster on the outer
37 // layer, closer than fMaxDistOnOuterLayer to the tracklet prolongation.
38 //   We can reject (potentially pathological) events with the muon track
39 // tangential to the layer by the requiring the radial position of
40 // the vertex to be smaller than fMaxVtxRadius.
41 //   Due to background clusters, more than one vertex per event can 
42 // be found. We consider the first found.
43 //   The errors on x,y,z of the vertex are calculated as errors on the mean
44 // of clusters coordinates. Non-diag elements of vertex cov. mat. are set to 0.
45 //   The number of contributors set in the AliESDVertex object is the
46 // number of the layer on which the tracklet was built; if this number is -1, 
47 // the procedure could not find a vertex position and by default 
48 // the vertex coordinates are set to (0,0,0) with large errors (100,100,100)
49 //
50 // Origin: A.Dainese, andrea.dainese@lnl.infn.it
51 //-------------------------------------------------------------------------
52
53 ClassImp(AliITSVertexerCosmics)
54
55 //-------------------------------------------------------------------------
56 AliITSVertexerCosmics::AliITSVertexerCosmics():AliITSVertexer(),
57 fMaxDistOnOuterLayer(0),
58 fMinDist2Vtxs(0)
59 {
60   // Default constructor
61   SetFirstLastModules(0,0,79);
62   SetFirstLastModules(1,80,239);
63   SetFirstLastModules(2,240,323);
64   SetFirstLastModules(3,324,499);
65   SetFirstLastModules(4,500,1247);
66   SetFirstLastModules(5,1248,2197);
67   /*
68   SetMaxVtxRadius(0,3.5);
69   SetMaxVtxRadius(1,6.5);
70   SetMaxVtxRadius(2,14.5);
71   SetMaxVtxRadius(3,23.5);
72   SetMaxVtxRadius(4,37.5);
73   SetMaxVtxRadius(5,42.5);
74   */  
75   SetMaxVtxRadius(0,5.5);
76   SetMaxVtxRadius(1,8.5);
77   SetMaxVtxRadius(2,18.5);
78   SetMaxVtxRadius(3,28.5);
79   SetMaxVtxRadius(4,39.5);
80   SetMaxVtxRadius(5,48.5);
81   
82   SetMaxDistOnOuterLayer();
83   SetMinDist2Vtxs();
84 }
85 //--------------------------------------------------------------------------
86 AliESDVertex* AliITSVertexerCosmics::FindVertexForCurrentEvent(Int_t evnumber) 
87 {
88   // Defines the AliESDVertex for the current event
89
90   fCurrentVertex = 0;
91   AliRunLoader *rl =AliRunLoader::GetRunLoader();
92   AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
93   itsLoader->LoadRecPoints();
94   rl->GetEvent(evnumber);
95
96   TTree *rpTree = itsLoader->TreeR();
97
98   TClonesArray *recpoints=new TClonesArray("AliITSRecPoint",10000);
99   rpTree->SetBranchAddress("ITSRecPoints",&recpoints);
100
101   Int_t lay,lad,det; 
102
103   Double_t pos[3]={0.,0.,0.};
104   Double_t err[3]={100.,100.,100.};
105   Int_t ncontributors = -1;
106
107   // Search for innermost layer with at least two clusters 
108   // on two different modules
109   Int_t ilayer=0,ilayer2=0;
110   while(ilayer<6) {
111     if(AliITSReconstructor::GetRecoParam()->GetLayersToSkip(ilayer)) {
112       ilayer++;
113       continue;
114     }
115     Int_t nHitModules=0;
116     for(Int_t imodule=fFirst[ilayer]; imodule<=fLast[ilayer]; imodule++) {
117       rpTree->GetEvent(imodule);
118       AliITSgeomTGeo::GetModuleId(imodule,lay,lad,det);
119       lay -= 1;  // AliITSgeomTGeo gives layer from 1 to 6, we want 0 to 5
120       if(lay!=ilayer) AliFatal("Layer mismatch!");
121       if(recpoints->GetEntriesFast()>0) nHitModules++;
122     }
123     if(nHitModules>=2) break;
124     ilayer++;
125   }
126
127   ilayer2=ilayer+1;
128   while(ilayer2<6) {
129     if(!AliITSReconstructor::GetRecoParam()->GetLayersToSkip(ilayer2)) break;
130     ilayer2++;
131   }
132
133   // try tracklet on SPD2 and point on SPD1
134   if(ilayer==1 && ilayer2>5 && !AliITSReconstructor::GetRecoParam()->GetLayersToSkip(0)) {ilayer=0; ilayer2=1;}
135
136   if(ilayer>4 || ilayer2>5) {
137     AliWarning("Not enough clusters");
138     delete recpoints;
139     itsLoader->UnloadRecPoints();
140     fCurrentVertex = new AliESDVertex(pos,err,"cosmics");
141     fCurrentVertex->SetTitle("cosmics fake vertex (failed)");
142     fCurrentVertex->SetNContributors(ncontributors);
143     return fCurrentVertex;
144   }
145
146
147   const Int_t arrSize = 1000;
148   Float_t xclInnLay[arrSize],yclInnLay[arrSize],zclInnLay[arrSize],modclInnLay[arrSize];
149   Float_t e2xclInnLay[arrSize],e2yclInnLay[arrSize],e2zclInnLay[arrSize];
150   Float_t e2xclOutLay[arrSize],e2yclOutLay[arrSize],e2zclOutLay[arrSize];
151   Int_t nclInnLayStored=0;
152   Float_t xclOutLay[arrSize],yclOutLay[arrSize],zclOutLay[arrSize],modclOutLay[arrSize];
153   Int_t nclOutLayStored=0;
154   Int_t nRecPoints,nRecPointsInnLay=0;
155
156   Float_t gc[3],gcov[5];
157
158   Float_t matchOutLayValue;
159   Float_t distxyInnLay,distxyInnLayBest=0.;
160   Double_t p1[3],p2[3],p3[3];
161
162   Float_t xvtx,yvtx,zvtx,rvtx;
163   Int_t i1InnLayBest=-1,i2InnLayBest=-1;
164
165
166   // Collect clusters in the selected layer and the outer one
167   for(Int_t imodule=fFirst[ilayer]; imodule<=fLast[ilayer2]; imodule++) {
168     rpTree->GetEvent(imodule);
169     AliITSgeomTGeo::GetModuleId(imodule,lay,lad,det);
170     lay -= 1; // AliITSgeomTGeo gives layer from 1 to 6, we want 0 to 5
171     nRecPoints=recpoints->GetEntriesFast();
172     if(imodule<=fLast[ilayer]) nRecPointsInnLay += nRecPoints;
173     //printf("cosmics: module %d clusters %d\n",imodule,nRecPoints);
174     for(Int_t irp=0; irp<nRecPoints; irp++) {
175       AliITSRecPoint *rp=(AliITSRecPoint*)recpoints->UncheckedAt(irp);
176       // Local coordinates of this recpoint
177       rp->GetGlobalXYZ(gc);
178       if(lay==ilayer) { // store InnLay clusters
179         xclInnLay[nclInnLayStored]=gc[0];
180         yclInnLay[nclInnLayStored]=gc[1];
181         zclInnLay[nclInnLayStored]=gc[2];
182         rp->GetGlobalCov(gcov);
183         e2xclInnLay[nclInnLayStored]=gcov[0];
184         e2yclInnLay[nclInnLayStored]=gcov[3];
185         e2zclInnLay[nclInnLayStored]=gcov[5];
186         modclInnLay[nclInnLayStored]=imodule;
187         nclInnLayStored++;
188       }
189       if(lay==ilayer2) { // store OutLay clusters
190         xclOutLay[nclOutLayStored]=gc[0];
191         yclOutLay[nclOutLayStored]=gc[1];
192         zclOutLay[nclOutLayStored]=gc[2];
193         rp->GetGlobalCov(gcov);
194         e2xclOutLay[nclOutLayStored]=gcov[0];
195         e2yclOutLay[nclOutLayStored]=gcov[3];
196         e2zclOutLay[nclOutLayStored]=gcov[5];
197         modclOutLay[nclOutLayStored]=imodule;
198         nclOutLayStored++;
199       }
200       if(nclInnLayStored>arrSize || nclOutLayStored>arrSize) {
201         //AliFatal("More than arrSize clusters per layer");
202         AliWarning("Too many clusters per layer");
203         delete recpoints;
204         itsLoader->UnloadRecPoints();
205         fCurrentVertex = new AliESDVertex(pos,err,"cosmics");
206         fCurrentVertex->SetTitle("cosmics fake vertex (failed)");
207         fCurrentVertex->SetNContributors(ncontributors);
208         return fCurrentVertex;
209       }
210     }// end clusters in a module
211   }// end modules
212
213
214   Int_t i1InnLay,i2InnLay,iOutLay;
215
216   // build fake vertices
217   printf("Building tracklets on layer %d\n",ilayer);
218
219   // InnLay - first cluster
220   for(i1InnLay=0; i1InnLay<nclInnLayStored; i1InnLay++) { 
221     p1[0]=xclInnLay[i1InnLay]; 
222     p1[1]=yclInnLay[i1InnLay]; 
223     p1[2]=zclInnLay[i1InnLay];
224     // InnLay - second cluster
225     for(i2InnLay=i1InnLay+1; i2InnLay<nclInnLayStored; i2InnLay++) { 
226       if(modclInnLay[i1InnLay]==modclInnLay[i2InnLay]) continue;
227       p2[0]=xclInnLay[i2InnLay]; 
228       p2[1]=yclInnLay[i2InnLay]; 
229       p2[2]=zclInnLay[i2InnLay];
230       // look for point on OutLay
231       AliStrLine innLayline(p1,p2,kTRUE);
232       for(iOutLay=0; iOutLay<nclOutLayStored; iOutLay++) {
233         p3[0]=xclOutLay[iOutLay]; 
234         p3[1]=yclOutLay[iOutLay]; 
235         p3[2]=zclOutLay[iOutLay];
236         //printf("(%f,%f) (%f,%f)     (%f,%f) %f\n",p1[0],p1[1],p2[0],p2[1],p3[0],p3[1],innLayline.GetDistFromPoint(p3));
237         matchOutLayValue=innLayline.GetDistFromPoint(p3);
238         distxyInnLay = (p1[0]-p2[0])*(p1[0]-p2[0])+(p1[1]-p2[1])*(p1[1]-p2[1]);
239         if(matchOutLayValue<fMaxDistOnOuterLayer &&
240            distxyInnLay>distxyInnLayBest) { 
241           //printf("found\n");
242           distxyInnLayBest=distxyInnLay;
243           i1InnLayBest=i1InnLay;
244           i2InnLayBest=i2InnLay;
245         }
246       }
247     } // InnLay - second cluster
248   } // InnLay - first cluster
249
250   if(i1InnLayBest>-1 && i2InnLayBest>-1) { 
251     xvtx = 0.5*(xclInnLay[i1InnLayBest]+xclInnLay[i2InnLayBest]);
252     yvtx = 0.5*(yclInnLay[i1InnLayBest]+yclInnLay[i2InnLayBest]);
253     zvtx = 0.5*(zclInnLay[i1InnLayBest]+zclInnLay[i2InnLayBest]);
254     rvtx = TMath::Sqrt(xvtx*xvtx+yvtx*yvtx);
255     if(rvtx<fMaxVtxRadius[ilayer]) {
256       ncontributors = ilayer;
257       pos[0] = xvtx;
258       pos[1] = yvtx;
259       pos[2] = zvtx;
260       err[0]=TMath::Sqrt(0.25*(e2xclInnLay[i1InnLayBest]+e2xclInnLay[i2InnLayBest])); 
261       err[1]=TMath::Sqrt(0.25*(e2yclInnLay[i1InnLayBest]+e2yclInnLay[i2InnLayBest])); 
262       err[2]=TMath::Sqrt(0.25*(e2zclInnLay[i1InnLayBest]+e2zclInnLay[i2InnLayBest]));
263     }
264
265   } else { // give it a try exchanging InnLay and OutLay
266
267     // OutLay - first cluster
268     for(i1InnLay=0; i1InnLay<nclOutLayStored; i1InnLay++) { 
269       p1[0]=xclOutLay[i1InnLay]; 
270       p1[1]=yclOutLay[i1InnLay]; 
271       p1[2]=zclOutLay[i1InnLay];
272       // OutLay - second cluster
273       for(i2InnLay=i1InnLay+1; i2InnLay<nclOutLayStored; i2InnLay++) { 
274         if(modclOutLay[i1InnLay]==modclOutLay[i2InnLay]) continue;
275         p2[0]=xclOutLay[i2InnLay]; 
276         p2[1]=yclOutLay[i2InnLay]; 
277         p2[2]=zclOutLay[i2InnLay];
278         // look for point on InnLay
279         AliStrLine outLayline(p1,p2,kTRUE);
280         for(iOutLay=0; iOutLay<nclInnLayStored; iOutLay++) {
281           p3[0]=xclInnLay[iOutLay]; 
282           p3[1]=yclInnLay[iOutLay]; 
283           p3[2]=zclInnLay[iOutLay];
284           //printf(" %f\n",InnLayline.GetDistFromPoint(p3));
285           matchOutLayValue=outLayline.GetDistFromPoint(p3);
286           distxyInnLay = (p1[0]-p2[0])*(p1[0]-p2[0])+(p1[1]-p2[1])*(p1[1]-p2[1]);
287           if(matchOutLayValue<fMaxDistOnOuterLayer &&
288              distxyInnLay>distxyInnLayBest) { 
289             distxyInnLayBest=distxyInnLay;
290             i1InnLayBest=i1InnLay;
291             i2InnLayBest=i2InnLay;
292           }
293         }
294       } // OutLay - second cluster
295     } // OutLay - first cluster
296
297     if(i1InnLayBest>-1 && i2InnLayBest>-1) { 
298       xvtx = 0.5*(xclOutLay[i1InnLayBest]+xclOutLay[i2InnLayBest]);
299       yvtx = 0.5*(yclOutLay[i1InnLayBest]+yclOutLay[i2InnLayBest]);
300       zvtx = 0.5*(zclOutLay[i1InnLayBest]+zclOutLay[i2InnLayBest]);
301       rvtx = TMath::Sqrt(xvtx*xvtx+yvtx*yvtx);
302       if(rvtx<fMaxVtxRadius[ilayer]) {
303         ncontributors = ilayer2;
304         pos[0] = xvtx;
305         pos[1] = yvtx;
306         pos[2] = zvtx;
307         err[0]=TMath::Sqrt(0.25*(e2xclOutLay[i1InnLayBest]+e2xclOutLay[i2InnLayBest])); 
308         err[1]=TMath::Sqrt(0.25*(e2yclOutLay[i1InnLayBest]+e2yclOutLay[i2InnLayBest])); 
309         err[2]=TMath::Sqrt(0.25*(e2zclOutLay[i1InnLayBest]+e2zclOutLay[i2InnLayBest]));
310       }
311     }
312   } // give it a try exchanging InnLay and OutLay
313
314   fCurrentVertex = new AliESDVertex(pos,err,"cosmics");
315   fCurrentVertex->SetTitle("cosmics fake vertex");
316   fCurrentVertex->SetNContributors(ncontributors);
317   //fCurrentVertex->Print();
318
319   delete recpoints;
320   itsLoader->UnloadRecPoints();
321
322   return fCurrentVertex;
323 }  
324 //-------------------------------------------------------------------------
325 void AliITSVertexerCosmics::FindVertices()
326 {
327   // computes the vertices of the events in the range FirstEvent - LastEvent
328   AliRunLoader *rl = AliRunLoader::GetRunLoader();
329   AliITSLoader* itsLoader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
330   itsLoader->ReloadRecPoints();
331   for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
332     //  printf("Processing event %d\n",i);
333     rl->GetEvent(i);
334     FindVertexForCurrentEvent(i);
335     if(fCurrentVertex){
336       WriteCurrentVertex();
337     }
338   }
339 }
340 //-------------------------------------------------------------------------
341 void AliITSVertexerCosmics::PrintStatus() const 
342 {
343   // Print current status
344   printf("=======================================================\n");
345   printf(" fMaxDistOnOuterLayer: %f\n",fMaxDistOnOuterLayer);
346   printf(" fMaxVtxRadius[0]:  %f\n",fMaxVtxRadius[0]);
347   printf(" fMinDist2Vtxs:  %f\n",fMinDist2Vtxs);
348   printf("=======================================================\n");
349 }
350 //-------------------------------------------------------------------------