]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerZ.cxx
5d274cb3498640d1de018ffbd140358ac67066b7
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerZ.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2003, 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 <AliITSVertexerZ.h>
16 #include <TString.h>
17 #include "AliITSLoader.h"
18 #include<TBranch.h>
19 #include<TClonesArray.h>
20 #include<TH1.h>
21 #include<TTree.h>
22 #include <AliITSgeom.h>
23 #include "AliITSDetTypeRec.h"
24 #include <AliITSRecPoint.h>
25
26 /////////////////////////////////////////////////////////////////
27 // this class implements a fast method to determine
28 // the Z coordinate of the primary vertex
29 // for p-p collisions it seems to give comparable or better results
30 // with respect to what obtained with AliITSVertexerPPZ
31 // It can be used successfully with Pb-Pb collisions
32 ////////////////////////////////////////////////////////////////
33
34 ClassImp(AliITSVertexerZ)
35
36
37
38 //______________________________________________________________________
39 AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer() {
40   // Default constructor
41   SetDiffPhiMax(0);
42   fX0 = 0.;
43   fY0 = 0.;
44   SetFirstLayerModules(0);
45   SetSecondLayerModules(0);
46   fZFound = 0;
47   fZsig = 0.;
48   fZCombc = 0;
49   fZCombf = 0;
50   fZCombv = 0;
51   SetLowLimit(0.);
52   SetHighLimit(0.);
53   SetBinWidthCoarse(0.);
54   SetBinWidthFine(0.);
55   SetTolerance(0.);
56   SetPPsetting(0.,0.);
57 }
58
59 //______________________________________________________________________
60 AliITSVertexerZ::AliITSVertexerZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn) {
61   // Standard Constructor
62   SetDiffPhiMax();
63   fX0 = x0;
64   fY0 = y0;
65   SetFirstLayerModules();
66   SetSecondLayerModules();
67   fZFound = 0;
68   fZsig = 0.;
69   fZCombc = 0;
70   fZCombf = 0;
71   fZCombv = 0;
72   SetLowLimit();
73   SetHighLimit();
74   SetBinWidthCoarse();
75   SetBinWidthFine();
76   SetTolerance();
77   SetPPsetting();
78
79 }
80
81 //______________________________________________________________________
82 AliITSVertexerZ::AliITSVertexerZ(const AliITSVertexerZ &vtxr) : AliITSVertexer(vtxr) {
83   // Copy constructor
84   // Copies are not allowed. The method is protected to avoid misuse.
85   Error("AliITSVertexerZ","Copy constructor not allowed\n");
86 }
87
88 //______________________________________________________________________
89 AliITSVertexerZ& AliITSVertexerZ::operator=(const AliITSVertexerZ& /* vtxr */){
90   // Assignment operator
91   // Assignment is not allowed. The method is protected to avoid misuse.
92   Error("= operator","Assignment operator not allowed\n");
93   return *this;
94 }
95
96 //______________________________________________________________________
97 AliITSVertexerZ::~AliITSVertexerZ() {
98   // Destructor
99   delete fZCombc;
100   delete fZCombf;
101   delete fZCombv;
102 }
103
104 //______________________________________________________________________
105 AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(Int_t evnumber){
106   // Defines the AliESDVertex for the current event
107   
108   fCurrentVertex = 0;
109   AliRunLoader *rl =AliRunLoader::GetRunLoader();
110   AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
111   TDirectory * olddir = gDirectory;
112   rl->CdGAFile();
113   AliITSgeom* geom = (AliITSgeom*)gDirectory->Get("AliITSgeom");
114   olddir->cd();
115
116   itsLoader->LoadRecPoints();
117   rl->GetEvent(evnumber);
118
119   AliITSDetTypeRec detTypeRec;
120
121   TTree *tR = itsLoader->TreeR();
122   detTypeRec.SetTreeAddressR(tR);
123   TClonesArray *itsRec  = 0;
124   // lc and gc are local and global coordinates for layer 1
125   Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
126   Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
127   // lc2 and gc2 are local and global coordinates for layer 2
128   Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
129   Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
130
131   itsRec = detTypeRec.RecPoints();
132   TBranch *branch;
133   branch = tR->GetBranch("ITSRecPoints");
134
135   Int_t nrpL1 = 0;
136   Int_t nrpL2 = 0;
137   // By default fFirstL1=0 and fLastL1=79
138   // This loop counts the number of recpoints on layer1 (central modules)
139   for(Int_t module= fFirstL1; module<=fLastL1;module++){
140     // Keep only central modules
141     if(module%4==0 || module%4==3)continue;
142     //   cout<<"Procesing module "<<module<<" ";
143     branch->GetEvent(module);
144     //    cout<<"Number of clusters "<<clusters->GetEntries()<<endl;
145     nrpL1+= itsRec->GetEntries();
146     detTypeRec.ResetRecPoints();
147   }
148   //By default fFirstL2=80 and fLastL2=239
149   //This loop counts the number of RP on layer 2
150   for(Int_t module= fFirstL2; module<=fLastL2;module++){
151     branch->GetEvent(module);
152     nrpL2+= itsRec->GetEntries();
153     detTypeRec.ResetRecPoints();
154   }
155   if(nrpL1 == 0 || nrpL2 == 0){
156     ResetHistograms();
157     return fCurrentVertex;
158   }
159   // The vertex finding is attempted only if the number of RP is !=0 on
160   // both layers
161   Float_t *xc1 = new Float_t [nrpL1]; // coordinates of the L1 Recpoints
162   Float_t *yc1 = new Float_t [nrpL1];
163   Float_t *zc1 = new Float_t [nrpL1];
164   Float_t *phi1 = new Float_t [nrpL1];
165   Float_t *xc2 = new Float_t [nrpL2]; // coordinates of the L1 Recpoints
166   Float_t *yc2 = new Float_t [nrpL2];
167   Float_t *zc2 = new Float_t [nrpL2];
168   Float_t *phi2 = new Float_t [nrpL2];
169   Int_t ind = 0;// running index for RP
170   // Force a coarse bin size of 200 microns if the number of clusters on layer 2
171   // is low
172   if(nrpL2<fPPsetting[0])SetBinWidthCoarse(fPPsetting[1]);
173   // By default nbinfine = (10+10)/0.0005=40000
174   Int_t nbinfine = static_cast<Int_t>((fHighLim-fLowLim)/fStepFine);
175   // By default nbincoarse=(10+10)/0.01=2000
176   Int_t nbincoarse = static_cast<Int_t>((fHighLim-fLowLim)/fStepCoarse);
177   // Set stepverycoarse = 3*fStepCoarse
178   Int_t nbinvcoarse = static_cast<Int_t>((fHighLim-fLowLim)/(3.*fStepCoarse));
179   if(fZCombc)delete fZCombc;
180   fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
181   if(fZCombv)delete fZCombv;
182   fZCombv = new TH1F("fZCombv","Z",nbinvcoarse,fLowLim,fLowLim+nbinvcoarse*3.*fStepCoarse);
183   if(fZCombf)delete fZCombf;
184   fZCombf = new TH1F("fZCombf","Z",nbinfine,fLowLim,fLowLim+nbinfine*fStepFine);
185   // Loop on modules of layer 1 (restricted to central modules)
186   for(Int_t module= fFirstL1; module<=fLastL1;module++){
187     if(module%4==0 || module%4==3)continue;
188     branch->GetEvent(module);
189     Int_t nrecp1 = itsRec->GetEntries();
190     for(Int_t j=0;j<nrecp1;j++){
191       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
192       // Local coordinates of this recpoint
193       lc[0]=recp->GetDetLocalX();
194       lc[2]=recp->GetDetLocalZ();
195       geom->LtoG(module,lc,gc);
196       // Global coordinates of this recpoints
197       gc[0]-=fX0; // Possible beam offset in the bending plane
198       gc[1]-=fY0; //   "               "
199       xc1[ind]=gc[0];
200       yc1[ind]=gc[1];
201       zc1[ind]=gc[2];
202       // azimuthal angle is computed in the interval 0 --> 2*pi
203       phi1[ind] = TMath::ATan2(gc[1],gc[0]);
204       if(phi1[ind]<0)phi1[ind]=2*TMath::Pi()+phi1[ind];
205       ind++;
206     }
207     detTypeRec.ResetRecPoints();
208   }
209   ind = 0; // the running index is reset for Layer 2
210   for(Int_t module= fFirstL2; module<=fLastL2;module++){
211     branch->GetEvent(module);
212     Int_t nrecp2 = itsRec->GetEntries();
213     for(Int_t j=0;j<nrecp2;j++){
214       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
215       lc[0]=recp->GetDetLocalX();
216       lc[2]=recp->GetDetLocalZ();
217       geom->LtoG(module,lc,gc);
218       gc[0]-=fX0;
219       gc[1]-=fY0;
220       xc2[ind]=gc[0];
221       yc2[ind]=gc[1];
222       zc2[ind]=gc[2];
223       phi2[ind] = TMath::ATan2(gc[1],gc[0]);
224       if(phi2[ind]<0)phi2[ind]=2*TMath::Pi()+phi2[ind];
225       ind++;
226     }
227     detTypeRec.ResetRecPoints();
228   }
229   for(Int_t i=0;i<nrpL1;i++){ // loop on L1 RP
230     Float_t r1=TMath::Sqrt(xc1[i]*xc1[i]+yc1[i]*yc1[i]); // radius L1 RP
231     for(Int_t j=0;j<nrpL2;j++){ // loop on L2 RP
232       Float_t diff = TMath::Abs(phi2[j]-phi1[i]); // diff in azimuth
233       if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; //diff<pi
234       if(diff<fDiffPhiMax){ // cut on 10 milliradians by def.
235         Float_t r2=TMath::Sqrt(xc2[j]*xc2[j]+yc2[j]*yc2[j]); // radius L2 RP
236         Float_t zr0=(r2*zc1[i]-r1*zc2[j])/(r2-r1); //Z @ null radius
237         fZCombf->Fill(zr0);
238         fZCombc->Fill(zr0);
239         fZCombv->Fill(zr0);
240       }
241     }
242   }
243   delete [] xc1;
244   delete [] yc1;
245   delete [] zc1;
246   delete [] phi1;
247   delete [] xc2;
248   delete [] yc2;
249   delete [] zc2;
250   delete [] phi2;
251   if(fZCombc->GetEntries()==0){
252     Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
253     ResetHistograms();
254     return fCurrentVertex;
255   }
256   /*
257   else {
258     if(fDebug>0)cout<<"Number of entries in hist. "<<fZCombc->GetEntries()<<endl;
259   }
260   */
261
262   TH1F *hc = fZCombc;
263   Bool_t goon = kFALSE;
264   Int_t cnt = 0;
265   Int_t bi;
266
267   do {
268     goon = kFALSE;
269     cnt++;
270     bi = hc->GetMaximumBin();   // bin with maximal content on coarse histogram
271     if(hc->GetBinContent(bi)<3){
272       if(cnt==1)goon = kTRUE;
273       hc = fZCombv;
274     }
275   } while(goon);
276
277
278   Float_t centre = hc->GetBinCenter(bi);  // z value of the bin with maximal content
279   
280   // evaluation of the centroid
281   Int_t ii1=TMath::Max(bi-3,1);
282   Int_t ii2=TMath::Min(bi+3,hc->GetNbinsX());
283   centre = 0.;
284   Int_t nn=0;
285   for(Int_t ii=ii1;ii<ii2;ii++){
286     centre+= hc->GetBinCenter(ii)*hc->GetBinContent(ii);
287     nn+=static_cast<Int_t>(hc->GetBinContent(ii));
288   }
289   // PH Sometimes nn is 0, so we need a protection
290   if (nn) centre/=nn;
291   /*
292   if(fDebug>0){
293     cout<<"Value of center "<<centre<<endl;
294     cout<<"Population of 3 central bins: "<<hc->GetBinContent(bi-1)<<", ";
295     cout<<hc->GetBinContent(bi)<<", ";
296     cout<<hc->GetBinContent(bi+1)<<endl;
297   }
298   */
299   // n1 is the bin number of fine histogram containing the point located 1 coarse bin less than "centre"
300   Int_t n1 = static_cast<Int_t>((centre-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
301   // n2 is the bin number of fine histogram containing the point located 1 coarse bin more than "centre"
302   Int_t n2 = static_cast<Int_t>((centre+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
303   if(n1<1)n1=1;
304   if(n2>nbinfine)n2=nbinfine;
305   Int_t niter = 0;
306   goon = kTRUE;
307   Int_t num=0;
308   Bool_t last = kFALSE;
309
310   while(goon || last){
311     fZFound = 0.;
312     fZsig = 0.;
313     num=0;
314     // at the end of the loop:
315     // fZFound = N*(Average Z) - where N is the number of tracklets
316     // num=N
317     // fZsig = N*Q - where Q is the average of Z*Z
318     for(Int_t n=n1;n<=n2;n++){
319       fZFound+=fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
320       num+=static_cast<Int_t>(fZCombf->GetBinContent(n));
321       fZsig+=fZCombf->GetBinCenter(n)*fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
322     }
323     if(num<2){
324       fZsig = 0.;
325     }
326     else{
327       Float_t radi =  fZsig/(num-1)-fZFound*fZFound/num/(num-1);
328       // radi = square root of sample variance of Z
329       if(radi>0.)fZsig=TMath::Sqrt(radi);
330       else fZsig=0.;
331       // fZfound - Average Z
332       fZFound/=num;
333     }
334     if(!last){
335       // goon is true if the distance between the found Z and the lower bin differs from the distance between the found Z and
336       // the upper bin by more than tolerance (0.002)
337       goon = TMath::Abs(TMath::Abs(fZFound-fZCombf->GetBinCenter(n1))-TMath::Abs(fZFound-fZCombf->GetBinCenter(n2)))>fTolerance;
338       // a window in the fine grained histogram is centered aroung the found Z. The width is 2 bins of
339       // the coarse grained histogram
340       n1 = static_cast<Int_t>((fZFound-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
341       if(n1<1)n1=1;
342       n2 = static_cast<Int_t>((fZFound+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
343       if(n2>nbinfine)n2=nbinfine;
344       /*
345       if(fDebug>0){
346         cout<<"Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
347         cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
348       }
349       */
350       niter++;
351       // no more than 10 adjusting iterations
352       if(niter>=10){
353         goon = kFALSE;
354         Warning("FindVertexForCurrentEvent","The procedure dows not converge\n");
355       }
356
357       if((fZsig> 0.0001) && !goon && num>5){
358         last = kTRUE;
359         n1 = static_cast<Int_t>((fZFound-fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
360         if(n1<1)n1=1;
361         n2 = static_cast<Int_t>((fZFound+fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
362         if(n2>nbinfine)n2=nbinfine;
363         /*
364         if(fDebug>0){
365           cout<<"FINAL Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
366           cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
367         }
368         */
369       }
370     }
371     else {
372       last = kFALSE;
373     }
374
375   }
376   //  if(fDebug>0)cout<<"Numer of Iterations "<<niter<<endl<<endl;
377   if(num!=0)fZsig/=TMath::Sqrt(num);
378   fCurrentVertex = new AliESDVertex(fZFound,fZsig,num);
379   fCurrentVertex->SetTitle("vertexer: B");
380   ResetHistograms();
381   return fCurrentVertex;
382 }
383
384 //_____________________________________________________________________
385 void AliITSVertexerZ::ResetHistograms(){
386   // delete TH1 data members
387   if(fZCombc)delete fZCombc;
388   if(fZCombf)delete fZCombf;
389   if(fZCombv)delete fZCombv;
390   fZCombc = 0;
391   fZCombf = 0;
392   fZCombv = 0;
393 }
394
395 //______________________________________________________________________
396 void AliITSVertexerZ::FindVertices(){
397   // computes the vertices of the events in the range FirstEvent - LastEvent
398   AliRunLoader *rl = AliRunLoader::GetRunLoader();
399   AliITSLoader* itsLoader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
400   itsLoader->ReloadRecPoints();
401   for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
402     cout<<"Processing event "<<i<<endl;
403     rl->GetEvent(i);
404     FindVertexForCurrentEvent(i);
405     if(fCurrentVertex){
406       WriteCurrentVertex();
407     }
408     /*
409     else {
410       if(fDebug>0){
411         cout<<"Vertex not found for event "<<i<<endl;
412         cout<<"fZFound = "<<fZFound<<", fZsig= "<<fZsig<<endl;
413       }
414     }
415     */
416   }
417 }
418
419 //________________________________________________________
420 void AliITSVertexerZ::PrintStatus() const {
421   // Print current status
422   cout <<"=======================================================\n";
423   cout <<" First layer first and last modules: "<<fFirstL1<<", ";
424   cout <<fLastL1<<endl;
425   cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
426   cout <<fLastL2<<endl;
427   cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
428   cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
429   cout <<"Bin sizes for coarse and fine z histos "<<fStepCoarse<<"; "<<fStepFine<<endl;
430   cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
431   cout <<" Debug flag: "<<fDebug<<endl;
432   cout <<"First event to be processed "<<fFirstEvent;
433   cout <<"\n Last event to be processed "<<fLastEvent<<endl;
434   if(fZCombc){
435     cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
436   }
437   else{
438     cout<<"fZCombc does not exist\n";
439   }
440   if(fZCombv){
441     cout<<"fZCombv exists - entries="<<fZCombv->GetEntries()<<endl;
442   }
443   else{
444     cout<<"fZCombv does not exist\n";
445   }
446    if(fZCombf){
447     cout<<"fZCombf exists - entries="<<fZCombv->GetEntries()<<endl;
448   }
449   else{
450     cout<<"fZCombf does not exist\n";
451   }
452  
453   cout <<"=======================================================\n";
454 }
455