1 /**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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"
20 #include "AliITSLoader.h"
21 #include "AliITSgeom.h"
22 #include "AliITSDetTypeRec.h"
23 #include "AliITSRecPoint.h"
25 /////////////////////////////////////////////////////////////////
26 // this class implements a fast method to determine
27 // the Z coordinate of the primary vertex
28 // for p-p collisions it seems to give comparable or better results
29 // with respect to what obtained with AliITSVertexerPPZ
30 // It can be used successfully with Pb-Pb collisions
31 ////////////////////////////////////////////////////////////////
33 ClassImp(AliITSVertexerZ)
37 //______________________________________________________________________
38 AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer() {
39 // Default constructor
43 SetFirstLayerModules(0);
44 SetSecondLayerModules(0);
52 SetBinWidthCoarse(0.);
59 //______________________________________________________________________
60 AliITSVertexerZ::AliITSVertexerZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn) {
61 // Standard Constructor
65 SetFirstLayerModules();
66 SetSecondLayerModules();
82 //______________________________________________________________________
83 AliITSVertexerZ::AliITSVertexerZ(const AliITSVertexerZ &vtxr) : AliITSVertexer(vtxr) {
85 // Copies are not allowed. The method is protected to avoid misuse.
86 Error("AliITSVertexerZ","Copy constructor not allowed\n");
89 //______________________________________________________________________
90 AliITSVertexerZ& AliITSVertexerZ::operator=(const AliITSVertexerZ& /* vtxr */){
91 // Assignment operator
92 // Assignment is not allowed. The method is protected to avoid misuse.
93 Error("= operator","Assignment operator not allowed\n");
97 //______________________________________________________________________
98 AliITSVertexerZ::~AliITSVertexerZ() {
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.05,0.1,0.2,0.3,0.5};
112 Error("ConfigIterations","Maximum number of iterations is 5\n");
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];
119 //______________________________________________________________________
120 AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(Int_t evnumber){
121 // Defines the AliESDVertex for the current event
122 VertexZFinder(evnumber);
124 for(Int_t iteraz=0;iteraz<fMaxIter;iteraz++){
125 if(fCurrentVertex) ntrackl=fCurrentVertex->GetNContributors();
126 if(!fCurrentVertex || ntrackl==0 || ntrackl==-1){
127 Float_t diffPhiMaxOrig=fDiffPhiMax;
128 fDiffPhiMax=GetPhiMaxIter(iteraz);
129 VertexZFinder(evnumber);
130 fDiffPhiMax=diffPhiMaxOrig;
133 FindMultiplicity(evnumber);
134 return fCurrentVertex;
140 //______________________________________________________________________
141 void AliITSVertexerZ::VertexZFinder(Int_t evnumber){
142 // Defines the AliESDVertex for the current event
144 AliRunLoader *rl =AliRunLoader::GetRunLoader();
145 AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
146 AliITSgeom* geom = itsLoader->GetITSgeom();
147 itsLoader->LoadRecPoints();
148 rl->GetEvent(evnumber);
150 AliITSDetTypeRec detTypeRec;
152 TTree *tR = itsLoader->TreeR();
153 detTypeRec.SetTreeAddressR(tR);
154 TClonesArray *itsRec = 0;
155 // lc and gc are local and global coordinates for layer 1
156 Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
157 Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
158 // lc2 and gc2 are local and global coordinates for layer 2
159 Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
160 Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
162 itsRec = detTypeRec.RecPoints();
164 branch = tR->GetBranch("ITSRecPoints");
168 // By default fFirstL1=0 and fLastL1=79
169 // This loop counts the number of recpoints on layer1 (central modules)
170 for(Int_t module= fFirstL1; module<=fLastL1;module++){
171 // Keep only central modules
172 if(module%4==0 || module%4==3)continue;
173 // cout<<"Procesing module "<<module<<" ";
174 branch->GetEvent(module);
175 // cout<<"Number of clusters "<<clusters->GetEntries()<<endl;
176 nrpL1+= itsRec->GetEntries();
177 detTypeRec.ResetRecPoints();
179 //By default fFirstL2=80 and fLastL2=239
180 //This loop counts the number of RP on layer 2
181 for(Int_t module= fFirstL2; module<=fLastL2;module++){
182 branch->GetEvent(module);
183 nrpL2+= itsRec->GetEntries();
184 detTypeRec.ResetRecPoints();
186 if(nrpL1 == 0 || nrpL2 == 0){
188 itsLoader->UnloadRecPoints();
189 fCurrentVertex = new AliESDVertex(0.,5.3,-2);
192 // The vertex finding is attempted only if the number of RP is !=0 on
194 Float_t *xc1 = new Float_t [nrpL1]; // coordinates of the L1 Recpoints
195 Float_t *yc1 = new Float_t [nrpL1];
196 Float_t *zc1 = new Float_t [nrpL1];
197 Float_t *phi1 = new Float_t [nrpL1];
198 Float_t *xc2 = new Float_t [nrpL2]; // coordinates of the L1 Recpoints
199 Float_t *yc2 = new Float_t [nrpL2];
200 Float_t *zc2 = new Float_t [nrpL2];
201 Float_t *phi2 = new Float_t [nrpL2];
202 Int_t ind = 0;// running index for RP
203 // Force a coarse bin size of 200 microns if the number of clusters on layer 2
205 if(nrpL2<fPPsetting[0])SetBinWidthCoarse(fPPsetting[1]);
206 // By default nbinfine = (10+10)/0.0005=40000
207 Int_t nbinfine = static_cast<Int_t>((fHighLim-fLowLim)/fStepFine);
208 // By default nbincoarse=(10+10)/0.01=2000
209 Int_t nbincoarse = static_cast<Int_t>((fHighLim-fLowLim)/fStepCoarse);
210 // Set stepverycoarse = 3*fStepCoarse
211 Int_t nbinvcoarse = static_cast<Int_t>((fHighLim-fLowLim)/(3.*fStepCoarse));
212 if(fZCombc)delete fZCombc;
213 fZCombc = new TH1F("fZCombc","Z",nbincoarse,fLowLim,fLowLim+nbincoarse*fStepCoarse);
214 if(fZCombv)delete fZCombv;
215 fZCombv = new TH1F("fZCombv","Z",nbinvcoarse,fLowLim,fLowLim+nbinvcoarse*3.*fStepCoarse);
216 if(fZCombf)delete fZCombf;
217 fZCombf = new TH1F("fZCombf","Z",nbinfine,fLowLim,fLowLim+nbinfine*fStepFine);
218 // Loop on modules of layer 1 (restricted to central modules)
219 for(Int_t module= fFirstL1; module<=fLastL1;module++){
220 if(module%4==0 || module%4==3)continue;
221 branch->GetEvent(module);
222 Int_t nrecp1 = itsRec->GetEntries();
223 for(Int_t j=0;j<nrecp1;j++){
224 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
225 // Local coordinates of this recpoint
226 lc[0]=recp->GetDetLocalX();
227 lc[2]=recp->GetDetLocalZ();
228 geom->LtoG(module,lc,gc);
229 // Global coordinates of this recpoints
230 gc[0]-=fX0; // Possible beam offset in the bending plane
235 // azimuthal angle is computed in the interval 0 --> 2*pi
236 phi1[ind] = TMath::ATan2(gc[1],gc[0]);
237 if(phi1[ind]<0)phi1[ind]=2*TMath::Pi()+phi1[ind];
240 detTypeRec.ResetRecPoints();
242 ind = 0; // the running index is reset for Layer 2
243 for(Int_t module= fFirstL2; module<=fLastL2;module++){
244 branch->GetEvent(module);
245 Int_t nrecp2 = itsRec->GetEntries();
246 for(Int_t j=0;j<nrecp2;j++){
247 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
248 lc[0]=recp->GetDetLocalX();
249 lc[2]=recp->GetDetLocalZ();
250 geom->LtoG(module,lc,gc);
256 phi2[ind] = TMath::ATan2(gc[1],gc[0]);
257 if(phi2[ind]<0)phi2[ind]=2*TMath::Pi()+phi2[ind];
260 detTypeRec.ResetRecPoints();
264 for(Int_t i=0;i<nrpL1;i++){ // loop on L1 RP
265 Float_t r1=TMath::Sqrt(xc1[i]*xc1[i]+yc1[i]*yc1[i]); // radius L1 RP
266 for(Int_t j=0;j<nrpL2;j++){ // loop on L2 RP
267 Float_t diff = TMath::Abs(phi2[j]-phi1[i]); // diff in azimuth
268 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff; //diff<pi
269 if(diff<fDiffPhiMax){ // cut on 10 milliradians by def.
270 Float_t r2=TMath::Sqrt(xc2[j]*xc2[j]+yc2[j]*yc2[j]); // radius L2 RP
271 Float_t zr0=(r2*zc1[i]-r1*zc2[j])/(r2-r1); //Z @ null radius
283 MakeTracklet(pA,pB,nolines);
295 Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
297 // Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
299 itsLoader->UnloadRecPoints();
300 fCurrentVertex = new AliESDVertex(0.,5.3,-1);
305 if(fDebug>0)cout<<"Number of entries in hist. "<<fZCombc->GetEntries()<<endl;
310 Bool_t goon = kFALSE;
317 bi = hc->GetMaximumBin(); // bin with maximal content on coarse histogram
318 if(hc->GetBinContent(bi)<3){
319 if(cnt==1)goon = kTRUE;
325 Float_t centre = hc->GetBinCenter(bi); // z value of the bin with maximal content
327 Int_t bifine=static_cast<Int_t>((hc->GetBinCenter(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
328 Int_t nbinsfine=static_cast<Int_t>(3*hc->GetBinWidth(bi)/fStepFine);
329 // evaluation of the centroid
330 Int_t ii1=TMath::Max(bifine-nbinsfine,1);
331 Int_t ii2=TMath::Min(bifine+nbinsfine,fZCombf->GetNbinsX());
334 for(Int_t ii=ii1;ii<ii2;ii++){
335 centre+= fZCombf->GetBinCenter(ii)*fZCombf->GetBinContent(ii);
336 nn+=static_cast<Int_t>(fZCombf->GetBinContent(ii));
341 cout<<"Value of center "<<centre<<endl;
342 cout<<"Population of 3 central bins: "<<hc->GetBinContent(bi-1)<<", ";
343 cout<<hc->GetBinContent(bi)<<", ";
344 cout<<hc->GetBinContent(bi+1)<<endl;
347 // n1 is the bin number of fine histogram containing the point located 1 coarse bin less than "centre"
348 Int_t n1 = static_cast<Int_t>((centre-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
349 // n2 is the bin number of fine histogram containing the point located 1 coarse bin more than "centre"
350 Int_t n2 = static_cast<Int_t>((centre+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
352 if(n2>nbinfine)n2=nbinfine;
356 Bool_t last = kFALSE;
362 // at the end of the loop:
363 // fZFound = N*(Average Z) - where N is the number of tracklets
365 // fZsig = N*Q - where Q is the average of Z*Z
366 for(Int_t n=n1;n<=n2;n++){
367 fZFound+=fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
368 num+=static_cast<Int_t>(fZCombf->GetBinContent(n));
369 fZsig+=fZCombf->GetBinCenter(n)*fZCombf->GetBinCenter(n)*fZCombf->GetBinContent(n);
372 fZsig = 5.3; // Default error from the beam sigmoid
375 Float_t radi = fZsig/(num-1)-fZFound*fZFound/num/(num-1);
376 // radi = square root of sample variance of Z
377 if(radi>0.)fZsig=TMath::Sqrt(radi);
378 else fZsig=5.3; // Default error from the beam sigmoid
379 // fZfound - Average Z
383 // goon is true if the distance between the found Z and the lower bin differs from the distance between the found Z and
384 // the upper bin by more than tolerance (0.002)
385 goon = TMath::Abs(TMath::Abs(fZFound-fZCombf->GetBinCenter(n1))-TMath::Abs(fZFound-fZCombf->GetBinCenter(n2)))>fTolerance;
386 // a window in the fine grained histogram is centered aroung the found Z. The width is 2 bins of
387 // the coarse grained histogram
389 n1 = static_cast<Int_t>((fZFound-hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
391 n2 = static_cast<Int_t>((fZFound+hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
392 if(n2>nbinfine)n2=nbinfine;
395 cout<<"Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
396 cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
401 n1 = static_cast<Int_t>((centre-(niter+2)*hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
402 n2 = static_cast<Int_t>((centre+(niter+2)*hc->GetBinWidth(bi)-fZCombf->GetBinLowEdge(0))/fStepFine);
404 if(n2>nbinfine)n2=nbinfine;
407 // no more than 10 adjusting iterations
410 // Warning("FindVertexForCurrentEvent","The procedure does not converge\n");
413 if((fZsig> 0.0001) && !goon && num>5){
415 n1 = static_cast<Int_t>((fZFound-fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
417 n2 = static_cast<Int_t>((fZFound+fZsig-fZCombf->GetBinLowEdge(0))/fZCombf->GetBinWidth(0));
418 if(n2>nbinfine)n2=nbinfine;
421 cout<<"FINAL Search restricted to n1 = "<<n1<<", n2= "<<n2<<endl;
422 cout<<"z1= "<<fZCombf->GetBinCenter(n1)<<", n2 = "<<fZCombf->GetBinCenter(n2)<<endl;
432 // if(fDebug>0)cout<<"Numer of Iterations "<<niter<<endl<<endl;
433 // if(num!=0)fZsig/=TMath::Sqrt(num);
434 if (fZsig<=0) fZsig=5.3; // Default error from the beam sigmoid
435 fCurrentVertex = new AliESDVertex(fZFound,fZsig,num);
436 fCurrentVertex->SetTitle("vertexer: B");
438 itsLoader->UnloadRecPoints();
442 //_____________________________________________________________________
443 void AliITSVertexerZ::ResetHistograms(){
444 // delete TH1 data members
445 if(fZCombc)delete fZCombc;
446 if(fZCombf)delete fZCombf;
447 if(fZCombv)delete fZCombv;
453 //______________________________________________________________________
454 void AliITSVertexerZ::FindVertices(){
455 // computes the vertices of the events in the range FirstEvent - LastEvent
456 AliRunLoader *rl = AliRunLoader::GetRunLoader();
457 AliITSLoader* itsLoader = (AliITSLoader*) rl->GetLoader("ITSLoader");
458 itsLoader->ReloadRecPoints();
459 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
460 // cout<<"Processing event "<<i<<endl;
462 FindVertexForCurrentEvent(i);
464 WriteCurrentVertex();
469 cout<<"Vertex not found for event "<<i<<endl;
470 cout<<"fZFound = "<<fZFound<<", fZsig= "<<fZsig<<endl;
477 //________________________________________________________
478 void AliITSVertexerZ::PrintStatus() const {
479 // Print current status
480 cout <<"=======================================================\n";
481 cout <<" First layer first and last modules: "<<fFirstL1<<", ";
482 cout <<fLastL1<<endl;
483 cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
484 cout <<fLastL2<<endl;
485 cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
486 cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
487 cout <<"Bin sizes for coarse and fine z histos "<<fStepCoarse<<"; "<<fStepFine<<endl;
488 cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
489 cout <<" Debug flag: "<<fDebug<<endl;
490 cout <<"First event to be processed "<<fFirstEvent;
491 cout <<"\n Last event to be processed "<<fLastEvent<<endl;
493 cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
496 cout<<"fZCombc does not exist\n";
499 cout<<"fZCombv exists - entries="<<fZCombv->GetEntries()<<endl;
502 cout<<"fZCombv does not exist\n";
505 cout<<"fZCombf exists - entries="<<fZCombv->GetEntries()<<endl;
508 cout<<"fZCombf does not exist\n";
511 cout <<"=======================================================\n";