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 <AliITSVertexerPPZ.h>
18 #include <TObjArray.h>
20 #include <TClonesArray.h>
21 #include "AliITSDetTypeRec.h"
22 #include "AliITSgeom.h"
23 #include "AliITSLoader.h"
24 #include "AliITSRecPoint.h"
25 /////////////////////////////////////////////////////////////////////////
27 // This class is intended to compute the Z coordinate of the //
28 // primary vertex for p-p interactions, or in general when the //
29 // number of secondaries is too slow to use the class //
30 // AliITSVertexerIons, which in turn is optimized for A-A collsions //
31 // Origin: masera@to.infn.it 9/12/2002 //
33 /////////////////////////////////////////////////////////////////////////
34 ClassImp(AliITSVertexerPPZ)
38 //______________________________________________________________________
39 AliITSVertexerPPZ::AliITSVertexerPPZ():AliITSVertexer(),
50 // Default Constructor
53 SetFirstLayerModules(0);
54 SetSecondLayerModules(0);
59 AliITSVertexerPPZ::AliITSVertexerPPZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn),
70 // Standard constructor
74 SetFirstLayerModules();
75 SetSecondLayerModules();
80 //______________________________________________________________________
81 AliITSVertexerPPZ::~AliITSVertexerPPZ() {
86 //________________________________________________________
87 void AliITSVertexerPPZ::EvalZ(TH1F *hist,Int_t sepa, Int_t ncoinc, TArrayF *zval) {
90 Float_t deltaVal = hist->GetBinWidth(1)*fWindow; // max window in Z for searching
98 for(Int_t i=1;i<=sepa;i++){
99 Float_t cont=hist->GetBinContent(i);
100 if(cont!=0)curz = hist->GetBinLowEdge(i)+0.5*hist->GetBinWidth(i);
104 if(cont!=0)nbinNotZero++;
106 if(nbinNotZero==0){fZFound=-100; fZsig=-100; return;}
110 fCurrentVertex = new AliESDVertex(fZFound,fZsig,nbinNotZero);
113 Float_t errsq = totst2/(nN-1)-totst*totst/nN/(nN-1);
115 totst2=TMath::Sqrt(errsq);
118 Error("EvalZ","Negative variance: %d - %12.7f - %12.7f",nN,totst2,totst);
124 Float_t cut = totst+totst2*2.;
125 if(fDebug>1)cout<<"totst, totst2, cut: "<<totst<<", "<<totst2<<", "<<cut<<endl;
126 Float_t val1=hist->GetBinLowEdge(sepa);
127 Float_t val2=hist->GetBinLowEdge(1);
130 for(Int_t i=1;i<=sepa;i++){
131 Float_t cont=hist->GetBinContent(i);
134 zmax = hist->GetBinLowEdge(i)+0.5*hist->GetBinWidth(1);
137 curz=hist->GetBinLowEdge(i);
138 if(curz<val1)val1=curz;
139 if(curz>val2)val2=curz;
142 val2+=hist->GetBinWidth(1);
143 if((val2-val1)>deltaVal){
144 val1 = zmax-deltaVal/2.;
145 val2 = zmax+deltaVal/2.;
146 if(fDebug>0)cout<<"val1 and val2 recomputed\n";
148 if(fDebug>0)cout<<"Values for Z finding: "<<val1<<" "<<val2<<" "<<val2-val1<<endl;
152 for(Int_t i=0; i<ncoinc; i++){
153 Float_t z=zval->At(i);
159 /* weights defined by the curvature
160 Float_t wei = 1./curv->At(i);
166 if(nN<1){fZFound=-110; fZsig=-110; return;}
169 fCurrentVertex = new AliESDVertex(fZFound,fZsig,nN);
172 errsq = (fZsig/(nN-1)-fZFound*fZFound/nN/(nN-1))/nN;
174 fZsig=TMath::Sqrt(errsq);
177 Error("evalZ","Negative variance: %d - %12.7f %12.7f",nN,fZsig,fZFound);
181 /* weights defined by the curvature
184 fZsig = TMath::Sqrt(fZsig);
186 fCurrentVertex = new AliESDVertex(fZFound,fZsig,nN);
189 //______________________________________________________________________
190 AliESDVertex* AliITSVertexerPPZ::FindVertexForCurrentEvent(Int_t evnumber){
191 // Defines the AliESDVertex for the current event
195 AliRunLoader *rl =AliRunLoader::GetRunLoader();
196 AliITSLoader* iTSloader = (AliITSLoader*)rl->GetLoader("ITSLoader");
197 TDirectory * olddir = gDirectory;
199 AliITSgeom* geom = (AliITSgeom*)gDirectory->Get("AliITSgeom");
202 AliITSDetTypeRec detTypeRec;
205 Error("FindVertexForCurrentEvent","ITS geometry is not defined");
206 return fCurrentVertex;
209 TClonesArray *itsRec = 0;
210 Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
211 Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
212 Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
213 Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
215 tR = iTSloader->TreeR();
217 Error("FindVertexForCurrentEvent","TreeR not found");
218 return fCurrentVertex;
220 detTypeRec.SetTreeAddressR(tR);
221 itsRec = detTypeRec.RecPoints();
223 // TClonesArray dummy("AliITSRecPoint",10000), *clusters=&dummy;
225 branch = tR->GetBranch("ITSRecPoints");
227 branch = tR->GetBranch("ITSRecPointsF");
231 Error("FindVertexForCurrentEvent","branch for ITS rec points not found");
232 return fCurrentVertex;
238 for(Int_t module= fFirstL1; module<=fLastL1;module++){
239 branch->GetEvent(module);
240 Int_t nrecp1 = itsRec->GetEntries();
241 for(Int_t i=0; i<nrecp1;i++){
242 AliITSRecPoint *current = (AliITSRecPoint*)itsRec->At(i);
243 lc[0]=current->GetDetLocalX();
244 lc[2]=current->GetDetLocalZ();
245 geom->LtoG(module,lc,gc);
250 detTypeRec.ResetRecPoints();
253 rmszav=TMath::Sqrt(zave2/(firipixe-1)-zave*zave/firipixe/(firipixe-1));
255 if(fDebug>1)cout<<"++++++++++++++++++++++++++++++++++++++++++++++++++++++\n Number of firing pixels: "<<firipixe<<endl;
260 Warning("FindVertexForCurrentEvent","No rec points on first layer for this event");
261 return fCurrentVertex;
263 Float_t zlim1=zave-rmszav;
264 Float_t zlim2=zave+rmszav;
265 Int_t sepa=(Int_t)((zlim2-zlim1)*10.+1.);
266 zlim2=zlim1 + sepa/10.;
267 TH1F *zvdis = new TH1F("z_ev","zv distr",sepa,zlim1,zlim2);
269 cout<<"Z limits: "<<zlim1<<" "<<zlim2<<"; Bins= "<<sepa<<endl;
270 cout<<"Bin width: "<<zvdis->GetBinWidth(1)*10000.<<" micron\n";
273 TArrayF *zval = new TArrayF(sizarr);
274 // TArrayF *curv = new TArrayF(sizarr);
276 for(Int_t module= fFirstL1; module<=fLastL1;module++){
277 if(fDebug>0)cout<<"processing module "<<module<<" \r";
278 branch->GetEvent(module);
279 Int_t nrecp1 = itsRec->GetEntries();
280 TObjArray *poiL1 = new TObjArray(nrecp1);
281 for(Int_t i=0; i<nrecp1;i++)poiL1->AddAt(itsRec->At(i),i);
282 detTypeRec.ResetRecPoints();
283 for(Int_t i=0; i<nrecp1;i++){
284 AliITSRecPoint *current = (AliITSRecPoint*)poiL1->At(i);
285 lc[0]=current->GetDetLocalX();
286 lc[2]=current->GetDetLocalZ();
287 geom->LtoG(module,lc,gc);
290 Float_t r1=TMath::Sqrt(gc[0]*gc[0]+gc[1]*gc[1]);
291 Float_t phi1 = TMath::ATan2(gc[1],gc[0]);
292 if(phi1<0)phi1=2*TMath::Pi()+phi1;
293 if(fDebug>1)cout<<"module "<<module<<" "<<gc[0]<<" "<<gc[1]<<" "<<gc[2]<<" "<<phi1<<" \n";
294 for(Int_t modul2=fFirstL2; modul2<=fLastL2; modul2++){
295 branch->GetEvent(modul2);
296 Int_t nrecp2 = itsRec->GetEntries();
297 for(Int_t j=0; j<nrecp2;j++){
298 AliITSRecPoint *recp = (AliITSRecPoint*)itsRec->At(j);
299 lc2[0]=recp->GetDetLocalX();
300 lc2[2]=recp->GetDetLocalZ();
301 geom->LtoG(modul2,lc2,gc2);
304 Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
305 Float_t zr0=(r2*gc[2]-r1*gc2[2])/(r2-r1);
306 Float_t phi2 = TMath::ATan2(gc2[1],gc2[0]);
307 if(phi2<0)phi2=2.*TMath::Pi()+phi2;
308 Float_t diff = TMath::Abs(phi2-phi1);
309 if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
310 if(zr0>zlim1 && zr0<zlim2){
311 if(diff<fDiffPhiMax ){
313 zval->AddAt(zr0,ncoinc);
314 /* uncomment these lines to use curvature as a weight
315 Float_t cu = Curv(0.,0.,gc[0],gc[1],gc2[0],gc[1]);
316 curv->AddAt(cu,ncoinc);
319 if(ncoinc==(sizarr-1)){
322 //uncomment next line to use curvature as weight
323 // curv->Set(sizarr);
328 detTypeRec.ResetRecPoints();
334 cout<<endl<<"Number of coincidences = "<<ncoinc<<endl;
336 // EvalZ(zvdis,sepa,ncoinc,zval,curv);
337 EvalZ(zvdis,sepa,ncoinc,zval);
343 sprintf(name,"Vertex_%d",evnumber);
344 // fCurrentVertex->SetName(name);
345 fCurrentVertex->SetTitle("vertexer: PPZ");
348 return fCurrentVertex;
351 //______________________________________________________________________
352 void AliITSVertexerPPZ::FindVertices(){
353 // computes the vertices of the events in the range FirstEvent - LastEvent
354 AliRunLoader *rl = AliRunLoader::GetRunLoader();
355 AliITSLoader* iTSloader = (AliITSLoader*) rl->GetLoader("ITSLoader");
356 iTSloader->ReloadRecPoints();
357 for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
358 cout<<"Processing event "<<i<<endl;
360 FindVertexForCurrentEvent(i);
362 WriteCurrentVertex();
366 cout<<"Vertex not found for event "<<i<<endl;
367 cout<<"fZFound = "<<fZFound<<", fZsig= "<<fZsig<<endl;
373 //________________________________________________________
374 void AliITSVertexerPPZ::PrintStatus() const {
375 // Print current status
376 cout <<"=======================================================\n";
377 cout <<" First layer first and last modules: "<<fFirstL1<<", ";
378 cout <<fLastL1<<endl;
379 cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
380 cout <<fLastL2<<endl;
381 cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
382 cout <<" Window for Z search: "<<fWindow<<endl;
383 cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
384 cout <<" Debug flag: "<<fDebug<<endl;
385 cout<<"First event to be processed "<<fFirstEvent;
386 cout<<"\n Last event to be processed "<<fLastEvent<<endl;
389 //________________________________________________________
390 Float_t AliITSVertexerPPZ::Curv(Double_t x1,Double_t y1,
391 Double_t x2,Double_t y2,
392 Double_t x3,Double_t y3)
394 //-----------------------------------------------------------------
395 // Initial approximation of the track curvature (Y. Belikov) squared
396 //-----------------------------------------------------------------
397 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
398 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
399 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
400 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
401 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
403 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
405 Float_t val = (Float_t)(( -xr*yr/TMath::Sqrt(xr*xr+yr*yr))
406 *( -xr*yr/TMath::Sqrt(xr*xr+yr*yr)));