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