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