]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerZ.cxx
Extend q_ij range for pp and pPb
[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 using std::endl;
36 using std::cout;
37 ClassImp(AliITSVertexerZ)
38
39
40
41 //______________________________________________________________________
42 AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer(),
43 fFirstL1(0),
44 fLastL1(0),
45 fFirstL2(0),
46 fLastL2(0),
47 fDiffPhiMax(0),
48 fZFound(0),
49 fZsig(0.),
50 fZCombc(0),
51 fLowLim(0.),
52 fHighLim(0.),
53 fStepCoarse(0),
54 fTolerance(0.),
55 fMaxIter(0),
56 fWindowWidth(0),
57 fSearchForPileup(kTRUE)
58 {
59   // Default constructor
60   SetDiffPhiMax();
61   SetFirstLayerModules();
62   SetSecondLayerModules();
63   SetLowLimit();
64   SetHighLimit();
65   SetBinWidthCoarse();
66   SetTolerance();
67   SetPPsetting();
68   ConfigIterations();
69   SetWindowWidth();
70 }
71
72 //______________________________________________________________________
73 AliITSVertexerZ::AliITSVertexerZ(Float_t x0, Float_t y0):AliITSVertexer(),
74 fFirstL1(0),
75 fLastL1(0),
76 fFirstL2(0),
77 fLastL2(0),
78 fDiffPhiMax(0),
79 fZFound(0),
80 fZsig(0.),
81 fZCombc(0),
82 fLowLim(0.),
83 fHighLim(0.),
84 fStepCoarse(0),
85 fTolerance(0.),
86 fMaxIter(0),
87 fWindowWidth(0),
88 fSearchForPileup(kTRUE)
89 {
90   // Standard Constructor
91   SetDiffPhiMax();
92   SetFirstLayerModules();
93   SetSecondLayerModules();
94   SetLowLimit();
95   SetHighLimit();
96   SetBinWidthCoarse();
97   SetTolerance();
98   SetPPsetting();
99   ConfigIterations();
100   SetWindowWidth();
101   SetVtxStart((Double_t)x0,(Double_t)y0,0.);
102
103 }
104
105 //______________________________________________________________________
106 AliITSVertexerZ::~AliITSVertexerZ() {
107   // Destructor
108   delete fZCombc;
109 }
110
111 //______________________________________________________________________
112 void AliITSVertexerZ::ConfigIterations(Int_t noiter,Float_t *ptr){
113   // configure the iterative procedure to gain efficiency for
114   // pp events with very low multiplicity
115   Float_t defaults[5]={0.02,0.05,0.1,0.2,0.3};
116   fMaxIter=noiter;
117   if(noiter>5){
118     Error("ConfigIterations","Maximum number of iterations is 5\n");
119     fMaxIter=5;
120   }
121   for(Int_t j=0;j<5;j++)fPhiDiffIter[j]=defaults[j];
122   if(ptr)for(Int_t j=0;j<fMaxIter;j++)fPhiDiffIter[j]=ptr[j];
123 }
124
125 //______________________________________________________________________
126 Int_t AliITSVertexerZ::GetPeakRegion(TH1F*h, Int_t &binmin, Int_t &binmax){
127   // Finds a region around a peak in the Z histogram
128   // Case of 2 peaks is treated 
129   Int_t imax=h->GetNbinsX();
130   Float_t maxval=0;
131   Int_t bi1=h->GetMaximumBin();
132   Int_t bi2=0;
133   for(Int_t i=imax;i>=1;i--){
134     if(h->GetBinContent(i)>maxval){
135       maxval=h->GetBinContent(i);
136       bi2=i;
137     }
138   }
139   Int_t npeaks=0;
140
141   if(bi1==bi2){
142     binmin=bi1-3;
143     binmax=bi1+3;
144     npeaks=1;
145   }else{
146     TH1F *copy = new TH1F(*h);
147     copy->SetBinContent(bi1,0.);
148     copy->SetBinContent(bi2,0.);
149     Int_t l1=TMath::Max(bi1-3,1);
150     Int_t l2=TMath::Min(bi1+3,h->GetNbinsX());
151     Float_t cont1=copy->Integral(l1,l2);
152     Int_t ll1=TMath::Max(bi2-3,1);
153     Int_t ll2=TMath::Min(bi2+3,h->GetNbinsX());
154     Float_t cont2=copy->Integral(ll1,ll2);
155     if(cont1>cont2){
156       binmin=l1;
157       binmax=l2;
158       npeaks=1;
159     }
160     if(cont2>cont1){
161       binmin=ll1;
162       binmax=ll2;
163       npeaks=1;
164     }
165     if(cont1==cont2){
166       binmin=l1;
167       binmax=ll2;
168       if(bi2-bi1==1) npeaks=1;
169       else npeaks=2;
170     }  
171     delete copy;
172   }    
173   return npeaks;
174 }
175 //______________________________________________________________________
176 AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(TTree *itsClusterTree){
177   // Defines the AliESDVertex for the current event
178   VertexZFinder(itsClusterTree);
179   Int_t ntrackl=0;
180   for(Int_t iteraz=0;iteraz<fMaxIter;iteraz++){
181     if(fCurrentVertex) ntrackl=fCurrentVertex->GetNContributors();
182     if(!fCurrentVertex || ntrackl==0 || ntrackl==-1){
183       Float_t diffPhiMaxOrig=fDiffPhiMax;
184       fDiffPhiMax=GetPhiMaxIter(iteraz);
185       VertexZFinder(itsClusterTree);
186       fDiffPhiMax=diffPhiMaxOrig;
187     }
188   }
189   if(fComputeMultiplicity) FindMultiplicity(itsClusterTree);
190   return fCurrentVertex;
191 }  
192
193 //______________________________________________________________________
194 void AliITSVertexerZ::VertexZFinder(TTree *itsClusterTree){
195   // Defines the AliESDVertex for the current event
196   delete fCurrentVertex;
197   fCurrentVertex = 0;
198   Double_t startPos[3]={GetNominalPos()[0],GetNominalPos()[1],GetNominalPos()[2]};
199   Double_t startCov[6]={GetNominalCov()[0],GetNominalCov()[1],GetNominalCov()[2],
200                         GetNominalCov()[3],GetNominalCov()[4],GetNominalCov()[5]};
201   ResetVertex();
202   TClonesArray *itsRec  = 0;
203   // lc1 and gc1 are local and global coordinates for layer 1
204   Float_t gc1[3]={0.,0.,0.}; // ; for(Int_t ii=0; ii<3; ii++) gc1[ii]=0.;
205   // lc2 and gc2 are local and global coordinates for layer 2
206   Float_t gc2[3]={0.,0.,0.}; //; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
207   AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
208   rpcont->FetchClusters(0,itsClusterTree);
209   if(!rpcont->IsSPDActive()){
210     AliWarning("Null pointer for RecPoints branch, vertex not calculated");
211     ResetHistograms();
212     fCurrentVertex = new AliESDVertex(startPos,startCov,99999.,-2);
213     return;
214   }
215
216   Int_t nrpL1 = 0;
217   Int_t nrpL2 = 0;
218   nrpL1=rpcont->GetNClustersInLayerFast(1);
219   nrpL2=rpcont->GetNClustersInLayerFast(2);
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(startPos,startCov,99999.,-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 nEntriesMod[kNSPDMod];
247   TClonesArray* recpArr[kNSPDMod];
248   for(Int_t modul=0; modul<kNSPDMod; ++modul) {
249     if(!fUseModule[modul]) {
250       nEntriesMod[modul]=0;
251       recpArr[modul]=0;
252     } else {
253       recpArr[modul]=rpcont->UncheckedGetClusters(modul);
254       nEntriesMod[modul]=recpArr[modul]->GetEntriesFast();
255     }
256   }
257       
258   Int_t maxdim=TMath::Min(nrpL1*nrpL2,50000);  // temporary; to limit the size in PbPb
259   static TClonesArray points("AliITSZPoint",maxdim);
260   Int_t nopoints =0;
261   for(Int_t modul1= fFirstL1; modul1<=fLastL1;modul1++){   // Loop on modules of layer 1
262     if(!fUseModule[modul1]) continue;
263     UShort_t ladder=int(modul1/4)+1;  // ladders are numbered starting from 1
264     TClonesArray *prpl1=recpArr[modul1]; //rpcont->UncheckedGetClusters(modul1);
265     Int_t nrecp1 = nEntriesMod[modul1]; //prpl1->GetEntries();
266     for(Int_t j1=0;j1<nrecp1;j1++){
267       AliITSRecPoint *recp1 = (AliITSRecPoint*)prpl1->At(j1);
268       recp1->GetGlobalXYZ(gc1);
269       gc1[0]-=GetNominalPos()[0]; // Possible beam offset in the bending plane
270       gc1[1]-=GetNominalPos()[1]; //   "               "
271       Float_t phi1 = TMath::ATan2(gc1[1],gc1[0]);
272       if(phi1<0)phi1+=TMath::TwoPi();
273       for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
274         for(Int_t k=0;k<4;k++){
275           Int_t ladmod=fLadders[ladder-1]+ladl2;
276           if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
277           Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
278           if(modul2<0)continue;
279           if(!fUseModule[modul2]) continue;
280           itsRec=recpArr[modul2]; // rpcont->UncheckedGetClusters(modul2);
281           Int_t nrecp2 = nEntriesMod[modul2]; // itsRec->GetEntries();
282           for(Int_t j2=0;j2<nrecp2;j2++){
283             AliITSRecPoint *recp2 = (AliITSRecPoint*)itsRec->At(j2);
284             recp2->GetGlobalXYZ(gc2);
285             gc2[0]-=GetNominalPos()[0];
286             gc2[1]-=GetNominalPos()[1];
287             Float_t phi2 = TMath::ATan2(gc2[1],gc2[0]);
288             if(phi2<0)phi2+=TMath::TwoPi();
289
290             Float_t diff = TMath::Abs(phi2-phi1); 
291             if(diff>TMath::Pi())diff=TMath::TwoPi()-diff;
292             if(diff<fDiffPhiMax){
293               Float_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
294               Float_t zc1=gc1[2];
295               Float_t erz1=recp1->GetSigmaZ2();
296               Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
297               Float_t zc2=gc2[2];
298               Float_t erz2=recp2->GetSigmaZ2();
299               //        Float_t tgth=(zc2[j]-zc1[i])/(r2-r1); // slope (used for multiple scattering)
300               Float_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
301               Float_t ezr0q=(r2*r2*erz1+r1*r1*erz2)/((r2-r1)*(r2-r1)); //error on Z @ null radius
302               /*
303               // Multiple scattering
304               ezr0q+=r1*r1*(1+tgth*tgth)*theta2L1/2; // multiple scattering in layer 1
305               ezr0q+=rBP*rBP*(1+tgth*tgth)*theta2BP/2; // multiple scattering in beam pipe
306               */
307               if(nopoints<maxdim) new(points[nopoints++])AliITSZPoint(zr0,ezr0q);             
308               fZCombc->Fill(zr0);
309             }
310           }
311         }
312       }
313     }
314   }
315
316   points.Sort();
317
318   Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
319   if(contents<1.){
320     //    Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
321     ResetHistograms();
322     fCurrentVertex = new AliESDVertex(startPos,startCov,99999.,-1);
323     points.Clear();
324     return;
325   }
326
327   TH1F *hc = fZCombc;
328
329   
330   if(hc->GetBinContent(hc->GetMaximumBin())<3)hc->Rebin(4);
331   Int_t binmin,binmax;
332   Int_t nPeaks=GetPeakRegion(hc,binmin,binmax);   
333   if(nPeaks==2)AliDebug(2,"2 peaks found");
334   Float_t zm =0.;
335   Float_t ezm =0.;
336   Float_t lim1 = hc->GetBinLowEdge(binmin);
337   Float_t lim2 = hc->GetBinLowEdge(binmax)+hc->GetBinWidth(binmax);
338   Float_t widthSR=lim2-lim1;
339
340   if(nPeaks ==1 && (lim2-lim1)<fWindowWidth){
341     Float_t c=(lim1+lim2)/2.;
342     lim1=c-fWindowWidth/2.;
343     lim2=c+fWindowWidth/2.;
344   }
345   Int_t niter = 0, ncontr=0;
346   do {
347     // symmetrization
348     if(zm  !=0.){
349       Float_t semilarg=TMath::Min((lim2-zm),(zm-lim1));
350       lim1=zm - semilarg;
351       lim2=zm + semilarg;
352     }
353
354     zm=0.;
355     ezm=0.;
356     ncontr=0;
357     for(Int_t i =0; i<points.GetEntriesFast(); i++){
358       AliITSZPoint* p=(AliITSZPoint*)points.UncheckedAt(i);
359       if(p->GetZ()>lim1 && p->GetZ()<lim2){
360         Float_t deno = p->GetErrZ();
361         zm+=p->GetZ()/deno;
362         ezm+=1./deno;
363         ncontr++;
364       }
365     }
366     if(ezm>0) {
367       zm/=ezm;
368       ezm=TMath::Sqrt(1./ezm);
369     }
370     niter++;
371   } while(niter<10 && TMath::Abs((zm-lim1)-(lim2-zm))>fTolerance);
372   if(nPeaks==2) ezm=widthSR;
373   Double_t position[3]={GetNominalPos()[0],GetNominalPos()[1],zm};
374   Double_t covmatrix[6]={GetNominalCov()[0],0.,GetNominalCov()[2],0.,0.,ezm};
375   fCurrentVertex = new AliESDVertex(position,covmatrix,99999.,ncontr);
376   fCurrentVertex->SetTitle("vertexer: Z");
377   fCurrentVertex->SetDispersion(fDiffPhiMax);
378   fNoVertices=1;
379   points.Clear();
380   if(fSearchForPileup && ncontr>fMinTrackletsForPilup){ 
381     Float_t secPeakPos;
382     Int_t ncontr2=FindSecondPeak(fZCombc,binmin,binmax,secPeakPos);
383     if(ncontr2>=fMinTrackletsForPilup){ 
384       fIsPileup=kTRUE;
385       fNoVertices=2;
386       fZpuv=secPeakPos;
387       fNTrpuv=ncontr2;
388       AliESDVertex secondVert(secPeakPos,0.1,ncontr2);
389       fVertArray = new AliESDVertex[2];
390       fVertArray[0]=(*fCurrentVertex);
391       fVertArray[1]=secondVert;
392     }
393   }
394   if(fNoVertices==1){
395     fVertArray = new AliESDVertex[1];
396     fVertArray[0]=(*fCurrentVertex);      
397   }
398   
399   ResetHistograms();
400   return;
401 }
402
403 //_____________________________________________________________________
404 Int_t AliITSVertexerZ::FindSecondPeak(TH1F* h, Int_t binmin,Int_t binmax, Float_t& secPeakPos){ 
405   // Resets bin contents between binmin and binmax and then search 
406   // for a second peak position 
407   for(Int_t i=binmin-1;i<=binmax+1;i++){
408     h->SetBinContent(i,0.);
409   }
410   Int_t secPeakBin=h->GetMaximumBin();
411   secPeakPos=h->GetBinCenter(secPeakBin);
412   Int_t secPeakCont=(Int_t)h->GetBinContent(secPeakBin);
413   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin-1);
414   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin+1);  
415   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin-2);
416   secPeakCont+=(Int_t)h->GetBinContent(secPeakBin+2);  
417   return secPeakCont;
418 }
419
420 //_____________________________________________________________________
421 void AliITSVertexerZ::ResetHistograms(){
422   // delete TH1 data members
423   if(fZCombc)delete fZCombc;
424   fZCombc = 0;
425 }
426
427 //________________________________________________________
428 void AliITSVertexerZ::PrintStatus() const {
429   // Print current status
430   cout <<"=======================================================\n";
431   cout <<" First layer first and last modules: "<<fFirstL1<<", ";
432   cout <<fLastL1<<endl;
433   cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
434   cout <<fLastL2<<endl;
435   cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
436   cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
437   cout <<"Bin sizes for coarse z histos "<<fStepCoarse<<endl;
438   cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
439   if(fZCombc){
440     cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
441   }
442   else{
443     cout<<"fZCombc does not exist\n";
444   }
445  
446   cout <<"=======================================================\n";
447 }
448