Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoTrack.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 2003/04/10
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#include "TBrowser.h"
13#include "TPoint.h"
14#include "TVirtualPad.h"
15#include "TVirtualX.h"
16#include "TView.h"
17#include "TGeoManager.h"
18#include "TVirtualGeoPainter.h"
19#include "TGeoTrack.h"
20
21/** \class TGeoTrack
22\ingroup Geometry_painter
23
24\deprecated
25Use of TGeoTrack is deprecated. For the event display please switch to TEve (using TEveTracks
26or TEvePointSet to display tracks specifically) or to [REve](https://root.cern/doc/master/tracks_8C_source.html).
27
28Class for user-defined tracks attached to a geometry.
29Tracks are 3D objects made of points and they store a
30pointer to a TParticle. The geometry manager holds a list
31of all tracks that will be deleted on destruction of
32gGeoManager.
33*/
34
36
37////////////////////////////////////////////////////////////////////////////////
38/// Default constructor.
39
41{
42 fPointsSize = 0;
43 fNpoints = 0;
44 fPoints = 0;
45}
46
47////////////////////////////////////////////////////////////////////////////////
48/// Constructor.
49
51 :TVirtualGeoTrack(id,pdgcode,parent,particle)
52{
53 fPointsSize = 0;
54 fNpoints = 0;
55 fPoints = 0;
56 if (fParent==0) {
59 SetMarkerSize(0.6);
60 SetLineColor(2);
61 SetLineWidth(2);
62 } else {
65 SetMarkerSize(0.6);
66 SetLineColor(4);
67 SetLineWidth(2);
68 }
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// Destructor.
73
75{
76 if (fPoints) delete [] fPoints;
77// if (gPad) gPad->GetListOfPrimitives()->Remove(this);
78}
79
80////////////////////////////////////////////////////////////////////////////////
81/// Add a daughter track to this.
82
84{
85 if (!fTracks) fTracks = new TObjArray(1);
86 Int_t index = fTracks->GetEntriesFast();
87 TGeoTrack *daughter = new TGeoTrack(id,pdgcode,this,particle);
88 fTracks->AddAtAndExpand(daughter,index);
89 return daughter;
90}
91
92////////////////////////////////////////////////////////////////////////////////
93/// Add a daughter and return its index.
94
96{
97 if (!fTracks) fTracks = new TObjArray(1);
98 Int_t index = fTracks->GetEntriesFast();
99 fTracks->AddAtAndExpand(other,index);
100 other->SetParent(this);
101 return index;
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// Draw animation of this track
106
108{
109 if (tmin<0 || tmin>=tmax || nframes<1) return;
112 if (!gPad) {
114 }
115 TList *list = gPad->GetListOfPrimitives();
116 TIter next(list);
117 TObject *obj;
118 while ((obj = next())) {
119 if (!strcmp(obj->ClassName(), "TGeoTrack")) list->Remove(obj);
120 }
121 Double_t dt = (tmax-tmin)/Double_t(nframes);
122 Double_t delt = 2E-9;
123 Double_t t = tmin;
124 Bool_t geomanim = kFALSE;
125 Bool_t issave = kFALSE;
126 TString fname;
127
128 TString opt(option);
129 if (opt.Contains("/G")) geomanim = kTRUE;
130 if (opt.Contains("/S")) issave = kTRUE;
131
133 Double_t *box = p->GetViewBox();
134 box[0] = box[1] = box[2] = 0;
135 box[3] = box[4] = box[5] = 100;
137 Draw(opt.Data());
138 Double_t start[6] = {0,0,0,0,0,0};
139 Double_t end[6] = {0,0,0,0,0,0};
140 Int_t i, j;
141 Double_t dlat=0, dlong=0, dpsi=0;
142 Double_t dd[6] = {0,0,0,0,0,0};
143 if (geomanim) {
144 p->EstimateCameraMove(tmin+5*dt, tmin+15*dt, start, end);
145 for (i=0; i<3; i++) {
146 start[i+3] = 20 + 1.3*start[i+3];
147 end[i+3] = 20 + 0.9*end[i+3];
148 }
149 for (i=0; i<6; i++) {
150 dd[i] = (end[i]-start[i])/10.;
151 }
152 memcpy(box, start, 6*sizeof(Double_t));
153 p->GetViewAngles(dlong,dlat,dpsi);
154 dlong = (-206-dlong)/Double_t(nframes);
155 dlat = (126-dlat)/Double_t(nframes);
156 dpsi = (75-dpsi)/Double_t(nframes);
157 p->GrabFocus();
158 }
159
160 for (i=0; i<nframes; i++) {
161 if (t-delt<0) gGeoManager->SetTminTmax(0,t);
162 else gGeoManager->SetTminTmax(t-delt,t);
163 if (geomanim) {
164 for (j=0; j<6; j++) box[j]+=dd[j];
165 p->GrabFocus(1,dlong,dlat,dpsi);
166 } else {
167 gPad->Modified();
168 gPad->Update();
169 }
170 if (issave) {
171 fname = TString::Format("anim%04d.gif", i);
172 gPad->Print(fname);
173 }
174 t += dt;
175 }
177}
178
179////////////////////////////////////////////////////////////////////////////////
180/// Add a point on the track.
181
183{
184 if (!fPoints) {
185 fPointsSize = 16;
187 } else {
188 if (fNpoints>=fPointsSize) {
189 Double_t *temp = new Double_t[2*fPointsSize];
190 memcpy(temp, fPoints, fNpoints*sizeof(Double_t));
191 fPointsSize *= 2;
192 delete [] fPoints;
193 fPoints = temp;
194 }
195 }
196 fPoints[fNpoints++] = x;
197 fPoints[fNpoints++] = y;
198 fPoints[fNpoints++] = z;
199 fPoints[fNpoints++] = t;
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// How-to-browse for a track.
204
206{
207 if (!b) return;
208 Int_t nd = GetNdaughters();
209 if (!nd) {
210 b->Add(this);
211 return;
212 }
213 for (Int_t i=0; i<nd; i++)
214 b->Add(GetDaughter(i));
215
216}
217
218////////////////////////////////////////////////////////////////////////////////
219/// Returns distance to track primitive for picking.
220
222{
223 const Int_t inaxis = 7;
224 const Int_t maxdist = 5;
225 Int_t dist = 9999;
226
227
228 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
229 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
230 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
231 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
232
233 // return if point is not in the user area
234 if (px < puxmin - inaxis) return dist;
235 if (py > puymin + inaxis) return dist;
236 if (px > puxmax + inaxis) return dist;
237 if (py < puymax - inaxis) return dist;
238
239 TView *view = gPad->GetView();
240 if (!view) return dist;
241 Int_t imin, imax;
242 if (TObject::TestBit(kGeoPDrawn) && Size(imin,imax)>=2) {
243 Int_t i, dsegment;
244 Double_t x1,y1,x2,y2;
245 Double_t xndc[3];
246 Int_t np = fNpoints>>2;
247 if (imin<0) imin=0;
248 if (imax>np-1) imax=np-1;
249 for (i=imin;i<imax;i++) {
250 view->WCtoNDC(&fPoints[i<<2], xndc);
251 x1 = xndc[0];
252 y1 = xndc[1];
253 view->WCtoNDC(&fPoints[(i+1)<<2], xndc);
254 x2 = xndc[0];
255 y2 = xndc[1];
256 dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
257// printf("%i: dseg=%i\n", i, dsegment);
258 if (dsegment < dist) {
259 dist = dsegment;
260 if (dist<maxdist) {
261 gPad->SetSelected(this);
262 return 0;
263 }
264 }
265 }
266 }
267 // check now daughters
268 Int_t nd = GetNdaughters();
269 if (!nd) return dist;
270 TGeoTrack *track;
271 for (Int_t id=0; id<nd; id++) {
272 track = (TGeoTrack*)GetDaughter(id);
273 dist = track->DistancetoPrimitive(px,py);
274 if (dist<maxdist) return 0;
275 }
276 return dist;
277}
278
279////////////////////////////////////////////////////////////////////////////////
280/// Draw this track over-imposed on a geometry, according to option.
281/// Options (case sensitive):
282/// - default : track without daughters
283/// - /D : track and first level descendents only
284/// - /* : track and all descendents
285/// - /Ntype : descendents of this track with particle name matching input type.
286///
287/// Options can appear only once but can be combined : e.g. Draw("/D /Npion-")
288///
289/// Time range for visible track segments can be set via TGeoManager::SetTminTmax()
290
292{
294 char *opt1 = Compress(option); // we will have to delete this ?
295 TString opt(opt1);
296 Bool_t is_default = kTRUE;
297 Bool_t is_onelevel = kFALSE;
298 Bool_t is_all = kFALSE;
299 Bool_t is_type = kFALSE;
300 if (opt.Contains("/D")) {
301 is_onelevel = kTRUE;
302 is_default = kFALSE;
303 }
304 if (opt.Contains("/*")) {
305 is_all = kTRUE;
306 is_default = kFALSE;
307 }
308 if (opt.Contains("/N")) {
309 is_type = kTRUE;
310 Int_t ist = opt.Index("/N")+2;
311 Int_t ilast = opt.Index("/",ist);
312 if (ilast<0) ilast=opt.Length();
313 TString type = opt(ist, ilast-ist);
315 }
316 SetBits(is_default, is_onelevel, is_all, is_type);
317 AppendPad("SAME");
319 gPad->Modified();
320 gPad->Update();
321 }
322 delete [] opt1;
323 return;
324}
325
326 ///////////////////////////////////////////////////////////////////////////////
327 /// Event treatment.
328
329void TGeoTrack::ExecuteEvent(Int_t /*event*/, Int_t /*px*/, Int_t /*py*/)
330{
331 if (!gPad) return;
332 gPad->SetCursor(kHand);
333}
334
335////////////////////////////////////////////////////////////////////////////////
336/// Get some info about the track.
337
338char *TGeoTrack::GetObjectInfo(Int_t /*px*/, Int_t /*py*/) const
339{
340 static TString info;
341 Double_t x=0,y=0,z=0,t=0;
342 GetPoint(0,x,y,z,t);
343 info = TString::Format("%s (%g, %g, %g) tof=%g", GetName(),x,y,z,t);
344 return (char*)info.Data();
345}
346
347////////////////////////////////////////////////////////////////////////////////
348/// Get coordinates for point I on the track.
349
351{
352 Int_t np = fNpoints>>2;
353 if (i<0 || i>=np) {
354 Error("GetPoint", "no point %i, indmax=%d", i, np-1);
355 return -1;
356 }
357 Int_t icrt = 4*i;
358 x = fPoints[icrt];
359 y = fPoints[icrt+1];
360 z = fPoints[icrt+2];
361 t = fPoints[icrt+3];
362 return i;
363}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Return the pointer to the array of points starting with index I.
367
369{
370 if (!fNpoints) return 0;
371 return (&fPoints[i<<2]);
372}
373
374////////////////////////////////////////////////////////////////////////////////
375/// Return the index of point on track having closest TOF smaller than
376/// the input value. Output POINT is filled with the interpolated value.
377
379{
380 Int_t np = fNpoints>>2;
381 if (istart>(np-2)) return (np-1);
382 Int_t ip = SearchPoint(tof, istart);
383 if (ip<0 || ip>(np-2)) return ip;
384 // point in segment (ip, ip+1) where 0<=ip<fNpoints-1
385 Int_t i;
386 Int_t j = ip<<2;
387 Int_t k = (ip+1)<<2;
388 Double_t dt = tof-fPoints[j+3];
389 Double_t ddt = fPoints[k+3]-fPoints[j+3];
390 for (i=0; i<3; i++) point[i] = fPoints[j+i] +(fPoints[k+i]-fPoints[j+i])*dt/ddt;
391 return ip;
392}
393
394////////////////////////////////////////////////////////////////////////////////
395/// Paint this track (and descendents) with current attributes.
396
398{
403 Bool_t match_type = kTRUE;
405 if (is_type) {
406 const char *type = gGeoManager->GetParticleName();
407 if (strlen(type) && strcmp(type, GetName())) match_type=kFALSE;
408 }
409 if (match_type) {
410 if (is_default || is_onelevel || is_all) PaintTrack(option);
411 }
412 // paint now daughters
413 Int_t nd = GetNdaughters();
414 if (!nd || is_default) return;
415 TGeoTrack *track;
416 for (Int_t i=0; i<nd; i++) {
417 track = (TGeoTrack*)GetDaughter(i);
418 if (track->IsInTimeRange()) {
419 track->SetBits(is_default,kFALSE,is_all,is_type);
420 track->Paint(option);
421 }
422 }
423}
424
425////////////////////////////////////////////////////////////////////////////////
426/// Paint track and daughters.
427
429{
434 Bool_t match_type = kTRUE;
435 if (is_type) {
436 const char *type = gGeoManager->GetParticleName();
437 if (strlen(type) && strcmp(type, GetName())) match_type=kFALSE;
438 }
439 if (match_type) {
440 if (is_default || is_onelevel || is_all) PaintCollectTrack(time, box);
441 }
442 // loop now daughters
443 Int_t nd = GetNdaughters();
444 if (!nd || is_default) return;
445 TGeoTrack *track;
446 for (Int_t i=0; i<nd; i++) {
447 track = (TGeoTrack*)GetDaughter(i);
448 if (track) track->PaintCollect(time, box);
449 }
450}
451
452////////////////////////////////////////////////////////////////////////////////
453/// Paint just this track.
454
456{
458 if (!painter) return;
459 Int_t np = fNpoints>>2;
460 Double_t point[3], local[3];
462 Int_t ip = GetPoint(time, point);
463 if (ip>=0 && ip<np-1) {
464 if (convert) gGeoManager->MasterToTop(point, local);
465 else memcpy(local, point, 3*sizeof(Double_t));
466 painter->AddTrackPoint(local, box);
467 }
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// Paint current point of the track as marker.
472
474{
475 TPoint p;
476 Double_t xndc[3];
477 TView *view = gPad->GetView();
478 if (!view) return;
479 view->WCtoNDC(point, xndc);
480 if (xndc[0] < gPad->GetX1() || xndc[0] > gPad->GetX2()) return;
481 if (xndc[1] < gPad->GetY1() || xndc[1] > gPad->GetY2()) return;
482 p.fX = gPad->XtoPixel(xndc[0]);
483 p.fY = gPad->YtoPixel(xndc[1]);
485 gVirtualX->DrawPolyMarker(1, &p);
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Paint this track with its current attributes.
490
492{
493 // Check whether there is some 3D view class for this TPad
494// TPadView3D *view3D = (TPadView3D*)gPad->GetView3D();
495// if (view3D) view3D->PaintGeoTrack(this,option); // to be implemented
496
497 // Check if option is 'x3d'. NOTE: This is a simple checking
498 // but since there is no other
499 // options yet, this works fine.
500 TString opt(option);
501 opt.ToLower();
503 if (opt.Contains("x")) return;
504 Int_t np = fNpoints>>2;
505 Int_t imin=0;
506 Int_t imax=np-1;
507 Int_t ip;
508 Double_t start[3] = {0.,0.,0.};
509 Double_t end[3] = {0.,0.,0.};
510 Double_t seg[6] = {0.,0.,0.,0.,0.,0.};
512 Double_t tmin=0.,tmax=0.;
513 Bool_t is_time = gGeoManager->GetTminTmax(tmin,tmax);
514 if (is_time) {
515 imin = GetPoint(tmin, start);
516 if (imin>=0 && imin<np-1) {
517 // we have a starting point -> find ending point
518 imax = GetPoint(tmax, end, imin);
519 if (imax<np-1) {
520 // we also have an ending point -> check if on the same segment with imin
521 if (imax==imin) {
522 // paint the virtual segment between the 2 points
524 if (convert) {
525 gGeoManager->MasterToTop(start, &seg[0]);
526 gGeoManager->MasterToTop(end, &seg[3]);
527 gPad->PaintLine3D(&seg[0], &seg[3]);
528 } else {
529 gPad->PaintLine3D(start, end);
530 }
531 } else {
532 // paint the starting, ending and connecting segments
534 if (convert) {
535 gGeoManager->MasterToTop(start, &seg[0]);
536 gGeoManager->MasterToTop(&fPoints[(imin+1)<<2], &seg[3]);
537 gPad->PaintLine3D(&seg[0], &seg[3]);
538 gGeoManager->MasterToTop(&fPoints[imax<<2], &seg[0]);
539 gGeoManager->MasterToTop(end, &seg[3]);
540 gPad->PaintLine3D(&seg[0], &seg[3]);
541 for (ip=imin+1; ip<imax; ip++) {
542 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
543 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
544 gPad->PaintLine3D(&seg[0], &seg[3]);
545 }
546 } else {
547 gPad->PaintLine3D(start, &fPoints[(imin+1)<<2]);
548 gPad->PaintLine3D(&fPoints[imax<<2], end);
549 for (ip=imin+1; ip<imax; ip++) {
550 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
551 }
552 }
553 }
554 if (convert) {
555 gGeoManager->MasterToTop(end, &seg[0]);
556 PaintMarker(&seg[0]);
557 } else {
558 PaintMarker(end);
559 }
560 } else {
562 if (convert) {
563 gGeoManager->MasterToTop(start, &seg[0]);
564 gGeoManager->MasterToTop(&fPoints[(imin+1)<<2], &seg[3]);
565 gPad->PaintLine3D(&seg[0], &seg[3]);
566 for (ip=imin+1; ip<np-2; ip++) {
567 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
568 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
569 gPad->PaintLine3D(&seg[0], &seg[3]);
570 }
571 } else {
572 gPad->PaintLine3D(start, &fPoints[(imin+1)<<2]);
573 for (ip=imin+1; ip<np-2; ip++) {
574 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
575 }
576 }
577 }
578 } else {
579 imax = GetPoint(tmax, end);
580 if (imax<0 || imax>=(np-1)) return;
581 // we have to draw just the end of the track
583 if (convert) {
584 for (ip=0; ip<imax-1; ip++) {
585 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
586 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
587 gPad->PaintLine3D(&seg[0], &seg[3]);
588 }
589 } else {
590 for (ip=0; ip<imax-1; ip++) {
591 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
592 }
593 }
594 if (convert) {
595 gGeoManager->MasterToTop(&fPoints[imax<<2], &seg[0]);
596 gGeoManager->MasterToTop(end, &seg[3]);
597 gPad->PaintLine3D(&seg[0], &seg[3]);
598 PaintMarker(&seg[3]);
599 } else {
600 gPad->PaintLine3D(&fPoints[imax<<2], end);
601 PaintMarker(end);
602 }
603 }
605 return;
606 }
607
608 // paint all segments from track
610 TAttLine::Modify(); // change attributes if necessary
611 for (ip=imin; ip<imax; ip++) {
612 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
613 }
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Print some info about the track.
618
619void TGeoTrack::Print(Option_t * /*option*/) const
620{
621 Int_t np = fNpoints>>2;
622 printf(" TGeoTrack%6i : %s ===============================\n", fId,GetName());
623 printf(" parent =%6i nd =%3i\n", (fParent)?fParent->GetId():-1, GetNdaughters());
624 Double_t x=0,y=0,z=0,t=0;
625 GetPoint(0,x,y,z,t);
626 printf(" production vertex : (%g, %g, %g) at tof=%g\n", x,y,z,t);
627 GetPoint(np-1,x,y,z,t);
628 printf(" Npoints =%6i, last : (%g, %g, %g) at tof=%g\n\n", np,x,y,z,t);
629}
630
631////////////////////////////////////////////////////////////////////////////////
632/// Return the number of points within the time interval specified by
633/// TGeoManager class and the corresponding indices.
634
636{
637 Double_t tmin, tmax;
638 Int_t np = fNpoints>>2;
639 imin = 0;
640 imax = np-1;
641 Int_t size = np;
642 if (!gGeoManager->GetTminTmax(tmin, tmax)) return size;
643 imin = SearchPoint(tmin);
644 imax = SearchPoint(tmax, imin);
645 return (imax-imin+1);
646}
647
648////////////////////////////////////////////////////////////////////////////////
649/// Search index of track point having the closest time tag smaller than
650/// TIME. Optional start index can be provided.
651
653{
654 Int_t nabove, nbelow, middle, midloc;
655 Int_t np = fNpoints>>2;
656 nabove = np+1;
657 nbelow = istart;
658 while (nabove-nbelow > 1) {
659 middle = (nabove+nbelow)/2;
660 midloc = ((middle-1)<<2)+3;
661 if (time == fPoints[midloc]) return middle-1;
662 if (time < fPoints[midloc]) nabove = middle;
663 else nbelow = middle;
664 }
665 return (nbelow-1);
666}
667
668////////////////////////////////////////////////////////////////////////////////
669/// Set drawing bits for this track
670
671void TGeoTrack::SetBits(Bool_t is_default, Bool_t is_onelevel,
672 Bool_t is_all, Bool_t is_type)
673{
674 TObject::SetBit(kGeoPDefault, is_default);
675 TObject::SetBit(kGeoPOnelevel, is_onelevel);
677 TObject::SetBit(kGeoPType, is_type);
678}
679
680////////////////////////////////////////////////////////////////////////////////
681/// Returns 3D size for the track.
682
684{
685}
686
687////////////////////////////////////////////////////////////////////////////////
688/// Reset data for this track.
689
691{
692 fNpoints = 0;
693 fPointsSize = 0;
694 if (fTracks) {fTracks->Delete(); delete fTracks;}
695 fTracks = 0;
696 if (fPoints) delete [] fPoints;
697 fPoints = 0;
698}
699
@ kHand
Definition GuiTypes.h:374
#define b(i)
Definition RSha256.hxx:100
static const double x2[5]
static const double x1[5]
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
XFontStruct * id
Definition TGX11.cxx:109
int type
Definition TGX11.cxx:121
R__EXTERN TGeoManager * gGeoManager
char * Compress(const char *str)
Remove all blanks from the string str.
Definition TString.cxx:2530
#define gPad
#define gVirtualX
Definition TVirtualX.h:338
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void Modify()
Change current line attributes if necessary.
Definition TAttLine.cxx:245
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition TAttLine.cxx:209
virtual void Modify()
Change current marker attributes if necessary.
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:41
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
void SetParticleName(const char *pname)
TGeoVolume * GetMasterVolume() const
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
void SetVisLevel(Int_t level=3)
set default level down to which visualization is performed
void SetAnimateTracks(Bool_t flag=kTRUE)
Bool_t GetTminTmax(Double_t &tmin, Double_t &tmax) const
Get time cut for drawing tracks.
const char * GetParticleName() const
void SetTminTmax(Double_t tmin=0, Double_t tmax=999)
Set time cut interval for drawing tracks.
Bool_t IsAnimatingTracks() const
TGeoVolume * GetTopVolume() const
void MasterToTop(const Double_t *master, Double_t *top) const
Convert coordinates from master volume frame to top.
Int_t fPointsSize
Definition TGeoTrack.h:39
virtual void Print(Option_t *option="") const
Print some info about the track.
virtual void Sizeof3D() const
Returns 3D size for the track.
virtual void Paint(Option_t *option="")
Paint this track (and descendents) with current attributes.
Int_t Size(Int_t &imin, Int_t &imax)
Return the number of points within the time interval specified by TGeoManager class and the correspon...
virtual void PaintCollectTrack(Double_t time, Double_t *box)
Paint just this track.
Int_t SearchPoint(Double_t time, Int_t istart=0) const
Search index of track point having the closest time tag smaller than TIME.
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Event treatment.
void SetBits(Bool_t is_default=kTRUE, Bool_t is_onelevel=kFALSE, Bool_t is_all=kFALSE, Bool_t is_type=kFALSE)
Set drawing bits for this track.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Returns distance to track primitive for picking.
Double_t * fPoints
Definition TGeoTrack.h:41
virtual ~TGeoTrack()
Destructor.
Definition TGeoTrack.cxx:74
virtual void Draw(Option_t *option="")
Draw this track over-imposed on a geometry, according to option.
TGeoTrack()
Default constructor.
Definition TGeoTrack.cxx:40
virtual char * GetObjectInfo(Int_t px, Int_t py) const
Get some info about the track.
void PaintMarker(Double_t *point, Option_t *option="")
Paint current point of the track as marker.
virtual TVirtualGeoTrack * AddDaughter(Int_t id, Int_t pdgcode, TObject *particle=0)
Add a daughter track to this.
Definition TGeoTrack.cxx:83
@ kGeoPDefault
Definition TGeoTrack.h:31
@ kGeoPOnelevel
Definition TGeoTrack.h:32
@ kGeoPAllDaughters
Definition TGeoTrack.h:33
void Browse(TBrowser *b)
How-to-browse for a track.
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y, Double_t &z, Double_t &t) const
Get coordinates for point I on the track.
virtual void PaintCollect(Double_t time, Double_t *box)
Paint track and daughters.
virtual void ResetTrack()
Reset data for this track.
virtual void AddPoint(Double_t x, Double_t y, Double_t z, Double_t t)
Add a point on the track.
virtual void PaintTrack(Option_t *option="")
Paint this track with its current attributes.
virtual void AnimateTrack(Double_t tmin=0, Double_t tmax=5E-8, Double_t nframes=200, Option_t *option="/*")
Draw animation of this track.
Int_t fNpoints
Definition TGeoTrack.h:40
virtual void Draw(Option_t *option="")
draw top volume according to option
A doubly linked list.
Definition TList.h:38
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition TList.cxx:822
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
virtual void Delete(Option_t *option="")
Remove all objects from the array AND delete all heap based objects.
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:177
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
SCoord_t fY
Definition TPoint.h:36
SCoord_t fX
Definition TPoint.h:35
Basic string class.
Definition TString.h:136
Ssiz_t Length() const
Definition TString.h:410
void ToLower()
Change string to lower-case.
Definition TString.cxx:1150
const char * Data() const
Definition TString.h:369
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2336
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:639
See TView3D.
Definition TView.h:25
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
Abstract class for geometry painters.
virtual void GetViewAngles(Double_t &, Double_t &, Double_t &)
virtual void GrabFocus(Int_t nfr=0, Double_t dlong=0, Double_t dlat=0, Double_t dpsi=0)=0
virtual void AddTrackPoint(Double_t *point, Double_t *box, Bool_t reset=kFALSE)=0
virtual Double_t * GetViewBox()=0
virtual void EstimateCameraMove(Double_t, Double_t, Double_t *, Double_t *)
Base class for user-defined tracks attached to a geometry.
Int_t GetId() const
TVirtualGeoTrack * fParent
void SetParent(TVirtualGeoTrack *parent)
virtual const char * GetName() const
Get the PDG name.
TVirtualGeoTrack * GetDaughter(Int_t index) const
Bool_t IsInTimeRange() const
True if track TOF range overlaps with time interval of TGeoManager.
Int_t GetNdaughters() const
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
th1 Draw()