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