]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexerPPZ.cxx
Using TMath::Abs instead of fabs
[u/mrichter/AliRoot.git] / ITS / AliITSVertexerPPZ.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 <AliITSVertexerPPZ.h>
16 #include <Riostream.h>
17 #include <TArrayF.h>
18 #include <TDirectory.h>
19 #include <TH1F.h>
20 #include <TMath.h>
21 #include <TObjArray.h>
22 #include <TTree.h>
23 #include <TClonesArray.h>
24 #include "AliRun.h"
25 #include "AliITS.h"
26 #include "AliITSgeom.h"
27 #include "AliITSLoader.h"
28 #include "AliITSRecPoint.h"
29 /////////////////////////////////////////////////////////////////////////
30 //                                                                     //
31 // This class is intended to compute the Z coordinate of the           //
32 // primary vertex for p-p interactions, or in general when the         //
33 // number of secondaries is too slow to use the class                  //
34 // AliITSVertexerIons, which in turn is optimized for A-A collsions    //
35 // Origin: masera@to.infn.it      9/12/2002                            //
36 //                                                                     //
37 /////////////////////////////////////////////////////////////////////////
38 ClassImp(AliITSVertexerPPZ)
39
40
41
42 //______________________________________________________________________
43 AliITSVertexerPPZ::AliITSVertexerPPZ():AliITSVertexer() {
44   // Default Constructor
45
46   SetDiffPhiMax(0);
47   fX0 = 0.;
48   fY0 = 0.;
49   SetFirstLayerModules(0);
50   SetSecondLayerModules(0);
51   fZFound = 0;
52   fZsig = 0.;
53   fITS = 0;
54   SetWindow(0);
55 }
56
57 AliITSVertexerPPZ::AliITSVertexerPPZ(TString fn, Float_t x0, Float_t y0):AliITSVertexer(fn) {
58   // Standard constructor
59   SetDiffPhiMax();
60   fX0 = x0;
61   fY0 = y0;
62   SetFirstLayerModules();
63   SetSecondLayerModules();
64   fZFound = 0;
65   fZsig = 0.;
66   fITS = 0;
67   SetWindow();
68 }
69
70 //______________________________________________________________________
71 AliITSVertexerPPZ::~AliITSVertexerPPZ() {
72   // Default Destructor
73   fITS = 0;
74 }
75
76 //________________________________________________________
77 void  AliITSVertexerPPZ::EvalZ(TH1F *hist,Int_t sepa, Int_t ncoinc, TArrayF *zval) {
78   Float_t DeltaVal = hist->GetBinWidth(1)*fWindow; // max window in Z for searching
79   fZFound=0;
80   fZsig=0;
81   Int_t N=0;
82   Int_t NbinNotZero=0;
83   Float_t totst = 0.;
84   Float_t totst2 = 0.;
85   Float_t curz = 0.;
86   for(Int_t i=1;i<=sepa;i++){
87     Float_t cont=hist->GetBinContent(i);
88     if(cont!=0)curz = hist->GetBinLowEdge(i)+0.5*hist->GetBinWidth(i);
89     totst+=cont;
90     totst2+=cont*cont;
91     N++;
92     if(cont!=0)NbinNotZero++;
93   }
94   if(NbinNotZero==0){fZFound=-100; fZsig=-100; return;}
95   if(NbinNotZero==1){
96     fZFound = curz;
97     fZsig=0;
98     fCurrentVertex = new AliITSVertex(fZFound,fZsig,NbinNotZero);
99     return;
100   }
101   Float_t errsq = totst2/(N-1)-totst*totst/N/(N-1);
102   if(errsq>=0){
103   totst2=TMath::Sqrt(errsq);
104   }
105   else {
106     Error("EvalZ","Negative variance: %d - %12.7f - %12.7f",N,totst2,totst);
107     fZFound=-100; 
108     fZsig=-100;
109     return;
110   }
111   totst /= N;
112   Float_t cut = totst+totst2*2.;
113   if(fDebug>1)cout<<"totst, totst2, cut: "<<totst<<", "<<totst2<<", "<<cut<<endl;
114   Float_t val1=hist->GetBinLowEdge(sepa); 
115   Float_t val2=hist->GetBinLowEdge(1);
116   Float_t valm = 0.;
117   Float_t zmax = 0.;
118   for(Int_t i=1;i<=sepa;i++){
119     Float_t cont=hist->GetBinContent(i);
120     if(cont>valm){
121       valm = cont;
122       zmax = hist->GetBinLowEdge(i)+0.5*hist->GetBinWidth(1);
123     }
124     if(cont>cut){
125       curz=hist->GetBinLowEdge(i);
126       if(curz<val1)val1=curz;
127       if(curz>val2)val2=curz;
128     }
129   }
130   val2+=hist->GetBinWidth(1);
131   if((val2-val1)>DeltaVal){
132     val1 = zmax-DeltaVal/2.;
133     val2 = zmax+DeltaVal/2.;
134     if(fDebug>0)cout<<"val1 and val2 recomputed\n";
135   }
136   if(fDebug>0)cout<<"Values for Z finding: "<<val1<<" "<<val2<<" "<<val2-val1<<endl;
137   fZFound=0;
138   fZsig=0;
139   N=0;
140   for(Int_t i=0; i<ncoinc; i++){
141     Float_t z=zval->At(i);
142     if(z<val1)continue;
143     if(z>val2)continue;
144    
145     fZFound+=z;
146     fZsig+=z*z;
147     /*   weights defined by the curvature
148     Float_t wei = 1./curv->At(i);
149     fZFound+=z*wei;
150     fZsig+=wei;
151     */
152     N++;
153   }
154   if(N<1){fZFound=-110; fZsig=-110; return;}
155   if(N==1){
156     fZsig = 0;
157     fCurrentVertex = new AliITSVertex(fZFound,fZsig,N);
158     return;
159   }
160   errsq = (fZsig/(N-1)-fZFound*fZFound/N/(N-1))/N;
161   if(errsq>=0.){
162     fZsig=TMath::Sqrt(errsq);
163   }
164   else {
165     Error("evalZ","Negative variance: %d - %12.7f %12.7f",N,fZsig,fZFound);
166     fZsig=0;
167   }
168   fZFound=fZFound/N;
169   /* weights defined by the curvature
170   fZsig=1./fZsig;
171   fZFound*=fZsig;
172   fZsig = TMath::Sqrt(fZsig);
173   */
174   fCurrentVertex = new AliITSVertex(fZFound,fZsig,N);
175 }
176
177 //______________________________________________________________________
178 AliITSVertex* AliITSVertexerPPZ::FindVertexForCurrentEvent(Int_t evnumber){
179   // Defines the AliITSVertex for the current event
180   fCurrentVertex = 0;
181   fZFound = -999;
182   fZsig = -999;
183   AliRunLoader *rl =AliRunLoader::GetRunLoader();
184   AliITSLoader* ITSloader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
185   if(!fITS)  {
186     fITS = (AliITS*)gAlice->GetModule("ITS");
187     if(!fITS) {
188       Error("FindVertexForCurrentEvent","AliITS object was not found");
189       return fCurrentVertex;
190     }
191   }
192   fITS->SetTreeAddress();
193   AliITSgeom *geom = fITS->GetITSgeom();
194   if(!geom) {
195     Error("FindVertexForCurrentEvent","ITS geometry is not defined");
196     return fCurrentVertex;
197   }
198   TTree *TR=0;
199   TClonesArray *ITSrec  = 0;
200   Float_t lc[3]; for(Int_t ii=0; ii<3; ii++) lc[ii]=0.;
201   Float_t gc[3]; for(Int_t ii=0; ii<3; ii++) gc[ii]=0.;
202   Float_t lc2[3]; for(Int_t ii=0; ii<3; ii++) lc2[ii]=0.;
203   Float_t gc2[3]; for(Int_t ii=0; ii<3; ii++) gc2[ii]=0.;
204
205   TR = ITSloader->TreeR();
206   if(!TR){
207     Error("FindVertexForCurrentEvent","TreeR not found");
208     return fCurrentVertex;
209   }
210   ITSrec  = fITS->RecPoints();  // uses slow points or fast points if slow are
211   // missing
212   TBranch *branch = TR->GetBranch("ITSRecPoints");
213   if(!branch){ 
214     branch = TR->GetBranch("ITSRecPointsF");
215     //   Warning("FindVertexForCurrentEvent","Using Fast points");
216   }
217   if(!branch){
218    Error("FindVertexForCurrentEvent","branch for ITS rec points not found");
219    return fCurrentVertex;
220   }
221   Float_t zave=0;
222   Float_t rmszav=0;
223   Float_t zave2=0;
224   Int_t firipixe=0;
225   for(Int_t module= fFirstL1; module<=fLastL1;module++){
226     branch->GetEvent(module);
227     Int_t nrecp1 = ITSrec->GetEntries();
228     for(Int_t i=0; i<nrecp1;i++){
229       AliITSRecPoint *current = (AliITSRecPoint*)ITSrec->At(i);
230       lc[0]=current->GetX();
231       lc[2]=current->GetZ();
232       geom->LtoG(module,lc,gc);
233       zave+=gc[2];
234       zave2+=gc[2]*gc[2];
235       firipixe++;
236     }
237     fITS->ResetRecPoints();
238   }
239   if(firipixe>1){
240     rmszav=TMath::Sqrt(zave2/(firipixe-1)-zave*zave/firipixe/(firipixe-1));
241     zave=zave/firipixe;
242     if(fDebug>1)cout<<"++++++++++++++++++++++++++++++++++++++++++++++++++++++\n Number of firing pixels: "<<firipixe<<endl;
243   }
244   else {
245     fZFound = -200;
246     fZsig = -200;
247     Warning("FindVertexForCurrentEvent","No rec points on first layer for this event");
248     return fCurrentVertex;
249   }
250   Float_t zlim1=zave-rmszav;
251   Float_t zlim2=zave+rmszav;
252   Int_t sepa=(Int_t)((zlim2-zlim1)*10.+1.);
253   zlim2=zlim1 + sepa/10.;
254   TH1F *zvdis = new TH1F("z_ev","zv distr",sepa,zlim1,zlim2);
255   if(fDebug>0){
256     cout<<"Z limits: "<<zlim1<<" "<<zlim2<<"; Bins= "<<sepa<<endl;
257     cout<<"Bin width: "<<zvdis->GetBinWidth(1)*10000.<<" micron\n";
258   }
259   Int_t sizarr=100;
260   TArrayF *zval = new TArrayF(sizarr);
261   //  TArrayF *curv = new TArrayF(sizarr);
262   Int_t ncoinc=0;
263   for(Int_t module= fFirstL1; module<=fLastL1;module++){
264     if(fDebug>0)cout<<"processing module   "<<module<<"                  \r";
265     branch->GetEvent(module);
266     Int_t nrecp1 = ITSrec->GetEntries();
267     TObjArray *poiL1 = new TObjArray(nrecp1);
268     for(Int_t i=0; i<nrecp1;i++)poiL1->AddAt(ITSrec->At(i),i);
269     fITS->ResetRecPoints();
270     for(Int_t i=0; i<nrecp1;i++){
271       AliITSRecPoint *current = (AliITSRecPoint*)poiL1->At(i);
272       lc[0]=current->GetX();
273       lc[2]=current->GetZ();
274       geom->LtoG(module,lc,gc);
275       gc[0]-=fX0;
276       gc[1]-=fY0;
277       Float_t r1=TMath::Sqrt(gc[0]*gc[0]+gc[1]*gc[1]);
278       Float_t phi1 = TMath::ATan2(gc[1],gc[0]);
279       if(phi1<0)phi1=2*TMath::Pi()+phi1;
280       if(fDebug>1)cout<<"module "<<module<<" "<<gc[0]<<" "<<gc[1]<<" "<<gc[2]<<" "<<phi1<<"     \n";
281       for(Int_t modul2=fFirstL2; modul2<=fLastL2; modul2++){
282         branch->GetEvent(modul2);
283         Int_t nrecp2 = ITSrec->GetEntries();
284         for(Int_t j=0; j<nrecp2;j++){
285           AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(j);
286           lc2[0]=recp->GetX();
287           lc2[2]=recp->GetZ();
288           geom->LtoG(modul2,lc2,gc2);
289           gc2[0]-=fX0;
290           gc2[1]-=fY0;
291           Float_t r2=TMath::Sqrt(gc2[0]*gc2[0]+gc2[1]*gc2[1]);
292           Float_t zr0=(r2*gc[2]-r1*gc2[2])/(r2-r1);
293           Float_t phi2 = TMath::ATan2(gc2[1],gc2[0]);
294           if(phi2<0)phi2=2.*TMath::Pi()+phi2;
295           Float_t diff = TMath::Abs(phi2-phi1);
296           if(diff>TMath::Pi())diff=2.*TMath::Pi()-diff;
297           if(zr0>zlim1 && zr0<zlim2){
298             if(diff<fDiffPhiMax ){
299               zvdis->Fill(zr0);
300               zval->AddAt(zr0,ncoinc);
301               /* uncomment these lines to use curvature as a weight
302               Float_t cu = Curv(0.,0.,gc[0],gc[1],gc2[0],gc[1]);
303               curv->AddAt(cu,ncoinc);
304               */
305               ncoinc++;
306               if(ncoinc==(sizarr-1)){
307                 sizarr+=100;
308                 zval->Set(sizarr);
309                 //uncomment next line to use curvature as weight
310                 //              curv->Set(sizarr);
311               }
312             }
313           }
314         }
315         fITS->ResetRecPoints();
316       }
317     }
318     delete poiL1;
319   }         // loop on modules
320   if(fDebug>0){
321     cout<<endl<<"Number of coincidences = "<<ncoinc<<endl;
322   }
323   //  EvalZ(zvdis,sepa,ncoinc,zval,curv);
324   EvalZ(zvdis,sepa,ncoinc,zval);
325   delete zvdis;
326   delete zval;
327   //  delete curv;
328   if(fCurrentVertex){
329     char name[30];
330     sprintf(name,"Vertex_%d",evnumber);
331     //    fCurrentVertex->SetName(name);
332     fCurrentVertex->SetTitle("vertexer: PPZ");
333   }
334   return fCurrentVertex;
335 }
336
337 //______________________________________________________________________
338 void AliITSVertexerPPZ::FindVertices(){
339   // computes the vertices of the events in the range FirstEvent - LastEvent
340   AliRunLoader *rl = AliRunLoader::GetRunLoader();
341   AliITSLoader* ITSloader =  (AliITSLoader*) rl->GetLoader("ITSLoader");
342   ITSloader->ReloadRecPoints();
343   for(Int_t i=fFirstEvent;i<=fLastEvent;i++){
344     cout<<"Processing event "<<i<<endl;
345     rl->GetEvent(i);
346     FindVertexForCurrentEvent(i);
347     if(fCurrentVertex){
348       WriteCurrentVertex();
349     }
350     else {
351       if(fDebug>0){
352         cout<<"Vertex not found for event "<<i<<endl;
353         cout<<"fZFound = "<<fZFound<<", fZsig= "<<fZsig<<endl;
354       }
355     }
356   }
357 }
358
359 //________________________________________________________
360 void AliITSVertexerPPZ::PrintStatus() const {
361   // Print current status
362   cout <<"=======================================================\n";
363   cout <<" First layer first and last modules: "<<fFirstL1<<", ";
364   cout <<fLastL1<<endl;
365   cout <<" Second layer first and last modules: "<<fFirstL2<<", ";
366   cout <<fLastL2<<endl;
367   cout <<" Max Phi difference: "<<fDiffPhiMax<<endl;
368   cout <<" Window for Z search: "<<fWindow<<endl;
369   cout <<" Current Z "<<fZFound<<"; Z sig "<<fZsig<<endl;
370   cout <<" Debug flag: "<<fDebug<<endl;
371   cout<<"First event to be processed "<<fFirstEvent;
372   cout<<"\n Last event to be processed "<<fLastEvent<<endl;
373 }
374
375 //________________________________________________________
376 Float_t AliITSVertexerPPZ::Curv(Double_t x1,Double_t y1,
377                    Double_t x2,Double_t y2,
378                    Double_t x3,Double_t y3) 
379 {
380   //-----------------------------------------------------------------
381   // Initial approximation of the track curvature (Y. Belikov) squared
382   //-----------------------------------------------------------------
383   Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
384   Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
385                   (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
386   Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
387                   (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
388
389   Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
390
391   Float_t val = (Float_t)(( -xr*yr/TMath::Sqrt(xr*xr+yr*yr))
392                           *( -xr*yr/TMath::Sqrt(xr*xr+yr*yr)));
393   return val; 
394 }