Logo ROOT  
Reference Guide
THelix.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Ping Yeh 19/12/97
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class THelix
13\ingroup g3d
14THelix has two different constructors.
15
16If a particle with charge q passes through a point (x,y,z)
17with momentum (px,py,pz) with magnetic field B along an axis (nx,ny,nz),
18this helix can be constructed like:
19
20~~~ {.cpp}
21 THelix p(x,y,z, px,py,pz, q*B, nx,ny,nz);
22 (nx,ny,nz) defaults to (0,0,1).
23~~~
24
25A helix in its own frame can be defined with a pivotal point
26(x0,y0,z0), the velocity at that point (vx0,vy0,vz0), and
27an angular frequency w. Combining vx0 and vy0 to a transverse
28velocity vt0 one can parametrize the helix as:
29
30~~~ {.cpp}
31 x(t) = x0 - vt0 / w * sin(-w * t + phi0)
32 y(t) = y0 + vt0 / w * cos(-w * t + phi0)
33 z(t) = z0 + vz0 * t
34~~~
35
36The second constructor has 6 parameters,
37
38Example:
39
40~~~ {.cpp}
41 THelix pl1(xyz, v, w, range, rtype, axis);
42~~~
43
44where:
45
46 - xyz : array of initial position
47 - v : array of initial velocity
48 - w : angular frequency
49 - range: helix range
50 - rtype: kHelixZ specifies allowed drawing range in helix Z direction, i.e., along B field.
51 kLabZ specifies drawing range in lab frame.
52 kHelixX, kHelixY, kLabX, kLabY, kUnchanged ... etc can also be specified
53 - axis : helix axis
54
55Example constructing a helix with several default values and drawing it:
56
57Begin_Macro(source)
58{
59 TCanvas* helix_example_c1 = new TCanvas("helix_example_c1");
60 TView *view = TView::CreateView(1);
61 view->SetRange(-1,-1,-1,1,1,1);
62 THelix *helix = new THelix(0., 0., 0., 1., 0., 0.3, 10.);
63 helix->Draw();
64}
65End_Macro
66
67This initializes a helix with its axis in Z direction (rtype=kHelixZ).
68*/
69
70#include "Riostream.h"
71#include "TBuffer.h"
72#include "TROOT.h"
73#include "TVirtualPad.h"
74#include "THelix.h"
75#include "TClass.h"
76#include "TMath.h"
77
78Int_t THelix::fgMinNSeg=5; // at least 5 line segments in TPolyLine3D
79
81
82////////////////////////////////////////////////////////////////////////////////
83/// Set all helix parameters.
84
85void THelix::SetHelix(Double_t const* xyz, Double_t const* v, Double_t w,
86 Double_t const* range, EHelixRangeType rType,
87 Double_t const* axis )
88{
89 // Define the helix frame by setting the helix axis and rotation matrix
90 SetAxis(axis);
91
92 // Calculate initial position and velocity in helix frame
93 fW = w;
95 Double_t vx0, vy0, vz0;
96 vx0 = v[0] * m[0] + v[1] * m[1] + v[2] * m[2];
97 vy0 = v[0] * m[3] + v[1] * m[4] + v[2] * m[5];
98 vz0 = v[0] * m[6] + v[1] * m[7] + v[2] * m[8];
99 fVt = TMath::Sqrt(vx0*vx0 + vy0*vy0);
100 fPhi0 = TMath::ATan2(vy0,vx0);
101 fVz = vz0;
102 fX0 = xyz[0] * m[0] + xyz[1] * m[1] + xyz[2] * m[2];
103 fY0 = xyz[0] * m[3] + xyz[1] * m[4] + xyz[2] * m[5];
104 fZ0 = xyz[0] * m[6] + xyz[1] * m[7] + xyz[2] * m[8];
105 if (fW != 0) {
106 fX0 += fVt / fW * TMath::Sin(fPhi0);
107 fY0 -= fVt / fW * TMath::Cos(fPhi0);
108 }
109
110 // Then calculate the range in t and set the polyline representation
111 Double_t r1 = 0;
112 Double_t r2 = 1;
113 if (range) {r1 = range[0]; r2 = range[1];}
114 if (rType != kUnchanged) {
115 fRange[0] = 0.0; fRange[1] = TMath::Pi(); // initialize to half round
116 SetRange(r1,r2,rType);
117 }
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// Helix default constructor.
122
124{
125 fX0 = fY0 = fZ0 = fVt = fPhi0 = fVz = fAxis[0] = fAxis[1] = 0.0;
126 fAxis[2] = 1.0;
127 fW = 1.5E7; // roughly the cyclon frequency of proton in AMS
128 fRange[0] = 0.0;
129 fRange[1] = 1.0;
130 fRotMat = 0;
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// Helix normal constructor.
135
137 Double_t vx, Double_t vy, Double_t vz,
138 Double_t w)
139 : TPolyLine3D()
140{
141 Double_t p[3], v[3];
142 p[0] = x;
143 p[1] = y;
144 p[2] = z;
145 v[0] = vx;
146 v[1] = vy;
147 v[2] = vz;
148 Double_t *range = 0;
149 fRotMat = 0;
150
151 SetHelix(p, v, w, range, kHelixZ);
152 fOption = "";
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// Helix normal constructor.
157
159 Double_t const* range, EHelixRangeType rType, Double_t const* axis)
160 : TPolyLine3D()
161{
162 Double_t r[2];
163 if ( range ) {
164 r[0] = range[0]; r[1] = range[1];
165 } else {
166 r[0] = 0.0; r[1] = 1.0;
167 }
168
169 fRotMat = 0;
170 if ( axis ) { // specify axis
171 SetHelix(xyz, v, w, r, rType, axis);
172 } else { // default axis
173 SetHelix(xyz, v, w, r, rType);
174 }
175
176 fOption = "";
177}
178
179
180#if 0
181////////////////////////////////////////////////////////////////////////////////
182/// Helix copy constructor.
183
185{
186 fX0 = h.fX0;
187 fY0 = h.fY0;
188 fZ0 = h.fZ0;
189 fVt = h.fVt;
190 fPhi0 = h.fPhi0;
191 fVz = h.fVz;
192 fW = h.fW;
193 for (Int_t i=0; i<3; i++) fAxis[i] = h.fAxis[i];
194 fRotMat = new TRotMatrix(*(h.fRotMat));
195 fRange[0] = h.fRange[0];
196 fRange[1] = h.fRange[1];
197
198 fOption = h.fOption;
199}
200#endif
201
202////////////////////////////////////////////////////////////////////////////////
203/// assignment operator
204
206{
207 if(this!=&hx) {
209 fX0=hx.fX0;
210 fY0=hx.fY0;
211 fZ0=hx.fZ0;
212 fVt=hx.fVt;
213 fPhi0=hx.fPhi0;
214 fVz=hx.fVz;
215 fW=hx.fW;
216 for(Int_t i=0; i<3; i++)
217 fAxis[i]=hx.fAxis[i];
218 fRotMat=hx.fRotMat;
219 for(Int_t i=0; i<2; i++)
220 fRange[i]=hx.fRange[i];
221 }
222 return *this;
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// Helix destructor.
227
229{
230 if (fRotMat) delete fRotMat;
231}
232
233
234////////////////////////////////////////////////////////////////////////////////
235/// Helix copy constructor.
236
237THelix::THelix(const THelix &helix) : TPolyLine3D(helix)
238{
239 fRotMat=0;
240 ((THelix&)helix).THelix::Copy(*this);
241}
242
243
244////////////////////////////////////////////////////////////////////////////////
245/// Copy this helix to obj.
246
247void THelix::Copy(TObject &obj) const
248{
249 TObject::Copy(obj);
250 TAttLine::Copy(((THelix&)obj));
251
252 ((THelix&)obj).fX0 = fX0;
253 ((THelix&)obj).fY0 = fY0;
254 ((THelix&)obj).fZ0 = fZ0;
255 ((THelix&)obj).fVt = fVt;
256 ((THelix&)obj).fPhi0 = fPhi0;
257 ((THelix&)obj).fVz = fVz;
258 ((THelix&)obj).fW = fW;
259 for (Int_t i=0; i<3; i++)
260 ((THelix&)obj).fAxis[i] = fAxis[i];
261
262 if (((THelix&)obj).fRotMat)
263 delete ((THelix&)obj).fRotMat;
264 ((THelix&)obj).fRotMat = new TRotMatrix(*fRotMat);
265
266 ((THelix&)obj).fRange[0] = fRange[0];
267 ((THelix&)obj).fRange[1] = fRange[1];
268
269 ((THelix&)obj).fOption = fOption;
270
271 //
272 // Set range and make the graphic representation
273 //
274 ((THelix&)obj).SetRange(fRange[0], fRange[1], kHelixT);
275}
276
277
278////////////////////////////////////////////////////////////////////////////////
279/// Draw this helix with its current attributes.
280
282{
283 AppendPad(option);
284}
285
286
287////////////////////////////////////////////////////////////////////////////////
288/// Dump this helix with its attributes.
289
290void THelix::Print(Option_t *option) const
291{
292 std::cout <<" THelix Printing N=" <<fN<<" Option="<<option<<std::endl;
293}
294
295
296////////////////////////////////////////////////////////////////////////////////
297/// Save primitive as a C++ statement(s) on output stream out.
298
299void THelix::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
300{
301 char quote = '"';
302 out<<" "<<std::endl;
303 if (gROOT->ClassSaved(THelix::Class())) {
304 out<<" ";
305 } else {
306 out<<" THelix *";
307 }
308 out<<"helix = new THelix("<<fX0<<","<<fY0<<","<<fZ0<<","
309 <<fVt*TMath::Cos(fPhi0)<<","<<fVt*TMath::Sin(fPhi0)<<","<<fVz<<","
310 <<fW<<","<<fRange[0]<<","<<fRange[1]<<","<<(Int_t)kHelixT<<","
311 <<fAxis[0]<<","<<fAxis[1]<<","<<fAxis[2]<<","
312 <<quote<<fOption<<quote<<");"<<std::endl;
313
314 SaveLineAttributes(out,"helix",1,1,1);
315
316 out<<" helix->Draw();"<<std::endl;
317}
318
319
320////////////////////////////////////////////////////////////////////////////////
321/// Set a new axis for the helix. This will make a new rotation matrix.
322
323void THelix::SetAxis(Double_t const* axis)
324{
325 if (axis) {
326 Double_t len = TMath::Sqrt(axis[0]*axis[0] + axis[1]*axis[1] + axis[2]*axis[2]);
327 if (len <= 0) {
328 Error("SetAxis()", "Impossible! axis length %lf <= 0!", len);
329 return;
330 }
331 fAxis[0] = axis[0]/len;
332 fAxis[1] = axis[1]/len;
333 fAxis[2] = axis[2]/len;
334 } else {
335 fAxis[0] = 0;
336 fAxis[1] = 0;
337 fAxis[2] = 1;
338 }
339
340 // Construct the rotational matrix from the axis
341 SetRotMatrix();
342}
343
344
345////////////////////////////////////////////////////////////////////////////////
346/// Set axis.
347
349{
350 Double_t axis[3]; axis[0] = x; axis[1] = y; axis[2] = z;
351 SetAxis(axis);
352}
353
354
355////////////////////////////////////////////////////////////////////////////////
356/// Set a new range for the helix. This will remake the polyline.
357
359{
360 Double_t a[2];
361 Double_t halfpi = TMath::Pi()/2.0;
362 Int_t i;
365 Double_t phase;
366
367 if ( fW != 0 && fVz != 0 ) { // general case
368 switch ( rType ) {
369 case kHelixT :
370 fRange[0] = range[0]; fRange[1] = range[1]; break;
371
372 case kHelixX :
373 for (i=0; i<2; i++ ) {
374 a[i] = fW / fVt * (range[i] - fX0);
375 if ( a[i] < -1 || a[i] > 1 ) {
376 Error("SetRange()",
377 "range out of bound (%lf:%lf): %lf. Default used: %lf",
378 fX0-fVt/fW, fX0+fVt/fW, range[i], fRange[i]);
379 return;
380 }
381 phase = FindClosestPhase(fPhi0+halfpi, a[i]);
382 fRange[i] = ( fPhi0 + halfpi - phase ) / fW;
383 }
384 break;
385
386 case kHelixY :
387 for (i=0; i<2; i++ ) {
388 a[i] = fW / fVt * (range[i] - fY0);
389 if ( a[i] < -1 || a[i] > 1 ) {
390 Error("SetRange()",
391 "range out of bound (%lf:%lf): %lf. Default used: %lf",
392 fY0-fVt/fW, fY0+fVt/fW, range[i], fRange[i]);
393 return;
394 }
395 phase = FindClosestPhase(fPhi0, a[i]);
396 fRange[i] = ( fPhi0 - phase ) / fW;
397 }
398 break;
399
400 case kHelixZ :
401 if ( fVz != 0 ) {
402 for (i=0; i<2; i++ ) {
403 fRange[i] = (range[i] - fZ0) / fVz;
404 }
405 } else { // fVz = 0, z = constant = fZ0
406 Error("SetRange()",
407 "Vz = 0 and attempts to set range along helix axis!");
408 return;
409 }
410 break;
411
412 case kLabX :
413 case kLabY :
414 case kLabZ :
415 printf("setting range in lab axes is not implemented yet\n");
416 break;
417 default:
418 Error("SetRange()","unknown range type %d", rType);
419 break;
420 }
421 } else if ( fW == 0 ) { // straight line: x = x0 + vx * t
422 switch ( rType ) {
423 case kHelixT :
424 fRange[0] = range[0]; fRange[1] = range[1];
425 break;
426 case kHelixX :
427 if ( vx != 0 ) {
428 fRange[0] = (range[0] - fX0) / vx;
429 fRange[1] = (range[1] - fX0) / vx;
430 } else {
431 Error("SetRange()",
432 "Vx = 0 and attempts to set range on helix x axis!");
433 return;
434 }
435 break;
436 case kHelixY :
437 if ( vy != 0 ) {
438 fRange[0] = (range[0] - fY0) / vy;
439 fRange[1] = (range[1] - fY0) / vy;
440 } else {
441 Error("SetRange()",
442 "Vy = 0 and attempts to set range on helix y axis!");
443 return;
444 }
445 break;
446 case kHelixZ :
447 if ( fVz != 0 ) {
448 fRange[0] = (range[0] - fZ0) / fVz;
449 fRange[1] = (range[1] - fZ0) / fVz;
450 } else {
451 Error("SetRange()",
452 "Vz = 0 and attempts to set range on helix z axis!");
453 return;
454 }
455 break;
456 case kLabX :
457 case kLabY :
458 case kLabZ :
459 printf("setting range in lab axes is not implemented yet\n");
460 break;
461 default :
462 Error("SetRange()","unknown range type %d", rType);
463 break;
464 }
465 } else if ( fVz == 0 ) { // a circle, not fully implemented yet
466 switch ( rType ) {
467 case kHelixT :
468 fRange[0] = range[0]; fRange[1] = range[1]; break;
469 case kHelixX :
470 if ( vx != 0 ) {
471 fRange[0] = (range[0] - fX0) / vx;
472 fRange[1] = (range[1] - fX0) / vx;
473 } else {
474 Error("SetRange()",
475 "Vx = 0 and attempts to set range on helix x axis!");
476 return;
477 }
478 break;
479 case kHelixY :
480 if ( vy != 0 ) {
481 fRange[0] = (range[0] - fY0) / vy;
482 fRange[1] = (range[1] - fY0) / vy;
483 } else {
484 Error("SetRange()",
485 "Vy = 0 and attempts to set range on helix y axis!");
486 return;
487 }
488 break;
489 case kHelixZ :
490 Error("SetRange()",
491 "Vz = 0 and attempts to set range on helix z axis!");
492 return;
493 case kLabX :
494 case kLabY :
495 case kLabZ :
496 printf("setting range in lab axes is not implemented yet\n");
497 break;
498 default :
499 Error("SetRange()","unknown range type %d", rType);
500 break;
501 }
502 }
503
504 if ( fRange[0] > fRange[1] ) {
505 Double_t temp = fRange[1]; fRange[1] = fRange[0]; fRange[0] = temp;
506 }
507
508 // Set the polylines in global coordinates
509 Double_t degrad = TMath::Pi() / 180.0;
510 Double_t segment = 5.0 * degrad; // 5 degree segments
511 Double_t dt = segment / TMath::Abs(fW); // parameter span on each segm.
512
513 Int_t nSeg = Int_t((fRange[1]-fRange[0]) / dt) + 1;
514 if (nSeg < THelix::fgMinNSeg) {
515 nSeg = THelix::fgMinNSeg;
516 dt = (fRange[1]-fRange[0])/nSeg;
517 }
518
519 Double_t * xl = new Double_t[nSeg+1]; // polyline in local coordinates
520 Double_t * yl = new Double_t[nSeg+1];
521 Double_t * zl = new Double_t[nSeg+1];
522
523 for (i=0; i<=nSeg; i++) { // calculate xl[], yl[], zl[];
524 Double_t t, phase2;
525 if (i==nSeg) t = fRange[1]; // the last point
526 else t = fRange[0] + dt * i;
527 phase2 = -fW * t + fPhi0;
528 xl[i] = fX0 - fVt/fW * TMath::Sin(phase2);
529 yl[i] = fY0 + fVt/fW * TMath::Cos(phase2);
530 zl[i] = fZ0 + fVz * t;
531 }
532
533 Float_t xg, yg,zg; // global coordinates
534 // must be Float_t to call TPolyLine3D::SetPoint()
537 for (i=0; i<=nSeg; i++) { // m^{-1} = transpose of m
538 xg = xl[i] * m[0] + yl[i] * m[3] + zl[i] * m[6] ;
539 yg = xl[i] * m[1] + yl[i] * m[4] + zl[i] * m[7] ;
540 zg = xl[i] * m[2] + yl[i] * m[5] + zl[i] * m[8] ;
541 TPolyLine3D::SetPoint(i,xg,yg,zg);
542 }
543
544 delete[] xl; delete[] yl; delete[] zl;
545}
546
547
548////////////////////////////////////////////////////////////////////////////////
549/// Set range.
550
552{
553 Double_t range[2];
554 range[0] = r1; range[1] = r2;
555 SetRange(range, rType);
556}
557
558
559////////////////////////////////////////////////////////////////////////////////
560// //
561// Protected Member Functions //
562// //
563////////////////////////////////////////////////////////////////////////////////
564
565
566////////////////////////////////////////////////////////////////////////////////
567/// Set the rotational matrix according to the helix axis.
568
570{
571 // Calculate all 6 angles.
572 // Note that TRotMatrix::TRotMatrix() expects angles in degrees.
573 Double_t raddeg = 180.0 / TMath::Pi();
574 Double_t halfpi = TMath::Pi()/2.0 * raddeg;
575 // (theta3,phi3) is the helix axis
576 Double_t theta3 = TMath::ACos(fAxis[2]) * raddeg;
577 Double_t phi3 = TMath::ATan2(fAxis[1], fAxis[0]) * raddeg;
578 // (theta1,phi1) is the x-axis in helix frame
579 Double_t theta1 = theta3 + halfpi;
580 Double_t phi1 = phi3;
581 // (theta2,phi2) is the y-axis in helix frame
582 Double_t theta2 = halfpi;
583 Double_t phi2 = phi1 + halfpi;
584
585 // Delete the old rotation matrix
586 if (fRotMat) delete fRotMat;
587
588 // Make a new rotation matrix
589 fRotMat = new TRotMatrix("HelixRotMat", "Master frame -> Helix frame",
590 theta1, phi1, theta2, phi2, theta3, phi3 );
591 return;
592}
593
594
595////////////////////////////////////////////////////////////////////////////////
596/// Finds the closest phase to phi0 that gives cos(phase) = cosine
597
599{
600 const Double_t pi = TMath::Pi();
601 const Double_t twopi = TMath::Pi() * 2.0;
602 Double_t phi1 = TMath::ACos(cosine);
603 Double_t phi2 = - phi1;
604
605 while ( phi1 - phi0 > pi ) phi1 -= twopi;
606 while ( phi1 - phi0 < -pi ) phi1 += twopi;
607
608 while ( phi2 - phi0 > pi ) phi2 -= twopi;
609 while ( phi2 - phi0 < -pi ) phi2 += twopi;
610
611 // Now phi1, phi2 and phi0 are within the same 2pi range
612 // and cos(phi1) = cos(phi2) = cosine
613 if ( TMath::Abs(phi1-phi0) < TMath::Abs(phi2-phi0) ) return phi1;
614 else return phi2;
615}
616
617
618////////////////////////////////////////////////////////////////////////////////
619/// Stream an object of class THelix.
620
621void THelix::Streamer(TBuffer &R__b)
622{
623 if (R__b.IsReading()) {
624 UInt_t R__s, R__c;
625 Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
626 if (R__v > 1) {
627 R__b.ReadClassBuffer(THelix::Class(), this, R__v, R__s, R__c);
628 return;
629 }
630 //====process old versions before automatic schema evolution
631 TPolyLine3D::Streamer(R__b);
632 R__b >> fX0;
633 R__b >> fY0;
634 R__b >> fZ0;
635 R__b >> fVt;
636 R__b >> fPhi0;
637 R__b >> fVz;
638 R__b >> fW;
640 R__b >> fRotMat;
642 R__b.CheckByteCount(R__s, R__c, THelix::IsA());
643 //====end of old versions
644
645 } else {
646 R__b.WriteClassBuffer(THelix::Class(),this);
647 }
648}
void Class()
Definition: Class.C:29
ROOT::R::TRInterface & r
Definition: Object.C:4
#define h(i)
Definition: RSha256.hxx:106
int Int_t
Definition: RtypesCore.h:43
short Version_t
Definition: RtypesCore.h:63
unsigned int UInt_t
Definition: RtypesCore.h:44
double Double_t
Definition: RtypesCore.h:57
float Float_t
Definition: RtypesCore.h:55
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
EHelixRangeType
Definition: THelix.h:18
@ kHelixY
Definition: THelix.h:19
@ kHelixX
Definition: THelix.h:19
@ kLabZ
Definition: THelix.h:19
@ kHelixT
Definition: THelix.h:19
@ kLabX
Definition: THelix.h:19
@ kHelixZ
Definition: THelix.h:19
@ kLabY
Definition: THelix.h:19
@ kUnchanged
Definition: THelix.h:19
#define gROOT
Definition: TROOT.h:406
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition: TAttLine.cxx:172
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:270
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition: TBuffer.h:85
virtual Int_t ReadStaticArray(Bool_t *b)=0
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
THelix has two different constructors.
Definition: THelix.h:23
Double_t fRange[2]
Definition: THelix.h:35
Double_t fZ0
Definition: THelix.h:28
virtual void Draw(Option_t *option="")
Draw this helix with its current attributes.
Definition: THelix.cxx:281
THelix & operator=(const THelix &)
assignment operator
Definition: THelix.cxx:205
virtual void SetAxis(Double_t const *axis)
Set a new axis for the helix. This will make a new rotation matrix.
Definition: THelix.cxx:323
virtual void Copy(TObject &helix) const
Copy this helix to obj.
Definition: THelix.cxx:247
Double_t fPhi0
Definition: THelix.h:30
void SetHelix(Double_t const *xyz, Double_t const *v, Double_t w, Double_t const *range=0, EHelixRangeType type=kUnchanged, Double_t const *axis=0)
Set all helix parameters.
Definition: THelix.cxx:85
virtual void SetRange(Double_t *range, EHelixRangeType rtype=kHelixZ)
Set a new range for the helix. This will remake the polyline.
Definition: THelix.cxx:358
THelix()
Helix default constructor.
Definition: THelix.cxx:123
Double_t fVz
Definition: THelix.h:31
Double_t fAxis[3]
Definition: THelix.h:33
Double_t fX0
Definition: THelix.h:26
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: THelix.cxx:299
virtual ~THelix()
Helix destructor.
Definition: THelix.cxx:228
Double_t FindClosestPhase(Double_t phi0, Double_t cosine)
Finds the closest phase to phi0 that gives cos(phase) = cosine.
Definition: THelix.cxx:598
void SetRotMatrix()
Set the rotational matrix according to the helix axis.
Definition: THelix.cxx:569
Double_t fY0
Definition: THelix.h:27
Double_t fVt
Definition: THelix.h:29
virtual void Print(Option_t *option="") const
Dump this helix with its attributes.
Definition: THelix.cxx:290
TRotMatrix * fRotMat
Definition: THelix.h:34
Double_t fW
Definition: THelix.h:32
static Int_t fgMinNSeg
Definition: THelix.h:42
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
virtual void Copy(TObject &object) const
Copy this to obj.
Definition: TObject.cxx:61
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:891
A 3-dimensional polyline.
Definition: TPolyLine3D.h:32
TPolyLine3D & operator=(const TPolyLine3D &polylin)
assignment operator
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
Int_t fN
Number of points.
Definition: TPolyLine3D.h:34
TString fOption
options
Definition: TPolyLine3D.h:36
virtual void SetPolyLine(Int_t n, Option_t *option="")
Re-initialize polyline with n points (0,0,0).
Manages a detector rotation matrix.
Definition: TRotMatrix.h:28
virtual Double_t * GetMatrix()
Definition: TRotMatrix.h:54
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
static constexpr double halfpi
static constexpr double pi
static constexpr double twopi
Double_t ACos(Double_t)
Definition: TMath.h:658
Double_t ATan2(Double_t y, Double_t x)
Definition: TMath.h:669
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Double_t Cos(Double_t)
Definition: TMath.h:631
constexpr Double_t Pi()
Definition: TMath.h:38
Double_t Sin(Double_t)
Definition: TMath.h:627
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * m
Definition: textangle.C:8
auto * a
Definition: textangle.C:12