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