]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerZ.cxx
Fixing bug #59311
[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 "AliITSRecPoint.h"
24 #include "AliITSRecPointContainer.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   TClonesArray *itsRec  = 0;
193   // lc1 and gc1 are local and global coordinates for layer 1
194   Float_t lc1[3]; for(Int_t ii=0; ii<3; ii++) lc1[ii]=0.;
195   Float_t gc1[3]; for(Int_t ii=0; ii<3; ii++) gc1[ii]=0.;
196   // lc2 and gc2 are local and global coordinates for layer 2
197   Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
198   Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
199   AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
200   itsRec=rpcont->FetchClusters(0,itsClusterTree);
201   if(!rpcont->IsSPDActive()){
202     AliWarning("Null pointer for RecPoints branch, vertex not calculated");
203     ResetHistograms();
204     fCurrentVertex = new AliESDVertex(0.,5.3,-2);
205     return;
206   }
207
208   Int_t nrpL1 = 0;
209   Int_t nrpL2 = 0;
210
211   // By default fFirstL1=0 and fLastL1=79
212   for(Int_t module= fFirstL1; module<=fLastL1;module++){
213     itsRec=rpcont->UncheckedGetClusters(module);
214     nrpL1+= itsRec->GetEntries();
215   }
216   //By default fFirstL2=80 and fLastL2=239
217   for(Int_t module= fFirstL2; module<=fLastL2;module++){
218     itsRec=rpcont->UncheckedGetClusters(module);
219     nrpL2+= itsRec->GetEntries();
220   }
221   if(nrpL1 == 0 || nrpL2 == 0){
222     AliDebug(1,Form("No RecPoints in at least one SPD layer (%d %d)",nrpL1,nrpL2));
223     ResetHistograms();
224     fCurrentVertex = new AliESDVertex(0.,5.3,-2);
225     return;
226   }
227   // Force a coarse bin size of 200 microns if the number of clusters on layer 2
228   // is low
229   if(nrpL2<fPPsetting[0])SetBinWidthCoarse(fPPsetting[1]);
230   // By default nbincoarse=(10+10)/0.01=2000
231   Int_t nbincoarse = static_cast<Int_t>((fHighLim-fLowLim)/fStepCoarse);
232   if(fZCombc)delete fZCombc;
233   fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
234
235  /* Test the ffect of mutiple scatternig on error. Negligible
236   // Multiple scattering
237   Float_t beta=1.,pmed=0.875; //pmed=875 MeV (for tracks with dphi<0.01 rad)
238   Float_t beta2=beta*beta;
239   Float_t p2=pmed*pmed;
240   Float_t rBP=3; //Beam Pipe radius = 3cm
241   Float_t dBP=0.08/35.3; // 800 um of Be
242   Float_t dL1=0.01; //approx. 1% of radiation length  
243   Float_t theta2BP=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dBP);
244   Float_t theta2L1=14.1*14.1/(beta2*p2*1e6)*TMath::Abs(dL1);
245 */
246   Int_t maxdim=TMath::Min(nrpL1*nrpL2,50000);  // temporary; to limit the size in PbPb
247   static TClonesArray points("AliITSZPoint",maxdim);
248   Int_t nopoints =0;
249   for(Int_t modul1= fFirstL1; modul1<=fLastL1;modul1++){   // Loop on modules of layer 1
250     if(!fUseModule[modul1]) continue;
251     UShort_t ladder=int(modul1/4)+1;  // ladders are numbered starting from 1
252     TClonesArray *prpl1=rpcont->UncheckedGetClusters(modul1);
253     Int_t nrecp1 = prpl1->GetEntries();
254     for(Int_t j1=0;j1<nrecp1;j1++){
255       AliITSRecPoint *recp = (AliITSRecPoint*)prpl1->At(j1);
256       recp->GetGlobalXYZ(gc1);
257       gc1[0]-=GetNominalPos()[0]; // Possible beam offset in the bending plane
258       gc1[1]-=GetNominalPos()[1]; //   "               "
259       Float_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
260       Float_t phi1 = TMath::ATan2(gc1[1],gc1[0]);
261       if(phi1<0)phi1+=2*TMath::Pi();
262       Float_t zc1=gc1[2];
263       Float_t erz1=recp->GetSigmaZ2();
264       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
265         for(Int_t k=0;k<4;k++){
266           Int_t ladmod=fLadders[ladder-1]+ladl2;
267           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
268           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
269           if(!fUseModule[modul2]) continue;
270           itsRec=rpcont->UncheckedGetClusters(modul2);
271           Int_t nrecp2 = itsRec->GetEntries();
272           for(Int_t j2=0;j2<nrecp2;j2++){
273             recp = (AliITSRecPoint*)itsRec->At(j2);
274             recp->GetGlobalXYZ(gc2);
275             gc2[0]-=GetNominalPos()[0];
276             gc2[1]-=GetNominalPos()[1];
277             Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
278             Float_t phi2 = TMath::ATan2(gc2[1],gc2[0]);
279             if(phi2<0)phi2+=2*TMath::Pi();
280             Float_t zc2=gc2[2];
281             Float_t erz2=recp->GetSigmaZ2();
282
283             Float_t diff = TMath::Abs(phi2-phi1); 
284             if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
285             if(diff<fDiffPhiMax){
286               //        Float_t tgth=(zc2[j]-zc1[i])/(r2-r1); // slope (used for multiple scattering)
287               Float_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
288               Float_t ezr0q=(r2*r2*erz1+r1*r1*erz2)/(r2-r1)/(r2-r1); //error on Z @ null radius
289          /*
290          // Multiple scattering
291         ezr0q+=r1*r1*(1+tgth*tgth)*theta2L1/2; // multiple scattering in layer 1
292         ezr0q+=rBP*rBP*(1+tgth*tgth)*theta2BP/2; // multiple scattering in beam pipe
293         */
294               if(nopoints<maxdim) new(points[nopoints++])AliITSZPoint(zr0,ezr0q);             
295               fZCombc->Fill(zr0);
296             }
297           }
298         }
299       }
300     }
301   }
302
303   points.Sort();
304
305   Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
306   if(contents<1.){
307     //    Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
308     ResetHistograms();
309     fCurrentVertex = new AliESDVertex(0.,5.3,-1);
310     points.Clear();
311     return;
312   }
313
314   TH1F *hc = fZCombc;
315
316   
317   if(hc->GetBinContent(hc->GetMaximumBin())<3)hc->Rebin(4);
318   Int_t binmin,binmax;
319   Int_t nPeaks=GetPeakRegion(hc,binmin,binmax);   
320   if(nPeaks==2)AliWarning("2 peaks found");
321   Float_t zm =0.;
322   Float_t ezm =0.;
323   Float_t lim1 = hc->GetBinLowEdge(binmin);
324   Float_t lim2 = hc->GetBinLowEdge(binmax)+hc->GetBinWidth(binmax);
325
326   if(nPeaks ==1 && (lim2-lim1)<fWindowWidth){
327     Float_t c=(lim1+lim2)/2.;
328     lim1=c-fWindowWidth/2.;
329     lim2=c+fWindowWidth/2.;
330   }
331   Int_t niter = 0, ncontr=0;
332   do {
333     // symmetrization
334     if(zm  !=0.){
335       Float_t semilarg=TMath::Min((lim2-zm),(zm-lim1));
336       lim1=zm - semilarg;
337       lim2=zm + semilarg;
338     }
339
340     zm=0.;
341     ezm=0.;
342     ncontr=0;
343     for(Int_t i =0; i<points.GetEntries(); i++){
344       AliITSZPoint* p=(AliITSZPoint*)points.UncheckedAt(i);
345       if(p->GetZ()>lim1 && p->GetZ()<lim2){
346         Float_t deno = p->GetErrZ();
347         zm+=p->GetZ()/deno;
348         ezm+=1./deno;
349         ncontr++;
350       }
351     }
352     if(ezm>0) {
353       zm/=ezm;
354       ezm=TMath::Sqrt(1./ezm);
355     }
356     niter++;
357   } while(niter<10 && TMath::Abs((zm-lim1)-(lim2-zm))>fTolerance);
358   fCurrentVertex = new AliESDVertex(zm,ezm,ncontr);
359   fCurrentVertex->SetTitle("vertexer: Z");
360   fNoVertices=1;
361   points.Clear();
362   if(ncontr>fMinTrackletsForPilup){ 
363     Float_t secPeakPos;
364     Int_t ncontr2=FindSecondPeak(fZCombc,binmin,binmax,secPeakPos);
365     if(ncontr2>=fMinTrackletsForPilup){ 
366       fIsPileup=kTRUE;
367       fNoVertices=2;
368       fZpuv=secPeakPos;
369       fNTrpuv=ncontr2;
370       AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
371       fVertArray = new AliESDVertex[2];
372       fVertArray[0]=(*fCurrentVertex);
373       fVertArray[1]=secondVert;
374     }
375   }
376   if(fNoVertices==1){
377     fVertArray = new AliESDVertex[1];
378     fVertArray[0]=(*fCurrentVertex);      
379   }
380   
381   ResetHistograms();
382   return;
383 }
384
385 //_____________________________________________________________________
386 Int_t AliITSVertexerZ::FindSecondPeak(TH1F* h, Int_t binmin,Int_t binmax, Float_t& secPeakPos){  
387   for(Int_t i=binmin-1;i<=binmax+1;i++){
388     h->SetBinContent(i,0.);
389   }
390   Int_t secPeakBin=h->GetMaximumBin();
391   secPeakPos=h->GetBinCenter(secPeakBin);
392   Int_t secPeakCont=(Int_t)h->GetBinContent(secPeakBin);
393   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin-1);
394   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin+1);  
395   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin-2);
396   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin+2);  
397   return secPeakCont;
398 }
399
400 //_____________________________________________________________________
401 void AliITSVertexerZ::ResetHistograms(){
402   // delete TH1 data members
403   if(fZCombc)delete fZCombc;
404   fZCombc = 0;
405 }
406
407 //________________________________________________________
408 void AliITSVertexerZ::PrintStatus() const {
409   // Print current status
410   cout <<"=======================================================\n";
411   cout <<" First layer first and last modules: "<<fFirstL1<<", ";
412   cout <<fLastL1<<endl;
413   cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
414   cout <<fLastL2<<endl;
415   cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
416   cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
417   cout <<"Bin sizes for coarse z histos "<<fStepCoarse<<endl;
418   cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
419   if(fZCombc){
420     cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
421   }
422   else{
423     cout<<"fZCombc does not exist\n";
424   }
425  
426   cout <<"=======================================================\n";
427 }
428