]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSVertexer3DTapan.cxx
Savannah bug 59091 (A. Dainese)
[u/mrichter/AliRoot.git] / ITS / AliITSVertexer3DTapan.cxx
1 /**************************************************************************
2  * Copyright(c) 2006-2008, 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
16 //-----------------------------------------------------------------
17 //           AliITSVertexer3DTapan class
18 //   This is a class for the 3d vertex finding
19 //    Origin: Tapan Nayak, VECC-CERN, Tapan.Nayak@cern.ch
20 //-----------------------------------------------------------------
21
22 #include <TH1.h>
23 #include <TTree.h>
24 #include <TClonesArray.h>
25
26 #include <AliITSVertexer3DTapan.h>
27 #include <AliITSRecPoint.h>
28 #include <AliITSgeomTGeo.h>
29 #include <AliESDVertex.h>
30
31 ClassImp(AliITSVertexer3DTapan)
32
33 void AliITSVertexer3DTapan::LoadClusters(TTree *cTree) {
34   //--------------------------------------------------------------------
35   //This function loads the SPD clusters
36   //--------------------------------------------------------------------
37    TClonesArray dummy("AliITSRecPoint",10000), *clusters=&dummy;
38    TBranch *branch=cTree->GetBranch("ITSRecPoints");
39    branch->SetAddress(&clusters);
40
41    Int_t nentr=cTree->GetEntries(),nc1=0,nc2=0;
42    for (Int_t i=0; i<nentr; i++) {
43        if (!cTree->GetEvent(i)) continue;
44        //
45        //   Below:
46        //   "alpha" is the angle from the global X-axis to the
47        //    local GEANT X'-axis  ( rot[0]=cos(alpha) and rot[1]=sin(alpha) )
48        //    "phi" is the angle from the global X-axis to the
49        //         local cluster X"-axis
50        //
51
52        Int_t    lay,lad,det; AliITSgeomTGeo::GetModuleId(i,lay,lad,det);
53
54        if (lay>2) break;  //load the SPD clusters only
55
56        Int_t ncl=clusters->GetEntriesFast();
57        Float_t hPhi;
58        while (ncl--) {
59           AliITSRecPoint *c=(AliITSRecPoint*)clusters->UncheckedAt(ncl);
60           Float_t pos[3];
61           c->GetGlobalXYZ(pos);
62           if (lay==1) {
63             /*             fX1[nc1]= r*cp - c->GetY()*sp;
64              fY1[nc1]= r*sp + c->GetY()*cp;
65              fZ1[nc1]= c->GetZ(); */
66             fX1[nc1] = pos[0]; fY1[nc1] = pos[1]; fZ1[nc1] = pos[2];
67              CalculatePhi(fX1[nc1], fY1[nc1], hPhi);
68              fPhi1[nc1]= hPhi;
69              nc1++;
70           } else {
71             /* fX2[nc2]= r*cp - c->GetY()*sp;
72              fY2[nc2]= r*sp + c->GetY()*cp;
73              fZ2[nc2]= c->GetZ(); */
74             fX2[nc2] = pos[0]; fY2[nc2] = pos[1]; fZ2[nc2] = pos[2];
75              CalculatePhi(fX2[nc2], fY2[nc2], hPhi);
76              fPhi2[nc2]= hPhi;
77              nc2++;
78           }
79        }
80    }
81    ficlu1 = nc1; ficlu2 = nc2;
82    AliInfo(Form("Number of clusters: %d (first layer) and %d (second layer)",ficlu1,ficlu2));
83 }
84
85 AliESDVertex *AliITSVertexer3DTapan::FindVertexForCurrentEvent(TTree *cTree) {
86   //
87   // This function reconstructs ....
88   //
89   //
90   LoadClusters(cTree);
91
92   Double_t pos[3], postemp[3], sigpos[3];
93   Int_t ncontr, ncontrtemp;
94   Float_t cuts[3];
95   Int_t vtxstatus=0;
96
97   //....
98   pos[0] = 0.;   pos[1] = 0.;   pos[2] = 0.;
99   cuts[0]=1.;  cuts[1]=1.;  cuts[2]=20.;
100   CalculateVertex3d1(pos, cuts, ncontr);
101   if(ncontr==0) {
102     pos[0] = 9999.;   pos[1] = 9999.;   pos[2] = 9999.;
103     vtxstatus = -1;
104   }
105   AliInfo(Form("1st step: %d %f %f %f st=%d",ncontr,pos[0],pos[1],pos[2],vtxstatus));
106
107   if(vtxstatus == 0) {
108     ncontrtemp = ncontr; postemp[0] = pos[0];   postemp[1] = pos[1];   postemp[2] = pos[2]; 
109     cuts[0]=0.3;  cuts[1]=0.3;  cuts[2]=1.;
110     CalculateVertex3d1(pos, cuts, ncontr);
111     if(ncontr==0) {
112       ncontr = ncontrtemp; pos[0] = postemp[0];   pos[1] = postemp[1];   pos[2] = postemp[2];
113       vtxstatus = 2;
114     }
115     AliInfo(Form("2nd step: %d %f %f %f st=%d",ncontr,pos[0],pos[1],pos[2],vtxstatus));
116   }
117
118   if(vtxstatus == 0) {
119     ncontrtemp = ncontr; postemp[0] = pos[0];   postemp[1] = pos[1];   postemp[2] = pos[2]; 
120     cuts[0]=0.25;  cuts[1]=0.25;  cuts[2]=1.0;
121     CalculateVertex3d2(pos, cuts, ncontr, sigpos);
122     if(ncontr==0) {
123       ncontr = ncontrtemp; pos[0] = postemp[0];   pos[1] = postemp[1];   pos[2] = postemp[2];
124       vtxstatus = 3;
125     }
126     AliInfo(Form("3rd step: %d %f %f %f st=%d",ncontr,pos[0],pos[1],pos[2],vtxstatus));
127   }
128
129   if(vtxstatus == 0) {
130     ncontrtemp = ncontr; postemp[0] = pos[0];   postemp[1] = pos[1];   postemp[2] = pos[2]; 
131     cuts[0]=0.1;  cuts[1]=0.1;  cuts[2]=0.2;
132     CalculateVertex3d2(pos, cuts, ncontr, sigpos);
133     if(ncontr==0) {
134       ncontr = ncontrtemp; pos[0] = postemp[0];   pos[1] = postemp[1];   pos[2] = postemp[2];
135       vtxstatus = 4;
136     }
137     AliInfo(Form("4th step: %d %f %f %f st=%d",ncontr,pos[0],pos[1],pos[2],vtxstatus));
138   }
139   AliInfo(Form("Final step: %d %f %f %f st=%d",ncontr,pos[0],pos[1],pos[2],vtxstatus));
140
141   return new AliESDVertex(pos,sigpos,(Double_t)vtxstatus,ncontr,"AliITSVertexer3DTapan");
142
143 }
144
145
146 void AliITSVertexer3DTapan::CalculateVertex3d1(Double_t pos[3], Float_t cuts[3], Int_t &ncontr) {
147   //
148   // This function reconstructs first two steps of vertex
149   //
150
151   Double_t p1[4], p2[4], p3[4], p4[4];
152   Double_t pa[3], pb[3];
153   Double_t hphi1, hphi2, hphi3, hphi4;
154
155   ncontr = 0;
156   Float_t  phicut = 1.0;
157   Double_t distance;  Float_t  distancecut = 1.0;
158   Int_t    ibin=20;  Float_t  ilow=-1.; Float_t  ihigh=1.; 
159   Int_t    ibinz=400; Float_t  ilowz=-20.; Float_t  ihighz=20.;
160
161   TH1F *hx = new TH1F("hx","",  ibin, ilow, ihigh);
162   TH1F *hy = new TH1F("hy","",  ibin, ilow, ihigh);
163   TH1F *hz = new TH1F("hz","",  ibinz,ilowz,ihighz);
164   
165   for (Int_t ip1=0; ip1<ficlu1; ip1++) {
166     // Two points on layer1: p1 and p3
167     p1[0] = fX1[ip1];   p1[1] = fY1[ip1];    p1[2] = fZ1[ip1];   
168     p3[0] = fX1[ip1+1]; p3[1] = fY1[ip1+1];  p3[2] = fZ1[ip1+1]; 
169     hphi1 = fPhi1[ip1]; hphi3 = fPhi1[ip1+1];
170     
171     for (Int_t ip2=0; ip2<ficlu2; ip2++) {
172       // Two points on layer 2: p2 and p4
173       p2[0] = fX2[ip2];   p2[1] = fY2[ip2];     p2[2] = fZ2[ip2]; 
174       p4[0] = fX2[ip2+1]; p4[1] = fY2[ip2+1];   p4[2] = fZ2[ip2+1]; 
175       hphi2 = fPhi2[ip2]; hphi4 = fPhi2[ip2+1];
176
177       // First line is formed by p1-p2 and second line by p3-p4
178       // We find two points on each line which form the closest distance of the two lines
179       // pa[0],pa[1],pa[2]: points on line 1 and pb[0],pb[1],pb[2]: points on line 2
180       // Next: Consider x, y and z to be less than cuts[0], cuts[1] and cuts[2], respectively
181
182       if(TMath::Abs(hphi1-hphi2)<phicut && TMath::Abs(hphi3-hphi4)<phicut){
183         CalculateLine(p1, p2, p3, p4, pa, pb);
184
185         if (pa[0]>pos[0]-cuts[0] && pa[0]<pos[0]+cuts[0] && pa[1]>pos[1]-cuts[1] && pa[1]<pos[1]+cuts[1] && pa[2]>pos[2]-cuts[2] && pa[2]<pos[2]+cuts[2]){
186           distance = (TMath::Sqrt(pow((pa[0]-pb[0]),2) + pow((pa[1]-pb[1]),2) + pow((pa[2]-pb[2]),2)));
187           if(distance<distancecut){
188             hx->Fill(pa[0]);            hy->Fill(pa[1]);                hz->Fill(pa[2]);
189             hx->Fill(pb[0]);            hy->Fill(pb[1]);                hz->Fill(pb[2]);
190             ncontr++;
191           }
192         }
193       }
194
195       // Third line is formed by p1-p4 and fourth line by p3-p2
196       // We find two points on each line which form the closest distance of the two lines
197       // pa[0],pa[1],pa[2]: points on line 3 and pb[0],pb[1],pb[2]: points on line 4
198       // Next: Consider x, y and z to be less than cuts[0], cuts[1] and cuts[2], respectively      
199       if(TMath::Abs(hphi1-hphi4)<phicut && TMath::Abs(hphi3-hphi2)<phicut) {
200         CalculateLine(p1, p4, p3, p2, pa, pb);
201         if (pa[0]>pos[0]-cuts[0] && pa[0]<pos[0]+cuts[0] && pa[1]>pos[1]-cuts[1] && pa[1]<pos[1]+cuts[1]){
202           distance = (TMath::Sqrt(pow((pa[0]-pb[0]),2) + pow((pa[1]-pb[1]),2) + pow((pa[2]-pb[2]),2)));
203           if(distance<distancecut){
204             hx->Fill(pa[0]);            hy->Fill(pa[1]);                hz->Fill(pa[2]);
205             hx->Fill(pb[0]);            hy->Fill(pb[1]);                hz->Fill(pb[2]);
206             ncontr++;
207           }
208         }
209       }
210     }
211   }
212
213   Int_t maxbinx = hx->GetMaximumBin(); 
214   Int_t maxbiny = hy->GetMaximumBin();
215   Int_t maxbinz = hz->GetMaximumBin();
216   pos[0] = ilow  + ((ihigh-ilow)/ibin)*maxbinx;
217   pos[1] = ilow  + ((ihigh-ilow)/ibin)*maxbiny;
218   pos[2] = ilowz + ((ihighz-ilowz)/ibinz)*maxbinz;
219   hx->Delete();
220   hy->Delete();
221   hz->Delete();
222 }
223
224 void AliITSVertexer3DTapan::CalculateVertex3d2(Double_t pos[3], Float_t cuts[3], Int_t &ncontr, Double_t sigpos[3]) {
225   //
226   // This function reconstructs second two steps of vertex
227   //
228
229   Double_t p1[4], p2[4], p3[4], p4[4];
230   Double_t pa[3], pb[3];
231   Double_t hphi1, hphi2, hphi3, hphi4;
232
233   ncontr = 0;
234   Float_t  phicut = 0.3;
235   Double_t distance; Float_t  distancecut = 1.0;
236
237   Double_t vertx  =0.;   Double_t verty  =0.;   Double_t vertz  =0.;     
238   Double_t vertx2 =0.;   Double_t verty2 =0.;   Double_t vertz2 =0.;     
239
240   for (Int_t ip1=0; ip1<ficlu1; ip1++) {
241     // Two points on layer1: p1 and p3
242     p1[0] = fX1[ip1];   p1[1] = fY1[ip1];    p1[2] = fZ1[ip1];   
243     p3[0] = fX1[ip1+1]; p3[1] = fY1[ip1+1];  p3[2] = fZ1[ip1+1]; 
244     hphi1 = fPhi1[ip1]; hphi3 = fPhi1[ip1+1];
245     
246     for (Int_t ip2=0; ip2<ficlu2; ip2++) {
247       // Two points on layer 2: p2 and p4
248       p2[0] = fX2[ip2];   p2[1] = fY2[ip2];     p2[2] = fZ2[ip2]; 
249       p4[0] = fX2[ip2+1]; p4[1] = fY2[ip2+1];   p4[2] = fZ2[ip2+1]; 
250       hphi2 = fPhi2[ip2]; hphi4 = fPhi2[ip2+1];
251
252       // First line is formed by p1-p2 and second line by p3-p4
253       // We find two points on each line which form the closest distance of the two lines
254       // pa[0],pa[1],pa[2] are the points on line 1 and pb[0],pb[1],pb[2] are the points on line 2
255       // Next: Consider x, y and z to be less than cuts[0], cuts[1] and cuts[2], respectively
256
257       if(TMath::Abs(hphi1-hphi2)<phicut && TMath::Abs(hphi3-hphi4)<phicut){
258         CalculateLine(p1, p2, p3, p4, pa, pb);
259
260         // We consider the points where x, y and z points are less than xcut, ycut and zcut, respectively
261         if (pa[0]>pos[0]-cuts[0] && pa[0]<pos[0]+cuts[0] && pa[1]>pos[1]-cuts[1] && pa[1]<pos[1]+cuts[1] && pa[2]>pos[2]-cuts[2] && pa[2]<pos[2]+cuts[2]){
262           distance = (TMath::Sqrt(pow((pa[0]-pb[0]),2) + pow((pa[1]-pb[1]),2) + pow((pa[2]-pb[2]),2)));
263           if(distance<distancecut){
264             ncontr++;
265             vertx  = vertx  + pa[0];       verty  = verty  + pa[1];       vertz  = vertz  + pa[2];
266             vertx2 = vertx2 + pa[0]*pa[0]; verty2 = verty2 + pa[1]*pa[1]; vertz2 = vertz2 + pa[2]*pa[2];
267             ncontr++;
268             vertx  = vertx  + pb[0];       verty  = verty  + pb[1];       vertz  = vertz  + pb[2];
269             vertx2 = vertx2 + pb[0]*pb[0]; verty2 = verty2 + pb[1]*pb[1]; vertz2 = vertz2 + pb[2]*pb[2];
270           }
271         }
272       }
273       
274       // Third line is formed by p1-p4 and fourth line by p3-p2
275       // We find two points on each line which form the closest distance of the two lines
276       // pa[0],pa[1],pa[2] are the points on line 3 and pb[0],pb[1],pb[2] are the points on line 4
277       // Next: Consider x, y and z to be less than cuts[0], cuts[1] and cuts[2], respectively      
278       if(TMath::Abs(hphi1-hphi4)<phicut && TMath::Abs(hphi3-hphi2)<phicut) {
279         
280         CalculateLine(p1, p4, p3, p2, pa, pb);
281         if (pa[0]>pos[0]-cuts[0] && pa[0]<pos[0]+cuts[0] && pa[1]>pos[1]-cuts[1] && pa[1]<pos[1]+cuts[1] && pa[2]>pos[2]-cuts[2] && pa[2]<pos[2]+cuts[2]){
282           distance = (TMath::Sqrt(pow((pa[0]-pb[0]),2) + pow((pa[1]-pb[1]),2) + pow((pa[2]-pb[2]),2)));
283           if(distance<distancecut){
284             ncontr++;
285             vertx  = vertx  + pa[0];       verty  = verty  + pa[1];       vertz  = vertz  + pa[2];
286             vertx2 = vertx2 + pa[0]*pa[0]; verty2 = verty2 + pa[1]*pa[1]; vertz2 = vertz2 + pa[2]*pa[2];
287             ncontr++;
288             vertx  = vertx  + pb[0];       verty  = verty  + pb[1];       vertz  = vertz  + pb[2];
289             vertx2 = vertx2 + pb[0]*pb[0]; verty2 = verty2 + pb[1]*pb[1]; vertz2 = vertz2 + pb[2]*pb[2];
290           }
291         }
292       }
293     }
294   }
295
296   if(ncontr>0){
297     pos[0]  = vertx/ncontr;  pos[1] = verty/ncontr;   pos[2] = vertz/ncontr;
298     vertx2 = vertx2/ncontr; verty2 = verty2/ncontr; vertz2 = vertz2/ncontr;
299     sigpos[0] = TMath::Sqrt(vertx2 - pos[0]*pos[0]); 
300     sigpos[1] = TMath::Sqrt(verty2 - pos[1]*pos[1]);
301     sigpos[2] = TMath::Sqrt(vertz2 - pos[2]*pos[2]);
302   }
303   ncontr = ncontr/2;
304 }
305
306 void AliITSVertexer3DTapan::CalculatePhi(Float_t fx, Float_t fy, Float_t & phi)
307 {
308   //calculates phi
309   Float_t ybyx, phi1;
310   const Float_t kradian = 180./3.141592654;
311   
312   if(fx==0.)
313     {
314       if(fy>0.) phi = 90.;
315       if(fy<0.) phi = 270.;
316     }
317   if(fx != 0.)
318     {
319       ybyx = fy/fx;
320       if(ybyx < 0) ybyx = - ybyx;
321       phi1 = TMath::ATan(ybyx)*kradian;
322       if(fx > 0 && fy > 0) phi = phi1;        // 1st Quadrant
323       if(fx < 0 && fy > 0) phi = 180 - phi1;  // 2nd Quadrant
324       if(fx < 0 && fy < 0) phi = 180 + phi1;  // 3rd Quadrant
325       if(fx > 0 && fy < 0) phi = 360 - phi1;  // 4th Quadrant
326       
327     }
328   phi = phi/kradian;
329 }
330
331 void AliITSVertexer3DTapan::CalculateLine(Double_t p1[4], Double_t p2[4], Double_t p3[4], Double_t p4[4], Double_t pa[3], Double_t pb[3]) const{
332   //calculates line
333   Double_t p13x, p13y, p13z;
334   Double_t p21x, p21y, p21z;
335   Double_t p43x, p43y, p43z;
336   Double_t d1343, d4321, d1321, d4343, d2121;
337   Double_t numer, denom;
338   Double_t mua,   mub;
339   mua = 0.; mub = 0.;
340   
341   p13x = p1[0] - p3[0];
342   p13y = p1[1] - p3[1];
343   p13z = p1[2] - p3[2];
344
345   p21x = p2[0] - p1[0];
346   p21y = p2[1] - p1[1];
347   p21z = p2[2] - p1[2];
348
349   p43x = p4[0] - p3[0];
350   p43y = p4[1] - p3[1];
351   p43z = p4[2] - p3[2];
352
353   d1343 = p13x * p43x + p13y * p43y + p13z * p43z;
354   d4321 = p43x * p21x + p43y * p21y + p43z * p21z;
355   d1321 = p13x * p21x + p13y * p21y + p13z * p21z;
356   d4343 = p43x * p43x + p43y * p43y + p43z * p43z;
357   d2121 = p21x * p21x + p21y * p21y + p21z * p21z;
358
359   denom = d2121 * d4343 - d4321 * d4321;
360   numer = d1343 * d4321 - d1321 * d4343;
361
362   if(denom>0) mua = numer / denom;
363   if(d4343>0) mub = (d1343 + d4321 * (mua)) / d4343;
364
365   pa[0] = p1[0] + mua * p21x;
366   pa[1] = p1[1] + mua * p21y;
367   pa[2] = p1[2] + mua * p21z;
368
369   pb[0] = p3[0] + mub * p43x;
370   pb[1] = p3[1] + mub * p43y;
371   pb[2] = p3[2] + mub * p43z;
372 }
373