]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCLaserTrack.cxx
Using the new dEdx algorithm + signal below threshold
[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   TObjArray *arrLaserTracks = 0x0;
181   
182   AliCDBManager *man=AliCDBManager::Instance();
183   if (!man->GetDefaultStorage() && gSystem->Getenv("ALICE_ROOT")) man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
184   if (man->GetDefaultStorage()){
185     if (man->GetRun()<0) man->SetRun(0);
186     AliCDBEntry *entry=man->Get(AliCDBPath("TPC/Calib/LaserTracks"));
187     arrLaserTracks = (TObjArray*)entry->GetObject();
188     entry->SetOwner(kTRUE);
189   } else {
190     if (!gSystem->AccessPathName("LaserTracks.root")){
191       TFile f("LaserTracks.root");
192       arrLaserTracks=(TObjArray*)f.Get("arrLaserTracks");
193       f.Close();
194     }
195   }
196   if ( !arrLaserTracks ) {
197 //      AliWarning(Form("Could not get laser position data from file: '%s'",fgkDataFileName));
198     return;
199   }
200   
201   arrLaserTracks->SetOwner();
202   
203   fgArrLaserTracks = new TObjArray(fgkNLaserTracks);
204   fgArrLaserTracks->SetOwner();
205   for (Int_t itrack=0; itrack<fgkNLaserTracks; itrack++){
206     AliTPCLaserTrack *ltr = (AliTPCLaserTrack*)arrLaserTracks->At(itrack);
207     if ( !ltr ){
208 //          AliWarning(Form("No informatino found for Track %d!",itrack));
209       continue;
210     }
211     ltr->UpdatePoints();
212     fgArrLaserTracks->AddAt(new AliTPCLaserTrack(*ltr),itrack);
213   }
214
215   delete arrLaserTracks;
216 }
217
218
219 void AliTPCLaserTrack::UpdatePoints(){
220   //
221   // update track points
222   //
223   const Double_t kMaxSnp=0.97;
224   AliTPCROC* roc = AliTPCROC::Instance();
225   //
226   //
227   if (!fVecSec){
228     fVecSec=new TVectorD(159);
229     fVecP2 =new TVectorD(159);       //                - P2  
230     fVecPhi=new TVectorD(159);       //                - Phi
231     fVecGX=new TVectorD(159);       // points vectors - globalX
232     fVecGY=new TVectorD(159);       // points vectors - globalY
233     fVecGZ=new TVectorD(159);       // points vectors - globalZ
234     fVecLX=new TVectorD(159);       // points vectors - localX
235     fVecLY=new TVectorD(159);       // points vectors - localY
236     fVecLZ=new TVectorD(159);       // points vectors - localZ
237
238   }
239   for (Int_t irow=158; irow>=0; irow--){
240     (*fVecSec)[irow]= -1;       //                -
241     (*fVecP2)[irow] = 0;       //                - P2  -snp
242     (*fVecPhi)[irow]= 0;       //                - global phi
243     (*fVecGX)[irow] = 0;       // points vectors - globalX
244     (*fVecGY)[irow] = 0;       // points vectors - globalY
245     (*fVecGZ)[irow] = 0;       // points vectors - globalZ
246     (*fVecLX)[irow] = 0;       // points vectors - localX
247     (*fVecLY)[irow] = 0;       // points vectors - localY
248     (*fVecLZ)[irow] = 0;       // points vectors - localZ
249
250   }
251   Double_t gxyz[3];
252   Double_t lxyz[3];
253   AliTPCLaserTrack*ltrp=new AliTPCLaserTrack(*this);  //make temporary track
254
255   for (Int_t irow=158; irow>=0; irow--){
256     UInt_t srow = irow;
257     Int_t sector=0;
258    
259     if (srow >=roc->GetNRows(0)) {
260       srow-=roc->GetNRows(0);
261       sector=36    ;
262     }
263     lxyz[0]= roc->GetPadRowRadii(sector,srow);
264     if (!ltrp->PropagateTo(lxyz[0],5)) break;
265     ltrp->GetXYZ(gxyz);
266     //
267     Double_t alpha=TMath::ATan2(gxyz[1],gxyz[0]);
268     if (alpha<0) alpha+=2*TMath::Pi();
269     sector      +=TMath::Nint(-0.5+9*alpha/TMath::Pi());
270     if (gxyz[2]<0) sector+=18;
271     Double_t salpha   = TMath::Pi()*(sector+0.5)/9.;    
272     if (!ltrp->Rotate(salpha)) break;
273     if (!ltrp->PropagateTo(lxyz[0],5)) break;
274     if (TMath::Abs(ltrp->GetSnp())>kMaxSnp) break;
275     ltrp->GetXYZ(gxyz);
276     lxyz[1]=ltrp->GetY();
277     lxyz[2]=ltrp->GetZ();
278     (*fVecSec)[irow]= sector;
279     (*fVecP2)[irow] = ltrp->GetSnp();                 //                - P2  -snp
280     (*fVecPhi)[irow]= TMath::ATan2(gxyz[1],gxyz[0]);  //                - global phi
281     (*fVecGX)[irow] = gxyz[0];       // points vectors - globalX
282     (*fVecGY)[irow] = gxyz[1];       // points vectors - globalY
283     (*fVecGZ)[irow] = gxyz[2];       // points vectors - globalZ
284     (*fVecLX)[irow] = lxyz[0];       // points vectors - localX
285     (*fVecLY)[irow] = lxyz[1];       // points vectors - localY
286     (*fVecLZ)[irow] = lxyz[2];       // points vectors - localZ
287
288   }
289   delete ltrp;  // delete temporary track
290 }
291
292 Int_t AliTPCLaserTrack::IdentifyTrack(AliExternalTrackParam *track, Int_t side)
293 {
294   //
295   // Find the laser track which is corresponding closest to 'track'
296   // return its id
297   //
298   // 
299   const  Float_t   kMaxdphi=0.2;
300   const  Float_t   kMaxdphiP=0.05;
301   const  Float_t   kMaxdz=40;
302
303   if ( !fgArrLaserTracks ) LoadTracks();
304   TObjArray *arrTracks = GetTracks();
305   Double_t lxyz0[3];
306   Double_t lxyz1[3];
307   Double_t pxyz0[3];
308   Double_t pxyz1[3];
309   track->GetXYZ(lxyz0);
310   track->GetDirection(pxyz0);
311   //
312   Float_t mindist=10; // maxima minimal distance
313   Int_t id = -1;
314   AliExternalTrackParam*  ltr0= (AliExternalTrackParam*)arrTracks->UncheckedAt(0);
315   for (Int_t itrack=0; itrack<fgkNLaserTracks; itrack++){    
316     AliTPCLaserTrack *ltr = (AliTPCLaserTrack*)arrTracks->UncheckedAt(itrack);
317     if (side>=0) if (ltr->GetSide()!=side) continue;
318     Double_t * kokot = (Double_t*)ltr->GetParameter();
319     kokot[4]=-0.0000000001;
320     //
321     ltr->GetXYZ(lxyz1);
322     if (TMath::Abs(lxyz1[2]-lxyz0[2])>kMaxdz) continue;
323     // phi position
324     Double_t phi0 = TMath::ATan2(lxyz0[1],lxyz0[0]);
325     Double_t phi1 = TMath::ATan2(lxyz1[1],lxyz1[0]);
326     if (TMath::Abs(phi0-phi1)>kMaxdphi) continue;
327     // phi direction
328     ltr->GetDirection(pxyz1);
329     Float_t direction= pxyz0[0]*pxyz1[0] + pxyz0[1]*pxyz1[1] + pxyz0[2]*pxyz1[2];
330     Float_t distdir = (1-TMath::Abs(direction))*90.; //distance at entrance
331     if (1-TMath::Abs(direction)>kMaxdphiP)
332       continue;
333     //
334     Float_t dist=0;
335     dist+=TMath::Abs(lxyz1[0]-lxyz0[0]);
336     dist+=TMath::Abs(lxyz1[1]-lxyz0[1]);
337     //    dist+=TMath::Abs(lxyz1[2]-lxyz0[2]); //z is not used for distance calculation
338     dist+=distdir;
339     //    
340     if (id<0)  {
341       id =itrack; 
342       mindist=dist; 
343       ltr0=ltr;
344       continue;
345     }
346     if (dist>mindist) continue;
347     id = itrack;
348     mindist=dist;
349     ltr0=ltr;
350   }
351   return id;
352 }
353