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