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