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