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"
17 #include<TClonesArray.h>
21 #include "AliITSLoader.h"
22 #include "AliITSgeomTGeo.h"
23 #include "AliITSDetTypeRec.h"
24 #include "AliITSRecPoint.h"
25 #include "AliITSZPoint.h"
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 ////////////////////////////////////////////////////////////////
35 ClassImp(AliITSVertexerZ)
39 //______________________________________________________________________
40 AliITSVertexerZ::AliITSVertexerZ():AliITSVertexer(),
55 // Default constructor
57 SetFirstLayerModules();
58 SetSecondLayerModules();
68 //______________________________________________________________________
69 AliITSVertexerZ::AliITSVertexerZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn),
84 // Standard Constructor
86 SetFirstLayerModules();
87 SetSecondLayerModules();
95 SetVtxStart((Double_t)x0,(Double_t)y0,0.);
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),
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){
119 //______________________________________________________________________
120 AliITSVertexerZ& AliITSVertexerZ::operator=(const AliITSVertexerZ& vtxr ){
121 // Assignment operator
123 this->~AliITSVertexerZ();
124 new(this) AliITSVertexerZ(vtxr);
128 //______________________________________________________________________
129 AliITSVertexerZ::~AliITSVertexerZ() {
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};
141 Error("ConfigIterations","Maximum number of iterations is 5\n");
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];
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();
154 Int_t bi1=h->GetMaximumBin();
156 for(Int_t i=imax;i>=1;i--){
157 if(h->GetBinContent(i)>maxval){
158 maxval=h->GetBinContent(i);
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);
191 if(bi2-bi1==1) npeaks=1;
198 //______________________________________________________________________
199 AliESDVertex* AliITSVertexerZ::FindVertexForCurrentEvent(Int_t evnumber){
200 // Defines the AliESDVertex for the current event
201 VertexZFinder(evnumber);
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;
212 FindMultiplicity(evnumber);
213 return fCurrentVertex;
216 //______________________________________________________________________
217 void AliITSVertexerZ::VertexZFinder(Int_t evnumber){
218 // Defines the AliESDVertex for the current event
220 AliRunLoader *rl =AliRunLoader::GetRunLoader();
221 AliITSLoader* itsLoader = (AliITSLoader*)rl->GetLoader("ITSLoader");
222 // AliITSgeom* geom = itsLoader->GetITSgeom();
223 itsLoader->LoadRecPoints();
224 rl->GetEvent(evnumber);
226 AliITSDetTypeRec detTypeRec;
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.;
238 itsRec = detTypeRec.RecPoints();
240 branch = tR->GetBranch("ITSRecPoints");
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();
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();
257 if(nrpL1 == 0 || nrpL2 == 0){
259 itsLoader->UnloadRecPoints();
260 fCurrentVertex = new AliESDVertex(0.,5.3,-2);
263 // Force a coarse bin size of 200 microns if the number of clusters on layer 2
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);
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);
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;
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);
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);
297 detTypeRec.ResetRecPoints();
298 for(Int_t j1=0;j1<nrecp1;j1++){
299 AliITSRecPoint *recp = (AliITSRecPoint*)prpl1->At(j1);
301 lc1[0]=recp->GetDetLocalX();
302 lc1[2]=recp->GetDetLocalZ();
303 geom->LtoG(modul1,lc1,gc1);
304 // Global coordinates of this recpoints
306 recp->GetGlobalXYZ(gc1);
307 gc1[0]-=fNominalPos[0]; // Possible beam offset in the bending plane
308 gc1[1]-=fNominalPos[1]; // " "
309 Float_t r1=TMath::Sqrt(gc1[0]*gc1[0]+gc1[1]*gc1[1]);
310 Float_t phi1 = TMath::ATan2(gc1[1],gc1[0]);
311 if(phi1<0)phi1+=2*TMath::Pi();
313 Float_t erz1=recp->GetSigmaZ2();
314 for(Int_t ladl2=0 ; ladl2<fLadOnLay2*2+1;ladl2++){
315 for(Int_t k=0;k<4;k++){
316 Int_t ladmod=fLadders[ladder-1]+ladl2;
317 if(ladmod>AliITSgeomTGeo::GetNLadders(2)) ladmod=ladmod-AliITSgeomTGeo::GetNLadders(2);
318 Int_t modul2=AliITSgeomTGeo::GetModuleIndex(2,ladmod,k+1);
319 branch->GetEvent(modul2);
320 Int_t nrecp2 = itsRec->GetEntries();
321 for(Int_t j2=0;j2<nrecp2;j2++){
322 recp = (AliITSRecPoint*)itsRec->At(j2);
324 lc2[0]=recp->GetDetLocalX();
325 lc2[2]=recp->GetDetLocalZ();
326 geom->LtoG(modul2,lc2,gc2);
328 recp->GetGlobalXYZ(gc2);
329 gc2[0]-=fNominalPos[0];
330 gc2[1]-=fNominalPos[1];
331 Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
332 Float_t phi2 = TMath::ATan2(gc2[1],gc2[0]);
333 if(phi2<0)phi2+=2*TMath::Pi();
335 Float_t erz2=recp->GetSigmaZ2();
337 Float_t diff = TMath::Abs(phi2-phi1);
338 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
339 if(diff<fDiffPhiMax){
340 // Float_t tgth=(zc2[j]-zc1[i])/(r2-r1); // slope (used for multiple scattering)
341 Float_t zr0=(r2*zc1-r1*zc2)/(r2-r1); //Z @ null radius
342 Float_t ezr0q=(r2*r2*erz1+r1*r1*erz2)/(r2-r1)/(r2-r1); //error on Z @ null radius
344 // Multiple scattering
345 ezr0q+=r1*r1*(1+tgth*tgth)*theta2L1/2; // multiple scattering in layer 1
346 ezr0q+=rBP*rBP*(1+tgth*tgth)*theta2BP/2; // multiple scattering in beam pipe
348 if(nopoints<maxdim) new(pts[nopoints++])AliITSZPoint(zr0,ezr0q);
352 detTypeRec.ResetRecPoints();
361 Double_t contents = fZCombc->GetEntries()- fZCombc->GetBinContent(0)-fZCombc->GetBinContent(nbincoarse+1);
363 // Warning("FindVertexForCurrentEvent","Insufficient number of rec. points\n");
365 itsLoader->UnloadRecPoints();
366 fCurrentVertex = new AliESDVertex(0.,5.3,-1);
373 if(hc->GetBinContent(hc->GetMaximumBin())<3)hc->Rebin(3);
375 Int_t nPeaks=GetPeakRegion(hc,binmin,binmax);
376 if(nPeaks==2)AliWarning("2 peaks found");
379 Float_t lim1 = hc->GetBinLowEdge(binmin);
380 Float_t lim2 = hc->GetBinLowEdge(binmax)+hc->GetBinWidth(binmax);
382 if(nPeaks ==1 && (lim2-lim1)<fWindowWidth){
383 Float_t c=(lim1+lim2)/2.;
384 lim1=c-fWindowWidth/2.;
385 lim2=c+fWindowWidth/2.;
387 Int_t niter = 0, ncontr=0;
391 Float_t semilarg=TMath::Min((lim2-zm),(zm-lim1));
399 for(Int_t i =0; i<points->GetEntriesFast(); i++){
400 AliITSZPoint* p=(AliITSZPoint*)points->UncheckedAt(i);
401 if(p->GetZ()>lim1 && p->GetZ()<lim2){
402 Float_t deno = p->GetErrZ();
409 ezm=TMath::Sqrt(1./ezm);
411 } while(niter<10 && TMath::Abs((zm-lim1)-(lim2-zm))>fTolerance);
412 fCurrentVertex = new AliESDVertex(zm,ezm,ncontr);
413 fCurrentVertex->SetTitle("vertexer: B");
416 itsLoader->UnloadRecPoints();
422 //_____________________________________________________________________
423 void AliITSVertexerZ::ResetHistograms(){
424 // delete TH1 data members
425 if(fZCombc)delete fZCombc;
429 //______________________________________________________________________
430 void AliITSVertexerZ::FindVertices(){
431 // computes the vertices of the events in the range FirstEvent - LastEvent
432 AliRunLoader *rl = AliRunLoader::GetRunLoader();
433 AliITSLoader* itsLoader = (AliITSLoader*) rl->GetLoader("ITSLoader");
434 itsLoader->ReloadRecPoints();
435 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
436 // cout<<"Processing event "<<i<<endl;
438 FindVertexForCurrentEvent(i);
440 WriteCurrentVertex();
445 //________________________________________________________
446 void AliITSVertexerZ::PrintStatus() const {
447 // Print current status
448 cout <<"=======================================================\n";
449 cout <<" First layer first and last modules: "<<fFirstL1<<", ";
450 cout <<fLastL1<<endl;
451 cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
452 cout <<fLastL2<<endl;
453 cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
454 cout <<"Limits for Z histograms: "<<fLowLim<<"; "<<fHighLim<<endl;
455 cout <<"Bin sizes for coarse z histos "<<fStepCoarse<<endl;
456 cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
457 cout <<"First event to be processed "<<fFirstEvent;
458 cout <<"\n Last event to be processed "<<fLastEvent<<endl;
460 cout<<"fZCombc exists - entries="<<fZCombc->GetEntries()<<endl;
463 cout<<"fZCombc does not exist\n";
466 cout <<"=======================================================\n";