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