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