]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCLaserTracks.cxx
Tree chaching implemented and update of merging of trees
[u/mrichter/AliRoot.git] / TPC / AliTPCLaserTracks.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 //                                                                           //
19 //  This class is a container class for the laser beam positions of          //
20 //  the TPC laser system (ALICE-INT-2002-022  v.1 ).                         //
21 //  The system consits of 6 Laser Rods, each contains four mirror bundels    //
22 //  on each side of the TPC. Each bundle has seven mirrors.                  //
23 //                                                                           //
24 //  The TPC side 0 corresponds to the "Shaft Side (A)", side 1 corresponds   //
25 //  to the "Muon Side (C)".
26 //  The laser rods are counted counter clockwise starting with the one       //
27 //  with the smalles phi angle.                                              //
28 //  The mirror bundles in one rod for each side are counted starting from    //
29 //  form the greatest z value to the smalles.                                //
30 //  The beams in each bundel are counted counter clockwise. Having all beams //
31 //  pointing downward, the first beam will be the leftmost.                  //
32 //                                                                           //
33 ///////////////////////////////////////////////////////////////////////////////
34
35 /*
36  Create File with design positions
37
38  .L AliTPCLaserTrack.cxx+
39  AliTPCLaserTracks ltracks
40  ltracks.WriteTreeDesignData()
41  TFile f("LaserTracksDesing.root");
42  TTree * tree =(TTree*)f.Get("LaserTracks");
43  
44 */
45
46 #include <Riostream.h>
47 #include <TString.h>
48 #include <TPolyLine3D.h>
49 #include <TMath.h>
50 #include <TTree.h>
51 #include <TFile.h>
52 #include <TEventList.h>
53 #include <TVector2.h>
54 #include <TVector3.h>
55 #include <TROOT.h>
56
57 #include <TGraph.h>
58
59
60 #include "AliTPCLaserTracks.h"
61
62 ////////////////////////////////////////////////////////////////////////
63 //              Class AliTPCLaserTracks
64 ////////////////////////////////////////////////////////////////////////
65
66 using std::cout;
67 using std::endl;
68 ClassImp(AliTPCLaserTracks)
69
70   AliTPCLaserTracks::AliTPCLaserTracks():TObject(),
71     fId(-1),
72     fSide(-1),
73     fRod(-1),
74     fBundle(-1),
75     fBeam(-1),
76     fX(0.),
77     fY(0.),
78     fZ(0.),
79     fTime(0.),
80     fPhi(0.),
81     fTheta(0.),
82     fMaxSize(0),
83     fNpoints(0),
84     fXarr(0x0),
85     fYarr(0x0),
86     fZarr(0x0)
87
88 {
89     //
90     //  AliTPCLaserTracks default constructor
91     //
92 }
93
94 //_______________________________________________________________________
95 AliTPCLaserTracks::AliTPCLaserTracks(Int_t npoints):TObject(),
96     fId(-1),
97     fSide(-1),
98     fRod(-1),
99     fBundle(-1),
100     fBeam(-1),
101     fX(0.),
102     fY(0.),
103     fZ(0.),
104     fTime(0.),
105     fPhi(0.),
106     fTheta(0.),
107     fMaxSize(0),
108     fNpoints(0),
109     fXarr(0x0),
110     fYarr(0x0),
111     fZarr(0x0)
112
113 {
114     //
115     //  AliTPCLaserTracks constructor to initialise array of points
116     //
117     fNpoints = npoints;
118     InitPoints();
119 }
120
121 //_______________________________________________________________________
122 Double_t AliTPCLaserTracks::FindBeamLength(TVector3 vF, TVector3 vP)
123 {
124     //
125     //  Find distance between mirror and the intersection between
126     //  inner and outer field cage respectively
127     //
128     //  use the pq formula to find the intersection point between
129     //  the beam and the inner and outer fieldcage cylinders
130     //  The formulae solved are
131     //
132     //                        line = zylinder
133     //  xPoint = vF.X() + n*vP.X() = r*cos(phi)
134     //  yPoint = vF.Y() + n*vP.Y() = r*sin(phi)
135     //  zPoint = vF.Y() + n*vP.Y()
136     //
137     //  get n from the first two equations, calculate x,y,zPoint
138     //
139     //  vF is the mirror position and vP the beam direction
140     Double_t r[2];  //smalles r[0] and largest r[1] pad radius
141     r[0] = 80.0;  //smalles pad radius
142     r[1] = 260.0;   //largest pad radius
143
144     Double_t fxpxfypy = vF.X()*vP.X()+vF.Y()*vP.Y();
145     Double_t px2py2   = vP.X()*vP.X()+vP.Y()*vP.Y();
146
147     for (Int_t i=0; i<2; i++){
148         Double_t rad = (fxpxfypy/px2py2)*(fxpxfypy/px2py2)-
149                       (vF.X()*vF.X()+vF.Y()*vF.Y()-r[i]*r[i])/px2py2;
150
151         if ( rad >= 0 ){
152             Double_t n1   = -(fxpxfypy)/px2py2+TMath::Sqrt(rad);
153             TVector3 vI1(vF.X()+n1*vP.X(),
154                          vF.Y()+n1*vP.Y(),
155                          vF.Z()+n1*vP.Z());
156
157             Double_t n2   = -(fxpxfypy)/px2py2-TMath::Sqrt(rad);
158             TVector3 vI2(vF.X()+n2*vP.X(),
159                          vF.Y()+n2*vP.Y(),
160                          vF.Z()+n2*vP.Z());
161
162
163             //if we cross two boarders on our way return the closer boarder
164             if ( (n1>0) && (n2>0) ) {
165                 if ( (vF-vI1).Mag() <= (vF-vI2).Mag() )
166                     return (vF-vI1).Mag();
167                 else
168                     return (vF-vI2).Mag();
169             }
170
171             if ( n1>0 )
172                 return (vF-vI1).Mag();
173
174             if ( n2>0 )
175                 return (vF-vI2).Mag();
176
177         }
178     }
179
180     return 3.4e38;  //not found!!
181 }
182
183 //_______________________________________________________________________
184 void AliTPCLaserTracks::WriteTreeDesignData()
185 {
186     //
187     //   Write a tree with the design data of the laser track positions
188     //
189
190     TFile *f = TFile::Open("LaserTracksDesing.root","recreate");
191
192     AliTPCLaserTracks *ltp = new AliTPCLaserTracks(2);
193
194     TTree *t = new TTree("LaserTracks","LaserTracks");
195     t->Branch("LaserTracks","AliTPCLaserTracks",&ltp);
196     t->AutoSave();
197
198     Double_t phiBeam[7]         = {-31.8,-16.,-9.2,2.5,9.2,16.,31.8};
199     Double_t phiRod[6]          = {40.,100.,160.,220.,280.,340.};     //-20. for the muon side
200 //    Double_t zBundleCorse[2][4] = {{10,79,163,241},{13,85,169,247}};
201     Double_t zBundleCorse[2][4] = {{241.,163.,79.,10.},{247.,169.,85.,13.}};
202     Double_t zBeamFine[7]       = {1.45,1.3,1.15,1.3,1.15,1.3,1.45};
203
204     Int_t id=0;
205     //loop over TPC side -- 0:A (Shaft) side -- 1:C (Muon) side
206     for (Int_t side=0; side<2; side++){
207         //loop over laser rods -- counterclockwise
208         for (Int_t rod=0; rod<6; rod++){
209             //center of the rod
210             TVector2 vRod;
211
212             Double_t phiRod2 = phiRod[rod]-side*20;
213             vRod.SetMagPhi(254.25,                   //center of the rod at 254.25cm
214                            TMath::DegToRad()*
215                            (phiRod2)  //-20 deg on C-Side
216                           );
217
218             //loop over bundle -- counted from large z to small z
219             for (Int_t bundle=0; bundle<4; bundle++){
220                 //center of the bundle; not yet rotated
221                 Int_t bundleS = bundle;
222                 //              if ( side == 1 ) bundleS = 4-bundle;
223                 if ( side == 1 ) bundleS = 3-bundle;
224
225                 TVector2 vBundle;
226                 vBundle.SetMagPhi(.65,
227                                   TMath::DegToRad()*
228                                 //?  TMath::Power(-1,side)*
229                                   (phiRod2+180.-90.*bundleS)
230                                  );
231
232                 //loop over beam
233                 Int_t i=0;
234                 for (Int_t beam=0; beam<7; beam++){
235                     TVector2 vBeam;
236                     if ( beam == 3 )
237                         vBeam.Set(0.,0.);
238                     else{
239                         vBeam.SetMagPhi(1.,TMath::DegToRad()*(phiRod2+i*60));
240                         i++;
241                     }
242
243                     TVector2 xyMirror = vRod+vBundle+vBeam;
244                     Double_t zMirror   = (zBundleCorse[rod%2][bundleS]+zBeamFine[beam])*TMath::Power(-1,side);
245                     TVector3 v3Mirror(xyMirror.X(),xyMirror.Y(),zMirror);
246
247                     Double_t phi      = TMath::DegToRad()*(phiRod2+180+phiBeam[beam]*TMath::Power(-1,side));
248                     Double_t theta    = TMath::DegToRad()*90;
249
250                     TVector3 v3Beam;
251                     v3Beam.SetMagThetaPhi(1,theta,phi);  //Direction Vector
252                     v3Beam.SetMagThetaPhi(FindBeamLength(v3Mirror,v3Beam), theta, phi);
253
254
255
256                     ltp->SetId(id);
257                     ltp->SetSide(side);
258                     ltp->SetRod(rod);
259                     ltp->SetBundle(bundleS);
260                     ltp->SetBeam(beam);
261
262                     ltp->SetX(xyMirror.X());
263                     ltp->SetY(xyMirror.Y());
264                     ltp->SetZ(zMirror);
265
266                     ltp->SetPhi(phi);
267                     ltp->SetTheta(theta);
268
269 //                  ltp->InitPoints(2);
270                     ltp->SetPoint(0,xyMirror.X(),xyMirror.Y(),zMirror);
271                     ltp->SetPoint(1,xyMirror.X()+v3Beam.X(),xyMirror.Y()+v3Beam.Y(),zMirror+v3Beam.Z());
272
273                     cout<< "Id: " << id
274                         << "  Side: " << side
275                         << "  Rod: " << rod
276                         << "  Bundle: " << bundleS
277                         << "  Beam: " << beam
278                         << endl;
279
280
281                     t->Fill();
282                     //delete line;
283                     id++;
284 //                  cout << "Filled!!!" << endl;
285                 }
286             }
287         }
288     }
289
290     t->Write();
291     delete f;
292 //    delete t;
293 //    delete ltp;
294 //    delete line;
295 }
296
297 //_______________________________________________________________________
298 TPolyLine3D* AliTPCLaserTracks::GetLine()
299 {
300   //
301   // Make a polilyne object oout of the points
302   //
303    if ( fNpoints ) return new TPolyLine3D(fNpoints,fXarr,fYarr,fZarr);
304
305     return 0x0;
306 }
307
308 //_______________________________________________________________________
309 TObjArray* AliTPCLaserTracks::GetLines(const Char_t* file, const Char_t *cuts)
310 {
311   //
312   // Read the input files with the laser track 
313   // Make a array of polylines for visualization
314   
315     TObjArray *array = new TObjArray;
316
317     TFile *f = TFile::Open(file,"read");
318     TTree *tree = (TTree*)f->Get("LaserTracks");
319
320     AliTPCLaserTracks *ltp = new AliTPCLaserTracks();
321 //    TEventList *evList = new TEventList("evList");
322     TEventList evList("evList");
323
324     tree->SetBranchAddress("LaserTracks",&ltp);
325
326     tree->Draw(">>evList",cuts,"goff");
327
328 //    cout << "N: " << evList.GetN() << endl;
329     gROOT->cd();
330     for (Int_t ev = 0; ev < evList.GetN(); ev++){
331 //        cout << ev << endl;
332         tree->GetEntry( evList.GetEntry(ev) );
333         array->Add(ltp->GetLine());
334     }
335 //    cout << "delte f"<< endl;
336 //    delete f;
337
338 //    cout << "evlist" << endl;
339
340     return array;
341 }
342
343 //_______________________________________________________________________
344 void AliTPCLaserTracks::InitPoints()
345 {
346     //
347     //  Init point arrays
348     //
349
350     fMaxSize = fNpoints;
351     fXarr = new Double_t[fMaxSize];
352     fYarr = new Double_t[fMaxSize];
353     fZarr = new Double_t[fMaxSize];
354 }
355
356 //_______________________________________________________________________
357 Int_t AliTPCLaserTracks::SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
358 {
359     //
360     // Set point to point arrays
361     //
362
363     if ( !fXarr || !fYarr || !fZarr ) return 1;
364     if ( point > fNpoints-1 ) return 1;
365
366     fXarr[point] = x;
367     fYarr[point] = y;
368     fZarr[point] = z;
369     return 0;
370 }
371
372 //_______________________________________________________________________
373 Int_t AliTPCLaserTracks::FindMirror(Char_t *file, Double_t x, Double_t y, Double_t z, Double_t phi)
374 {
375     //
376     //  return the id of the mirror in 'file' that maches best the given parameters
377     //
378
379     TFile *f = TFile::Open(file,"read");
380     TTree *tree = (TTree*)f->Get("LaserTracks");
381
382     AliTPCLaserTracks *ltp = new AliTPCLaserTracks();
383     TEventList evList("evList");
384
385     tree->SetBranchAddress("LaserTracks",&ltp);
386
387     TString s;
388     s = "abs(fX-"; s+= x;
389     s+= ")<3&&abs(fY-"; s+=y;
390     s+= ")<3&&abs(fZ-"; s+= z;
391     s+= ")<1.5";
392 //    s+= "&&((abs(fPhi-";s+= phi;
393 //    s+= ")<.06)||(abs(fPhi-";s+= (phi-TMath::DegToRad()*180);
394 //    s+= ")<.06))";
395     cout << s.Data() << endl;
396
397     tree->Draw(">>evList",s.Data());
398
399     Double_t dphiMin=TMath::Pi();
400     Int_t id=-1;
401
402     cout << "nev: " << evList.GetN() << endl;
403     for (Int_t i=0; i<evList.GetN(); i++){
404         tree->GetEntry( evList.GetEntry(i) );
405         Double_t dphi = TMath::Abs(ltp->GetPhi() - phi);
406         if ( dphi > TMath::DegToRad()*180 ) dphi-=TMath::DegToRad()*180;
407         if ( dphi<dphiMin ) {
408             dphiMin = dphi;
409             id = ltp->GetId();
410         }
411     }
412     delete f;
413     return id;
414 }
415
416 //_______________________________________________________________________
417 AliTPCLaserTracks::~AliTPCLaserTracks()
418 {
419     //
420     //  standard destructor
421     //
422
423 //    if ( fP ) delete fP;
424     if ( fXarr ) delete [] fXarr;
425     if ( fYarr ) delete [] fYarr;
426     if ( fZarr ) delete [] fZarr;
427 }