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