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