]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerZ.cxx
Using AliDebug (Massimo)
[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<TBranch.h>
17 #include<TClonesArray.h>
18 #include<TFile.h>
19 #include<TH1.h>
20 #include <TString.h>
21 #include<TTree.h>
22 #include "AliITSLoader.h"
23 #include "AliITSgeom.h"
24 #include "AliITSDetTypeRec.h"
25 #include "AliITSRecPoint.h"
26 #include "AliITSZPoint.h"
27 #include "AliHeader.h"
28 #include "AliGenEventHeader.h"
29
30 /////////////////////////////////////////////////////////////////
31 // this class implements a fast method to determine
32 // the Z coordinate of the primary vertex
33 // for p-p collisions it seems to give comparable or better results
34 // with respect to what obtained with AliITSVertexerPPZ
35 // It can be used successfully with Pb-Pb collisions
36 ////////////////////////////////////////////////////////////////
37
38 ClassImp(AliITSVertexerZ)
39
40
41
42 //______________________________________________________________________
43 AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer(),
44 fFirstL1(0),
45 fLastL1(0),
46 fFirstL2(0),
47 fLastL2(0),
48 fDiffPhiMax(0),
49 fZFound(0),
50 fZsig(0.),
51 fZCombc(0),
52 fLowLim(0.),
53 fHighLim(0.),
54 fStepCoarse(0),
55 fTolerance(0.),
56 fMaxIter(0),
57 fWindowWidth(0) {
58   // Default constructor
59   SetDiffPhiMax();
60   SetFirstLayerModules();
61   SetSecondLayerModules();
62   SetLowLimit();
63   SetHighLimit();
64   SetBinWidthCoarse();
65   SetTolerance();
66   SetPPsetting();
67   ConfigIterations();
68   SetWindowWidth();
69 }
70
71 //______________________________________________________________________
72 AliITSVertexerZ::AliITSVertexerZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn),
73 fFirstL1(0),
74 fLastL1(0),
75 fFirstL2(0),
76 fLastL2(0),
77 fDiffPhiMax(0),
78 fZFound(0),
79 fZsig(0.),
80 fZCombc(0),
81 fLowLim(0.),
82 fHighLim(0.),
83 fStepCoarse(0),
84 fTolerance(0.),
85 fMaxIter(0),
86 fWindowWidth(0) {
87   // Standard Constructor
88   SetDiffPhiMax();
89   SetFirstLayerModules();
90   SetSecondLayerModules();
91   SetLowLimit();
92   SetHighLimit();
93   SetBinWidthCoarse();
94   SetTolerance();
95   SetPPsetting();
96   ConfigIterations();
97   SetWindowWidth();
98
99 }
100
101 //______________________________________________________________________
102 AliITSVertexerZ::AliITSVertexerZ(const AliITSVertexerZ &vtxr) : AliITSVertexer(vtxr),
103 fFirstL1(vtxr.fFirstL1),
104 fLastL1(vtxr.fLastL1),
105 fFirstL2(vtxr.fFirstL2),
106 fLastL2(vtxr.fLastL2),
107 fDiffPhiMax(vtxr.fDiffPhiMax),
108 fZFound(vtxr.fZFound),
109 fZsig(vtxr.fZsig),
110 fZCombc(vtxr.fZCombc),
111 fLowLim(vtxr.fLowLim),
112 fHighLim(vtxr.fHighLim),
113 fStepCoarse(vtxr.fStepCoarse),
114 fTolerance(vtxr.fTolerance),
115 fMaxIter(vtxr.fMaxIter),
116 fWindowWidth(vtxr.fWindowWidth){
117   // Copy constructor
118
119 }
120
121 //______________________________________________________________________
122 AliITSVertexerZ& AliITSVertexerZ::operator=(const AliITSVertexerZ&  vtxr ){
123   // Assignment operator
124
125   this->~AliITSVertexerZ();
126   new(this) AliITSVertexerZ(vtxr);
127   return *this;
128 }
129
130 //______________________________________________________________________
131 AliITSVertexerZ::~AliITSVertexerZ() {
132   // Destructor
133   delete fZCombc;
134 }
135
136 //______________________________________________________________________
137 void AliITSVertexerZ::ConfigIterations(Int_t noiter,Float_t *ptr){
138   // configure the iterative procedure to gain efficiency for
139   // pp events with very low multiplicity
140   Float_t defaults[5]={0.05,0.1,0.2,0.3,0.5};
141   fMaxIter=noiter;
142   if(noiter>5){
143     Error("ConfigIterations","Maximum number of iterations is 5\n");
144     fMaxIter=5;
145   }
146   for(Int_t j=0;j<5;j++)fPhiDiffIter[j]=defaults[j];
147   if(ptr)for(Int_t j=0;j<fMaxIter;j++)fPhiDiffIter[j]=ptr[j];
148 }
149
150 //______________________________________________________________________
151 Int_t AliITSVertexerZ::GetPeakRegion(TH1F*h, Int_t &binmin, Int_t &binmax) const {
152   // Finds a region around a peak in the Z histogram
153   // Case of 2 peaks is treated 
154   Int_t imax=h->GetNbinsX();
155   Float_t maxval=0;
156   Int_t bi1=h->GetMaximumBin();
157   Int_t bi2=0;
158   for(Int_t i=imax;i>=1;i--){
159     if(h->GetBinContent(i)>maxval){
160       maxval=h->GetBinContent(i);
161       bi2=i;
162     }
163   }
164   Int_t npeaks=0;
165
166   if(bi1==bi2){
167     binmin=bi1-3;
168     binmax=bi1+3;
169     npeaks=1;
170   }else{
171     TH1F *copy = new TH1F(*h);
172     copy->SetBinContent(bi1,0.);
173     copy->SetBinContent(bi2,0.);
174     Int_t l1=TMath::Max(bi1-3,1);
175     Int_t l2=TMath::Min(bi1+3,h->GetNbinsX());
176     Float_t cont1=copy->Integral(l1,l2);
177     Int_t ll1=TMath::Max(bi2-3,1);
178     Int_t ll2=TMath::Min(bi2+3,h->GetNbinsX());
179     Float_t cont2=copy->Integral(ll1,ll2);
180     if(cont1>cont2){
181       binmin=l1;
182       binmax=l2;
183       npeaks=1;
184     }
185     if(cont2>cont1){
186       binmin=ll1;
187       binmax=ll2;
188       npeaks=1;
189     }
190     if(cont1==cont2){
191       binmin=l1;
192       binmax=ll2;
193       if(bi2-bi1==1) npeaks=1;
194       else npeaks=2;
195     }  
196     delete copy;
197   }    
198   return npeaks;
199 }
200 //______________________________________________________________________
201 AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(Int_t evnumber){
202   // Defines the AliESDVertex for the current event
203   VertexZFinder(evnumber);
204   Int_t ntrackl=0;
205   for(Int_t iteraz=0;iteraz<fMaxIter;iteraz++){
206     if(fCurrentVertex) ntrackl=fCurrentVertex->GetNContributors();
207     if(!fCurrentVertex || ntrackl==0 || ntrackl==-1){
208       Float_t diffPhiMaxOrig=fDiffPhiMax;
209       fDiffPhiMax=GetPhiMaxIter(iteraz);
210       VertexZFinder(evnumber);
211       fDiffPhiMax=diffPhiMaxOrig;
212     }
213   }
214   FindMultiplicity(evnumber);
215   return fCurrentVertex;
216 }  
217
218
219
220
221 //______________________________________________________________________
222 void AliITSVertexerZ::VertexZFinder(Int_t evnumber){
223   // Defines the AliESDVertex for the current event
224   fCurrentVertex = 0;
225   AliRunLoader *rl =AliRunLoader::GetRunLoader();
226   AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
227   AliITSgeom* geom = itsLoader->GetITSgeom();
228   itsLoader->LoadRecPoints();
229   rl->GetEvent(evnumber);
230
231   AliITSDetTypeRec detTypeRec;
232
233   TTree *tR = itsLoader->TreeR();
234   detTypeRec.SetTreeAddressR(tR);
235   TClonesArray *itsRec  = 0;
236   // lc and gc are local and global coordinates for layer 1
237   Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
238   Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
239   // lc2 and gc2 are local and global coordinates for layer 2
240   Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
241   Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
242
243   itsRec = detTypeRec.RecPoints();
244   TBranch *branch;
245   branch = tR->GetBranch("ITSRecPoints");
246
247   Int_t nrpL1 = 0;
248   Int_t nrpL2 = 0;
249   // By default fFirstL1=0 and fLastL1=79
250   // This loop counts the number of recpoints on layer1 (central modules)
251   for(Int_t module= fFirstL1; module<=fLastL1;module++){
252     // Keep only central modules
253     //    if(module%4==0 || module%4==3)continue;
254     //   cout<<"Procesing module "<<module<<" ";
255     branch->GetEvent(module);
256     //    cout<<"Number of clusters "<<clusters->GetEntries()<<endl;
257     nrpL1+= itsRec->GetEntries();
258     detTypeRec.ResetRecPoints();
259   }
260   //By default fFirstL2=80 and fLastL2=239
261   //This loop counts the number of RP on layer 2
262   for(Int_t module= fFirstL2; module<=fLastL2;module++){
263     branch->GetEvent(module);
264     nrpL2+= itsRec->GetEntries();
265     detTypeRec.ResetRecPoints();
266   }
267   if(nrpL1 == 0 || nrpL2 == 0){
268     ResetHistograms();
269     itsLoader->UnloadRecPoints();
270     fCurrentVertex = new AliESDVertex(0.,5.3,-2);
271     return;
272   }
273   // The vertex finding is attempted only if the number of RP is !=0 on
274   // both layers
275   Float_t *xc1 = new Float_t [nrpL1]; // coordinates of the L1 Recpoints
276   Float_t *yc1 = new Float_t [nrpL1];
277   Float_t *zc1 = new Float_t [nrpL1];
278   Float_t *phi1 = new Float_t [nrpL1];
279   Float_t *err1 = new Float_t [nrpL1];
280   Float_t *xc2 = new Float_t [nrpL2]; // coordinates of the L1 Recpoints
281   Float_t *yc2 = new Float_t [nrpL2];
282   Float_t *zc2 = new Float_t [nrpL2];
283   Float_t *phi2 = new Float_t [nrpL2];
284   Float_t *err2 = new Float_t [nrpL2];
285   Int_t ind = 0;// running index for RP
286   // Force a coarse bin size of 200 microns if the number of clusters on layer 2
287   // is low
288   if(nrpL2<fPPsetting[0])SetBinWidthCoarse(fPPsetting[1]);
289   // By default nbincoarse=(10+10)/0.01=2000
290   Int_t nbincoarse = static_cast<Int_t>((fHighLim-fLowLim)/fStepCoarse);
291   if(fZCombc)delete fZCombc;
292   fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
293
294   // Loop on modules of layer 1 
295
296   for(Int_t module= fFirstL1; module<=fLastL1;module++){
297     //    if(module%4==0 || module%4==3)continue;
298     branch->GetEvent(module);
299     Int_t nrecp1 = itsRec->GetEntries();
300     for(Int_t j=0;j<nrecp1;j++){
301       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
302       // Local coordinates of this recpoint
303       lc[0]=recp->GetDetLocalX();
304       lc[2]=recp->GetDetLocalZ();
305       geom->LtoG(module,lc,gc);
306       // Global coordinates of this recpoints
307       gc[0]-=fNominalPos[0]; // Possible beam offset in the bending plane
308       gc[1]-=fNominalPos[1]; //   "               "
309       xc1[ind]=gc[0];
310       yc1[ind]=gc[1];
311       zc1[ind]=gc[2];
312       // azimuthal angle is computed in the interval 0 --> 2*pi
313       phi1[ind] = TMath::ATan2(gc[1],gc[0]);
314       if(phi1[ind]<0)phi1[ind]=2*TMath::Pi()+phi1[ind];
315       err1[ind]=recp->GetSigmaZ2();
316       ind++;
317     }
318     detTypeRec.ResetRecPoints();
319   }
320   ind = 0; // the running index is reset for Layer 2
321   for(Int_t module= fFirstL2; module<=fLastL2;module++){
322     branch->GetEvent(module);
323     Int_t nrecp2 = itsRec->GetEntries();
324     for(Int_t j=0;j<nrecp2;j++){
325       AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
326       lc[0]=recp->GetDetLocalX();
327       lc[2]=recp->GetDetLocalZ();
328       geom->LtoG(module,lc,gc);
329       gc[0]-=fNominalPos[0];
330       gc[1]-=fNominalPos[1];
331       xc2[ind]=gc[0];
332       yc2[ind]=gc[1];
333       zc2[ind]=gc[2];
334       phi2[ind] = TMath::ATan2(gc[1],gc[0]);
335       if(phi2[ind]<0)phi2[ind]=2*TMath::Pi()+phi2[ind];
336       err2[ind]=recp->GetSigmaZ2();
337       ind++;
338     }
339     detTypeRec.ResetRecPoints();
340   }
341  
342 /* Test the ffect of mutiple scatternig on error. Negligible
343   // Multiple scattering
344   Float_t beta=1.,pmed=0.875; //pmed=875 MeV (for tracks with dphi<0.01 rad)
345   Float_t beta2=beta*beta;
346   Float_t p2=pmed*pmed;
347   Float_t rBP=3; //Beam Pipe radius = 3cm
348   Float_t dBP=0.08/35.3; // 800 um of Be
349   Float_t dL1=0.01; //approx. 1% of radiation length  
350   Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
351   Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
352 */
353   TClonesArray *points = new TClonesArray("AliITSZPoint",nrpL1*nrpL2);
354   TClonesArray &pts = *points;
355   Int_t nopoints =0;
356   for(Int_t i=0;i<nrpL1;i++){ // loop on L1 RP
357     Float_t r1=TMath::Sqrt(xc1[i]*xc1[i]+yc1[i]*yc1[i]); // radius L1 RP
358     for(Int_t j=0;j<nrpL2;j++){ // loop on L2 RP
359       Float_t diff = TMath::Abs(phi2[j]-phi1[i]); // diff in azimuth
360       if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; //diff<pi
361       if(diff<fDiffPhiMax){ // cut on 10 milliradians by def.
362         Float_t r2=TMath::Sqrt(xc2[j]*xc2[j]+yc2[j]*yc2[j]); // radius L2 RP
363 //      Float_t tgth=(zc2[j]-zc1[i])/(r2-r1); // slope
364         Float_t zr0=(r2*zc1[i]-r1*zc2[j])/(r2-r1); //Z @ null radius
365         Float_t ezr0q=(r2*r2*err1[i]+r1*r1*err2[j])/(r2-r1)/(r2-r1); //error on Z @ null radius
366         /*
367         ezr0q+=r1*r1*(1+tgth*tgth)*theta2L1/2; // multiple scattering in layer 1
368         ezr0q+=rBP*rBP*(1+tgth*tgth)*theta2BP/2; // multiple scattering in beam pipe
369         */
370         new(pts[nopoints++])AliITSZPoint(zr0,ezr0q);
371
372         fZCombc->Fill(zr0);
373       }
374     }
375   }
376   delete [] xc1;
377   delete [] yc1;
378   delete [] zc1;
379   delete [] phi1;
380   delete [] err1;
381   delete [] xc2;
382   delete [] yc2;
383   delete [] zc2;
384   delete [] phi2;
385   delete [] err2;
386
387   points->Sort();
388
389   Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
390   if(contents<1.){
391     //    Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
392     ResetHistograms();
393     itsLoader->UnloadRecPoints();
394     fCurrentVertex = new AliESDVertex(0.,5.3,-1);
395     return;
396   }
397
398   TH1F *hc = fZCombc;
399
400   
401   if(hc->GetBinContent(hc->GetMaximumBin())<3)hc->Rebin(3);
402   Int_t binmin,binmax;
403   Int_t nPeaks=GetPeakRegion(hc,binmin,binmax);   
404   if(nPeaks==2)AliWarning("2 peaks found");
405   Float_t zm =0.;
406   Float_t ezm =0.;
407   Float_t lim1 = hc->GetBinLowEdge(binmin);
408   Float_t lim2 = hc->GetBinLowEdge(binmax)+hc->GetBinWidth(binmax);
409
410   if(nPeaks ==1 && (lim2-lim1)<fWindowWidth){
411     Float_t c=(lim1+lim2)/2.;
412     lim1=c-fWindowWidth/2.;
413     lim2=c+fWindowWidth/2.;
414   }
415   Int_t niter = 0, ncontr=0;
416   do {
417     // symmetrization
418     if(zm  !=0.){
419       Float_t semilarg=TMath::Min((lim2-zm),(zm-lim1));
420       lim1=zm - semilarg;
421       lim2=zm + semilarg;
422     }
423
424     zm=0.;
425     ezm=0.;
426     ncontr=0;
427     for(Int_t i =0; i<points->GetEntriesFast(); i++){
428       AliITSZPoint* p=(AliITSZPoint*)points->UncheckedAt(i);
429       if(p->GetZ()>lim1 && p->GetZ()<lim2){
430         Float_t deno = p->GetErrZ();
431         zm+=p->GetZ()/deno;
432         ezm+=1./deno;
433         ncontr++;
434       }
435     }
436     zm/=ezm;
437     ezm=TMath::Sqrt(1./ezm);
438     niter++;
439   } while(niter<10 && TMath::Abs((zm-lim1)-(lim2-zm))>fTolerance);
440   fCurrentVertex = new AliESDVertex(zm,ezm,ncontr);
441   fCurrentVertex->SetTitle("vertexer: B");
442   ResetHistograms();
443   itsLoader->UnloadRecPoints();
444   return;
445 }
446
447 //_____________________________________________________________________
448 void AliITSVertexerZ::ResetHistograms(){
449   // delete TH1 data members
450   if(fZCombc)delete fZCombc;
451   fZCombc = 0;
452 }
453
454 //______________________________________________________________________
455 void AliITSVertexerZ::FindVertices(){
456   // computes the vertices of the events in the range FirstEvent - LastEvent
457   AliRunLoader *rl = AliRunLoader::GetRunLoader();
458   AliITSLoader* itsLoader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
459   itsLoader->ReloadRecPoints();
460   for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
461     //  cout<<"Processing event "<<i<<endl;
462     rl->GetEvent(i);
463     FindVertexForCurrentEvent(i);
464     if(fCurrentVertex){
465       WriteCurrentVertex();
466     }
467   }
468 }
469
470 //________________________________________________________
471 void AliITSVertexerZ::PrintStatus() const {
472   // Print current status
473   cout <<"=======================================================\n";
474   cout <<" First layer first and last modules: "<<fFirstL1<<", ";
475   cout <<fLastL1<<endl;
476   cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
477   cout <<fLastL2<<endl;
478   cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
479   cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
480   cout <<"Bin sizes for coarse z histos "<<fStepCoarse<<endl;
481   cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
482   cout <<"First event to be processed "<<fFirstEvent;
483   cout <<"\n Last event to be processed "<<fLastEvent<<endl;
484   if(fZCombc){
485     cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
486   }
487   else{
488     cout<<"fZCombc does not exist\n";
489   }
490  
491   cout <<"=======================================================\n";
492 }
493