]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCLaserTrack.cxx
AliTPCLaserTrack.cxx - max snp cut added
[u/mrichter/AliRoot.git] / TPC / AliTPCLaserTrack.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 #include <TObjArray.h>
18 #include <TFile.h>
19 #include <TString.h>
20 #include <TSystem.h>
21
22 #include "AliLog.h"
23 #include "AliTPCLaserTrack.h"
24 #include "AliTPCROC.h"
25
26 ClassImp(AliTPCLaserTrack)
27
28 TObjArray *AliTPCLaserTrack::fgArrLaserTracks=0x0;
29
30 AliTPCLaserTrack::AliTPCLaserTrack() :
31   AliExternalTrackParam(),
32   fId(-1),
33   fSide(-1),
34   fRod(-1),
35   fBundle(-1),
36   fBeam(-1),
37   fRayLength(0),
38   fVecSec(0),       // points vectors - sector
39   fVecP2(0),       // points vectors - snp
40   fVecPhi(0),       // points vectors - global phi
41   fVecGX(0),       // points vectors - globalX
42   fVecGY(0),       // points vectors - globalY
43   fVecGZ(0),       // points vectors - globalZ
44   fVecLX(0),       // points vectors - localX
45   fVecLY(0),       // points vectors - localY
46   fVecLZ(0)        // points vectors - localZ
47 {
48   //
49   // Default constructor
50   //
51
52 }
53
54 AliTPCLaserTrack::AliTPCLaserTrack(const AliTPCLaserTrack &ltr) :
55   AliExternalTrackParam(ltr),
56   fId(ltr.fId),
57   fSide(ltr.fSide),
58   fRod(ltr.fRod),
59   fBundle(ltr.fBundle),
60   fBeam(ltr.fBeam),
61   fRayLength(ltr.fRayLength),
62   fVecSec(0),       // points vectors - sector
63   fVecP2(0),       // points vectors - snp
64   fVecPhi(0),       // points vectors - global phi
65   fVecGX(0),       // points vectors - globalX
66   fVecGY(0),       // points vectors - globalY
67   fVecGZ(0),       // points vectors - globalZ
68   fVecLX(0),       // points vectors - localX
69   fVecLY(0),       // points vectors - localY
70   fVecLZ(0)        // points vectors - localZ
71 {
72   //
73   // Default constructor
74   //
75   fVecSec=new TVectorD(*ltr.fVecSec);       // points vectors - sector
76   fVecP2 =new TVectorD(*ltr.fVecP2);       // points vectors - snp
77   fVecPhi=new TVectorD(*ltr.fVecPhi);       // points vectors - global phi
78   fVecGX =new TVectorD(*ltr.fVecGX);       // points vectors - globalX
79   fVecGY =new TVectorD(*ltr.fVecGY);       // points vectors - globalY
80   fVecGZ =new TVectorD(*ltr.fVecGZ);       // points vectors - globalZ
81   fVecLX =new TVectorD(*ltr.fVecLX);       // points vectors - localX
82   fVecLY =new TVectorD(*ltr.fVecLY);       // points vectors - localY
83   fVecLZ =new TVectorD(*ltr.fVecLZ);       // points vectors - localY
84
85 }
86
87 AliTPCLaserTrack::AliTPCLaserTrack(const Int_t id, const Int_t side, const Int_t rod,
88                      const Int_t bundle, const Int_t beam,
89                      Double_t x, Double_t alpha,
90                      const Double_t param[5],
91                      const Double_t covar[15], const Float_t rayLength) :
92   AliExternalTrackParam(x,alpha,param,covar),
93   fId(id),
94   fSide(side),
95   fRod(rod),
96   fBundle(bundle),
97   fBeam(beam),
98   fRayLength(rayLength),
99   fVecSec(new TVectorD(159)),       // points vectors - sector
100   fVecP2(new TVectorD(159)),       // points vectors - snp
101   fVecPhi(new TVectorD(159)),       // points vectors - global phi
102   fVecGX(new TVectorD(159)),       // points vectors - globalX
103   fVecGY(new TVectorD(159)),       // points vectors - globalY
104   fVecGZ(new TVectorD(159)),       // points vectors - globalZ
105   fVecLX(new TVectorD(159)),       // points vectors - localX
106   fVecLY(new TVectorD(159)),       // points vectors - localY
107   fVecLZ(new TVectorD(159))        // points vectors - localZ
108
109 {
110   //
111   // create laser track from arguments
112   //
113   
114 }
115 //_____________________________________________________________________
116 AliTPCLaserTrack& AliTPCLaserTrack::operator = (const  AliTPCLaserTrack &source)
117 {
118   //
119   // assignment operator
120   //
121   if (&source == this) return *this;
122   new (this) AliTPCLaserTrack(source);
123   
124   return *this;
125 }
126
127
128 AliTPCLaserTrack::~AliTPCLaserTrack(){
129   //
130   // destructor
131   //
132   delete fVecSec;      //                - sector numbers  
133   delete fVecP2;       //                - P2  
134   delete fVecPhi;       // points vectors - global phi
135   delete fVecGX;       // points vectors - globalX
136   delete fVecGY;       // points vectors - globalY
137   delete fVecGZ;       // points vectors - globalZ
138   delete fVecLX;       // points vectors - localX
139   delete fVecLY;       // points vectors - localY
140   delete fVecLZ;       // points vectors - localZ
141 }
142
143 void AliTPCLaserTrack::LoadTracks()
144 {
145     //
146     // Load all design positions from file into the static array fgArrLaserTracks
147     //
148
149     if ( fgArrLaserTracks ) return;
150
151     TString dataFileName("$ALICE_ROOT/TPC/Calib/LaserTracks.root");  //Path to the Data File
152
153     TFile *f=TFile::Open(gSystem->ExpandPathName(dataFileName.Data()));
154     if ( !f || !f->IsOpen() ){
155 //      AliWarning(Form("Could not open laser data file: '%s'",dataFileName.Data()));
156 //      AliWarning("Could not open laser data file");
157         return;
158     }
159     TObjArray *arrLaserTracks = (TObjArray*)f->Get("arrLaserTracks");
160     if ( !arrLaserTracks ) {
161 //      AliWarning(Form("Could not get laser position data from file: '%s'",fgkDataFileName));
162         return;
163     }
164
165     fgArrLaserTracks = new TObjArray(fgkNLaserTracks);
166     for (Int_t itrack=0; itrack<fgkNLaserTracks; itrack++){
167         AliTPCLaserTrack *ltr = (AliTPCLaserTrack*)arrLaserTracks->At(itrack);
168         if ( !ltr ){
169 //          AliWarning(Form("No informatino found for Track %d!",itrack));
170             continue;
171         }
172         ltr->UpdatePoints();
173         fgArrLaserTracks->AddAt(new AliTPCLaserTrack(*ltr),itrack);
174     }
175     delete f;
176 }
177
178
179 void AliTPCLaserTrack::UpdatePoints(){
180   //
181   // update track points
182   //
183   const Double_t kMaxSnp=0.99;
184   AliTPCROC* roc = AliTPCROC::Instance();
185   //
186   //
187   if (!fVecSec){
188     fVecSec=new TVectorD(159);
189     fVecP2 =new TVectorD(159);       //                - P2  
190     fVecPhi=new TVectorD(159);       //                - Phi
191     fVecGX=new TVectorD(159);       // points vectors - globalX
192     fVecGY=new TVectorD(159);       // points vectors - globalY
193     fVecGZ=new TVectorD(159);       // points vectors - globalZ
194     fVecLX=new TVectorD(159);       // points vectors - localX
195     fVecLY=new TVectorD(159);       // points vectors - localY
196     fVecLZ=new TVectorD(159);       // points vectors - localZ
197
198   }
199   for (Int_t irow=158; irow>=0; irow--){
200     (*fVecSec)[irow]= 0;       //                - 
201     (*fVecP2)[irow] = 0;       //                - P2  -snp
202     (*fVecPhi)[irow]= 0;       //                - global phi
203     (*fVecGX)[irow] = 0;       // points vectors - globalX
204     (*fVecGY)[irow] = 0;       // points vectors - globalY
205     (*fVecGZ)[irow] = 0;       // points vectors - globalZ
206     (*fVecLX)[irow] = 0;       // points vectors - localX
207     (*fVecLY)[irow] = 0;       // points vectors - localY
208     (*fVecLZ)[irow] = 0;       // points vectors - localZ
209
210   }
211   Double_t gxyz[3];
212   Double_t lxyz[3];
213   AliTPCLaserTrack*ltrp=new AliTPCLaserTrack(*this);  //make temporary track
214
215   for (Int_t irow=158; irow>=0; irow--){
216     UInt_t srow = irow;
217     Int_t sector=0;
218    
219     if (srow >=roc->GetNRows(0)) {
220       srow-=roc->GetNRows(0);
221       sector=36    ;
222     }
223     lxyz[0]= roc->GetPadRowRadii(sector,srow);
224     if (!ltrp->PropagateTo(lxyz[0],5)) break;
225     ltrp->GetXYZ(gxyz);
226     //
227     Double_t alpha=TMath::ATan2(gxyz[1],gxyz[0]);
228     if (alpha<0) alpha+=2*TMath::Pi();
229     sector      +=TMath::Nint(-0.5+9*alpha/TMath::Pi());
230     if (gxyz[2]<0) sector+=18;
231     Double_t salpha   = TMath::Pi()*(sector+0.5)/9.;    
232     if (!ltrp->Rotate(salpha)) break;
233     if (!ltrp->PropagateTo(lxyz[0],5)) break;
234     if (TMath::Abs(ltrp->GetSnp())>kMaxSnp) break;
235     ltrp->GetXYZ(gxyz);
236     lxyz[1]=ltrp->GetY();
237     lxyz[2]=ltrp->GetZ();
238     (*fVecSec)[irow]= sector;
239     (*fVecP2)[irow] = ltrp->GetSnp();                 //                - P2  -snp
240     (*fVecPhi)[irow]= TMath::ATan2(gxyz[1],gxyz[0]);  //                - global phi
241     (*fVecGX)[irow] = gxyz[0];       // points vectors - globalX
242     (*fVecGY)[irow] = gxyz[1];       // points vectors - globalY
243     (*fVecGZ)[irow] = gxyz[2];       // points vectors - globalZ
244     (*fVecLX)[irow] = lxyz[0];       // points vectors - localX
245     (*fVecLY)[irow] = lxyz[1];       // points vectors - localY
246     (*fVecLZ)[irow] = lxyz[2];       // points vectors - localZ
247
248   }
249   delete ltrp;  // delete temporary track
250 }
251
252 Int_t AliTPCLaserTrack::IdentifyTrack(AliExternalTrackParam *track, Int_t side)
253 {
254   //
255   // Find the laser track which is corresponding closest to 'track'
256   // return its id
257   //
258   // 
259   const  Float_t   kMaxdphi=0.2;
260   const  Float_t   kMaxdphiP=0.05;
261   const  Float_t   kMaxdz=40;
262
263   if ( !fgArrLaserTracks ) LoadTracks();
264   TObjArray *arrTracks = GetTracks();
265   Double_t lxyz0[3];
266   Double_t lxyz1[3];
267   Double_t pxyz0[3];
268   Double_t pxyz1[3];
269   track->GetXYZ(lxyz0);
270   track->GetDirection(pxyz0);
271   //
272   Float_t mindist=10; // maxima minimal distance
273   Int_t id = -1;
274   AliExternalTrackParam*  ltr0= (AliExternalTrackParam*)arrTracks->UncheckedAt(0);
275   for (Int_t itrack=0; itrack<fgkNLaserTracks; itrack++){    
276     AliTPCLaserTrack *ltr = (AliTPCLaserTrack*)arrTracks->UncheckedAt(itrack);
277     if (side>=0) if (ltr->GetSide()!=side) continue;
278     Double_t * kokot = (Double_t*)ltr->GetParameter();
279     kokot[4]=-0.0000000001;
280     //
281     ltr->GetXYZ(lxyz1);
282     if (TMath::Abs(lxyz1[2]-lxyz0[2])>kMaxdz) continue;
283     // phi position
284     Double_t phi0 = TMath::ATan2(lxyz0[1],lxyz0[0]);
285     Double_t phi1 = TMath::ATan2(lxyz1[1],lxyz1[0]);
286     if (TMath::Abs(phi0-phi1)>kMaxdphi) continue;
287     // phi direction
288     ltr->GetDirection(pxyz1);
289     Float_t direction= pxyz0[0]*pxyz1[0] + pxyz0[1]*pxyz1[1] + pxyz0[2]*pxyz1[2];
290     Float_t distdir = (1-TMath::Abs(direction))*90.; //distance at entrance
291     if (1-TMath::Abs(direction)>kMaxdphiP)
292       continue;
293     //
294     Float_t dist=0;
295     dist+=TMath::Abs(lxyz1[0]-lxyz0[0]);
296     dist+=TMath::Abs(lxyz1[1]-lxyz0[1]);
297     //    dist+=TMath::Abs(lxyz1[2]-lxyz0[2]); //z is not used for distance calculation
298     dist+=distdir;
299     //    
300     if (id<0)  {
301       id =itrack; 
302       mindist=dist; 
303       ltr0=ltr;
304       continue;
305     }
306     if (dist>mindist) continue;
307     id = itrack;
308     mindist=dist;
309     ltr0=ltr;
310   }
311   return id;
312 }
313