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