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_classes
23
24Class for user-defined tracks attached to a geometry.
25Tracks are 3D objects made of points and they store a
26pointer to a TParticle. The geometry manager holds a list
27of all tracks that will be deleted on destruction of
28gGeoManager.
29*/
30
32
33////////////////////////////////////////////////////////////////////////////////
34/// Default constructor.
35
37{
38 fPointsSize = 0;
39 fNpoints = 0;
40 fPoints = 0;
41}
42
43////////////////////////////////////////////////////////////////////////////////
44/// Constructor.
45
47 :TVirtualGeoTrack(id,pdgcode,parent,particle)
48{
49 fPointsSize = 0;
50 fNpoints = 0;
51 fPoints = 0;
52 if (fParent==0) {
55 SetMarkerSize(0.6);
56 SetLineColor(2);
57 SetLineWidth(2);
58 } else {
61 SetMarkerSize(0.6);
62 SetLineColor(4);
63 SetLineWidth(2);
64 }
65}
66
67////////////////////////////////////////////////////////////////////////////////
68/// Destructor.
69
71{
72 if (fPoints) delete [] fPoints;
73// if (gPad) gPad->GetListOfPrimitives()->Remove(this);
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Add a daughter track to this.
78
80{
81 if (!fTracks) fTracks = new TObjArray(1);
82 Int_t index = fTracks->GetEntriesFast();
83 TGeoTrack *daughter = new TGeoTrack(id,pdgcode,this,particle);
84 fTracks->AddAtAndExpand(daughter,index);
85 return daughter;
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// Add a daughter and return its index.
90
92{
93 if (!fTracks) fTracks = new TObjArray(1);
94 Int_t index = fTracks->GetEntriesFast();
95 fTracks->AddAtAndExpand(other,index);
96 other->SetParent(this);
97 return index;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101/// Draw animation of this track
102
104{
105 if (tmin<0 || tmin>=tmax || nframes<1) return;
108 if (!gPad) {
110 }
111 TList *list = gPad->GetListOfPrimitives();
112 TIter next(list);
113 TObject *obj;
114 while ((obj = next())) {
115 if (!strcmp(obj->ClassName(), "TGeoTrack")) list->Remove(obj);
116 }
117 Double_t dt = (tmax-tmin)/Double_t(nframes);
118 Double_t delt = 2E-9;
119 Double_t t = tmin;
120 Bool_t geomanim = kFALSE;
121 Bool_t issave = kFALSE;
122 TString fname;
123
124 TString opt(option);
125 if (opt.Contains("/G")) geomanim = kTRUE;
126 if (opt.Contains("/S")) issave = kTRUE;
127
129 Double_t *box = p->GetViewBox();
130 box[0] = box[1] = box[2] = 0;
131 box[3] = box[4] = box[5] = 100;
133 Draw(opt.Data());
134 Double_t start[6], end[6];
135 Int_t i, j;
136 Double_t dlat=0, dlong=0, dpsi=0;
137 Double_t dd[6] = {0,0,0,0,0,0};
138 if (geomanim) {
139 p->EstimateCameraMove(tmin+5*dt, tmin+15*dt, start, end);
140 for (i=0; i<3; i++) {
141 start[i+3] = 20 + 1.3*start[i+3];
142 end[i+3] = 20 + 0.9*end[i+3];
143 }
144 for (i=0; i<6; i++) {
145 dd[i] = (end[i]-start[i])/10.;
146 }
147 memcpy(box, start, 6*sizeof(Double_t));
148 p->GetViewAngles(dlong,dlat,dpsi);
149 dlong = (-206-dlong)/Double_t(nframes);
150 dlat = (126-dlat)/Double_t(nframes);
151 dpsi = (75-dpsi)/Double_t(nframes);
152 p->GrabFocus();
153 }
154
155 for (i=0; i<nframes; i++) {
156 if (t-delt<0) gGeoManager->SetTminTmax(0,t);
157 else gGeoManager->SetTminTmax(t-delt,t);
158 if (geomanim) {
159 for (j=0; j<6; j++) box[j]+=dd[j];
160 p->GrabFocus(1,dlong,dlat,dpsi);
161 } else {
162 gPad->Modified();
163 gPad->Update();
164 }
165 if (issave) {
166 fname = TString::Format("anim%04d.gif", i);
167 gPad->Print(fname);
168 }
169 t += dt;
170 }
172}
173
174////////////////////////////////////////////////////////////////////////////////
175/// Add a point on the track.
176
178{
179 if (!fPoints) {
180 fPointsSize = 16;
182 } else {
183 if (fNpoints>=fPointsSize) {
184 Double_t *temp = new Double_t[2*fPointsSize];
185 memcpy(temp, fPoints, fNpoints*sizeof(Double_t));
186 fPointsSize *= 2;
187 delete [] fPoints;
188 fPoints = temp;
189 }
190 }
191 fPoints[fNpoints++] = x;
192 fPoints[fNpoints++] = y;
193 fPoints[fNpoints++] = z;
194 fPoints[fNpoints++] = t;
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// How-to-browse for a track.
199
201{
202 if (!b) return;
203 Int_t nd = GetNdaughters();
204 if (!nd) {
205 b->Add(this);
206 return;
207 }
208 for (Int_t i=0; i<nd; i++)
209 b->Add(GetDaughter(i));
210
211}
212
213////////////////////////////////////////////////////////////////////////////////
214/// Returns distance to track primitive for picking.
215
217{
218 const Int_t inaxis = 7;
219 const Int_t maxdist = 5;
220 Int_t dist = 9999;
221
222
223 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
224 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
225 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
226 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
227
228 // return if point is not in the user area
229 if (px < puxmin - inaxis) return dist;
230 if (py > puymin + inaxis) return dist;
231 if (px > puxmax + inaxis) return dist;
232 if (py < puymax - inaxis) return dist;
233
234 TView *view = gPad->GetView();
235 if (!view) return dist;
236 Int_t imin, imax;
237 if (TObject::TestBit(kGeoPDrawn) && Size(imin,imax)>=2) {
238 Int_t i, dsegment;
239 Double_t x1,y1,x2,y2;
240 Double_t xndc[3];
241 Int_t np = fNpoints>>2;
242 if (imin<0) imin=0;
243 if (imax>np-1) imax=np-1;
244 for (i=imin;i<imax;i++) {
245 view->WCtoNDC(&fPoints[i<<2], xndc);
246 x1 = xndc[0];
247 y1 = xndc[1];
248 view->WCtoNDC(&fPoints[(i+1)<<2], xndc);
249 x2 = xndc[0];
250 y2 = xndc[1];
251 dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
252// printf("%i: dseg=%i\n", i, dsegment);
253 if (dsegment < dist) {
254 dist = dsegment;
255 if (dist<maxdist) {
256 gPad->SetSelected(this);
257 return 0;
258 }
259 }
260 }
261 }
262 // check now daughters
263 Int_t nd = GetNdaughters();
264 if (!nd) return dist;
265 TGeoTrack *track;
266 for (Int_t id=0; id<nd; id++) {
267 track = (TGeoTrack*)GetDaughter(id);
268 dist = track->DistancetoPrimitive(px,py);
269 if (dist<maxdist) return 0;
270 }
271 return dist;
272}
273
274////////////////////////////////////////////////////////////////////////////////
275/// Draw this track over-imposed on a geometry, according to option.
276/// Options (case sensitive):
277/// - default : track without daughters
278/// - /D : track and first level descendents only
279/// - /* : track and all descendents
280/// - /Ntype : descendents of this track with particle name matching input type.
281///
282/// Options can appear only once but can be combined : e.g. Draw("/D /Npion-")
283///
284/// Time range for visible track segments can be set via TGeoManager::SetTminTmax()
285
287{
289 char *opt1 = Compress(option); // we will have to delete this ?
290 TString opt(opt1);
291 Bool_t is_default = kTRUE;
292 Bool_t is_onelevel = kFALSE;
293 Bool_t is_all = kFALSE;
294 Bool_t is_type = kFALSE;
295 if (opt.Contains("/D")) {
296 is_onelevel = kTRUE;
297 is_default = kFALSE;
298 }
299 if (opt.Contains("/*")) {
300 is_all = kTRUE;
301 is_default = kFALSE;
302 }
303 if (opt.Contains("/N")) {
304 is_type = kTRUE;
305 Int_t ist = opt.Index("/N")+2;
306 Int_t ilast = opt.Index("/",ist);
307 if (ilast<0) ilast=opt.Length();
308 TString type = opt(ist, ilast-ist);
310 }
311 SetBits(is_default, is_onelevel, is_all, is_type);
312 AppendPad("SAME");
314 gPad->Modified();
315 gPad->Update();
316 }
317 delete [] opt1;
318 return;
319}
320
321 ///////////////////////////////////////////////////////////////////////////////
322 /// Event treatment.
323
324void TGeoTrack::ExecuteEvent(Int_t /*event*/, Int_t /*px*/, Int_t /*py*/)
325{
326 if (!gPad) return;
327 gPad->SetCursor(kHand);
328}
329
330////////////////////////////////////////////////////////////////////////////////
331/// Get some info about the track.
332
333char *TGeoTrack::GetObjectInfo(Int_t /*px*/, Int_t /*py*/) const
334{
335 static TString info;
336 Double_t x=0,y=0,z=0,t=0;
337 GetPoint(0,x,y,z,t);
338 info = TString::Format("%s (%g, %g, %g) tof=%g", GetName(),x,y,z,t);
339 return (char*)info.Data();
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Get coordinates for point I on the track.
344
346{
347 Int_t np = fNpoints>>2;
348 if (i<0 || i>=np) {
349 Error("GetPoint", "no point %i, indmax=%d", i, np-1);
350 return -1;
351 }
352 Int_t icrt = 4*i;
353 x = fPoints[icrt];
354 y = fPoints[icrt+1];
355 z = fPoints[icrt+2];
356 t = fPoints[icrt+3];
357 return i;
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// Return the pointer to the array of points starting with index I.
362
364{
365 if (!fNpoints) return 0;
366 return (&fPoints[i<<2]);
367}
368
369////////////////////////////////////////////////////////////////////////////////
370/// Return the index of point on track having closest TOF smaller than
371/// the input value. Output POINT is filled with the interpolated value.
372
374{
375 Int_t np = fNpoints>>2;
376 if (istart>(np-2)) return (np-1);
377 Int_t ip = SearchPoint(tof, istart);
378 if (ip<0 || ip>(np-2)) return ip;
379 // point in segment (ip, ip+1) where 0<=ip<fNpoints-1
380 Int_t i;
381 Int_t j = ip<<2;
382 Int_t k = (ip+1)<<2;
383 Double_t dt = tof-fPoints[j+3];
384 Double_t ddt = fPoints[k+3]-fPoints[j+3];
385 for (i=0; i<3; i++) point[i] = fPoints[j+i] +(fPoints[k+i]-fPoints[j+i])*dt/ddt;
386 return ip;
387}
388
389////////////////////////////////////////////////////////////////////////////////
390/// Paint this track (and descendents) with current attributes.
391
393{
398 Bool_t match_type = kTRUE;
400 if (is_type) {
401 const char *type = gGeoManager->GetParticleName();
402 if (strlen(type) && strcmp(type, GetName())) match_type=kFALSE;
403 }
404 if (match_type) {
405 if (is_default || is_onelevel || is_all) PaintTrack(option);
406 }
407 // paint now daughters
408 Int_t nd = GetNdaughters();
409 if (!nd || is_default) return;
410 TGeoTrack *track;
411 for (Int_t i=0; i<nd; i++) {
412 track = (TGeoTrack*)GetDaughter(i);
413 if (track->IsInTimeRange()) {
414 track->SetBits(is_default,kFALSE,is_all,is_type);
415 track->Paint(option);
416 }
417 }
418}
419
420////////////////////////////////////////////////////////////////////////////////
421/// Paint track and daughters.
422
424{
429 Bool_t match_type = kTRUE;
430 if (is_type) {
431 const char *type = gGeoManager->GetParticleName();
432 if (strlen(type) && strcmp(type, GetName())) match_type=kFALSE;
433 }
434 if (match_type) {
435 if (is_default || is_onelevel || is_all) PaintCollectTrack(time, box);
436 }
437 // loop now daughters
438 Int_t nd = GetNdaughters();
439 if (!nd || is_default) return;
440 TGeoTrack *track;
441 for (Int_t i=0; i<nd; i++) {
442 track = (TGeoTrack*)GetDaughter(i);
443 if (track) track->PaintCollect(time, box);
444 }
445}
446
447////////////////////////////////////////////////////////////////////////////////
448/// Paint just this track.
449
451{
453 if (!painter) return;
454 Int_t np = fNpoints>>2;
455 Double_t point[3], local[3];
457 Int_t ip = GetPoint(time, point);
458 if (ip>=0 && ip<np-1) {
459 if (convert) gGeoManager->MasterToTop(point, local);
460 else memcpy(local, point, 3*sizeof(Double_t));
461 painter->AddTrackPoint(local, box);
462 }
463}
464
465////////////////////////////////////////////////////////////////////////////////
466/// Paint current point of the track as marker.
467
469{
470 TPoint p;
471 Double_t xndc[3];
472 TView *view = gPad->GetView();
473 if (!view) return;
474 view->WCtoNDC(point, xndc);
475 if (xndc[0] < gPad->GetX1() || xndc[0] > gPad->GetX2()) return;
476 if (xndc[1] < gPad->GetY1() || xndc[1] > gPad->GetY2()) return;
477 p.fX = gPad->XtoPixel(xndc[0]);
478 p.fY = gPad->YtoPixel(xndc[1]);
480 gVirtualX->DrawPolyMarker(1, &p);
481}
482
483////////////////////////////////////////////////////////////////////////////////
484/// Paint this track with its current attributes.
485
487{
488 // Check whether there is some 3D view class for this TPad
489// TPadView3D *view3D = (TPadView3D*)gPad->GetView3D();
490// if (view3D) view3D->PaintGeoTrack(this,option); // to be implemented
491
492 // Check if option is 'x3d'. NOTE: This is a simple checking
493 // but since there is no other
494 // options yet, this works fine.
495 TString opt(option);
496 opt.ToLower();
498 if (opt.Contains("x")) return;
499 Int_t np = fNpoints>>2;
500 Int_t imin=0;
501 Int_t imax=np-1;
502 Int_t ip;
503 Double_t start[3] = {0.,0.,0.};
504 Double_t end[3] = {0.,0.,0.};
505 Double_t seg[6] = {0.,0.,0.,0.,0.,0.};
507 Double_t tmin=0.,tmax=0.;
508 Bool_t is_time = gGeoManager->GetTminTmax(tmin,tmax);
509 if (is_time) {
510 imin = GetPoint(tmin, start);
511 if (imin>=0 && imin<np-1) {
512 // we have a starting point -> find ending point
513 imax = GetPoint(tmax, end, imin);
514 if (imax<np-1) {
515 // we also have an ending point -> check if on the same segment with imin
516 if (imax==imin) {
517 // paint the virtual segment between the 2 points
519 if (convert) {
520 gGeoManager->MasterToTop(start, &seg[0]);
521 gGeoManager->MasterToTop(end, &seg[3]);
522 gPad->PaintLine3D(&seg[0], &seg[3]);
523 } else {
524 gPad->PaintLine3D(start, end);
525 }
526 } else {
527 // paint the starting, ending and connecting segments
529 if (convert) {
530 gGeoManager->MasterToTop(start, &seg[0]);
531 gGeoManager->MasterToTop(&fPoints[(imin+1)<<2], &seg[3]);
532 gPad->PaintLine3D(&seg[0], &seg[3]);
533 gGeoManager->MasterToTop(&fPoints[imax<<2], &seg[0]);
534 gGeoManager->MasterToTop(end, &seg[3]);
535 gPad->PaintLine3D(&seg[0], &seg[3]);
536 for (ip=imin+1; ip<imax; ip++) {
537 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
538 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
539 gPad->PaintLine3D(&seg[0], &seg[3]);
540 }
541 } else {
542 gPad->PaintLine3D(start, &fPoints[(imin+1)<<2]);
543 gPad->PaintLine3D(&fPoints[imax<<2], end);
544 for (ip=imin+1; ip<imax; ip++) {
545 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
546 }
547 }
548 }
549 if (convert) {
550 gGeoManager->MasterToTop(end, &seg[0]);
551 PaintMarker(&seg[0]);
552 } else {
553 PaintMarker(end);
554 }
555 } else {
557 if (convert) {
558 gGeoManager->MasterToTop(start, &seg[0]);
559 gGeoManager->MasterToTop(&fPoints[(imin+1)<<2], &seg[3]);
560 gPad->PaintLine3D(&seg[0], &seg[3]);
561 for (ip=imin+1; ip<np-2; ip++) {
562 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
563 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
564 gPad->PaintLine3D(&seg[0], &seg[3]);
565 }
566 } else {
567 gPad->PaintLine3D(start, &fPoints[(imin+1)<<2]);
568 for (ip=imin+1; ip<np-2; ip++) {
569 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
570 }
571 }
572 }
573 } else {
574 imax = GetPoint(tmax, end);
575 if (imax<0 || imax>=(np-1)) return;
576 // we have to draw just the end of the track
578 if (convert) {
579 for (ip=0; ip<imax-1; ip++) {
580 gGeoManager->MasterToTop(&fPoints[ip<<2], &seg[0]);
581 gGeoManager->MasterToTop(&fPoints[(ip+1)<<2], &seg[3]);
582 gPad->PaintLine3D(&seg[0], &seg[3]);
583 }
584 } else {
585 for (ip=0; ip<imax-1; ip++) {
586 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
587 }
588 }
589 if (convert) {
590 gGeoManager->MasterToTop(&fPoints[imax<<2], &seg[0]);
591 gGeoManager->MasterToTop(end, &seg[3]);
592 gPad->PaintLine3D(&seg[0], &seg[3]);
593 PaintMarker(&seg[3]);
594 } else {
595 gPad->PaintLine3D(&fPoints[imax<<2], end);
596 PaintMarker(end);
597 }
598 }
600 return;
601 }
602
603 // paint all segments from track
605 TAttLine::Modify(); // change attributes if necessary
606 for (ip=imin; ip<imax; ip++) {
607 gPad->PaintLine3D(&fPoints[ip<<2], &fPoints[(ip+1)<<2]);
608 }
609}
610
611////////////////////////////////////////////////////////////////////////////////
612/// Print some info about the track.
613
614void TGeoTrack::Print(Option_t * /*option*/) const
615{
616 Int_t np = fNpoints>>2;
617 printf(" TGeoTrack%6i : %s ===============================\n", fId,GetName());
618 printf(" parent =%6i nd =%3i\n", (fParent)?fParent->GetId():-1, GetNdaughters());
619 Double_t x=0,y=0,z=0,t=0;
620 GetPoint(0,x,y,z,t);
621 printf(" production vertex : (%g, %g, %g) at tof=%g\n", x,y,z,t);
622 GetPoint(np-1,x,y,z,t);
623 printf(" Npoints =%6i, last : (%g, %g, %g) at tof=%g\n\n", np,x,y,z,t);
624}
625
626////////////////////////////////////////////////////////////////////////////////
627/// Return the number of points within the time interval specified by
628/// TGeoManager class and the corresponding indices.
629
631{
632 Double_t tmin, tmax;
633 Int_t np = fNpoints>>2;
634 imin = 0;
635 imax = np-1;
636 Int_t size = np;
637 if (!gGeoManager->GetTminTmax(tmin, tmax)) return size;
638 imin = SearchPoint(tmin);
639 imax = SearchPoint(tmax, imin);
640 return (imax-imin+1);
641}
642
643////////////////////////////////////////////////////////////////////////////////
644/// Search index of track point having the closest time tag smaller than
645/// TIME. Optional start index can be provided.
646
648{
649 Int_t nabove, nbelow, middle, midloc;
650 Int_t np = fNpoints>>2;
651 nabove = np+1;
652 nbelow = istart;
653 while (nabove-nbelow > 1) {
654 middle = (nabove+nbelow)/2;
655 midloc = ((middle-1)<<2)+3;
656 if (time == fPoints[midloc]) return middle-1;
657 if (time < fPoints[midloc]) nabove = middle;
658 else nbelow = middle;
659 }
660 return (nbelow-1);
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// Set drawing bits for this track
665
666void TGeoTrack::SetBits(Bool_t is_default, Bool_t is_onelevel,
667 Bool_t is_all, Bool_t is_type)
668{
669 TObject::SetBit(kGeoPDefault, is_default);
670 TObject::SetBit(kGeoPOnelevel, is_onelevel);
672 TObject::SetBit(kGeoPType, is_type);
673}
674
675////////////////////////////////////////////////////////////////////////////////
676/// Returns 3D size for the track.
677
679{
680}
681
682////////////////////////////////////////////////////////////////////////////////
683/// Reset data for this track.
684
686{
687 fNpoints = 0;
688 fPointsSize = 0;
689 if (fTracks) {fTracks->Delete(); delete fTracks;}
690 fTracks = 0;
691 if (fPoints) delete [] fPoints;
692 fPoints = 0;
693}
694
@ kHand
Definition GuiTypes.h:374
#define b(i)
Definition RSha256.hxx:100
static const double x2[5]
static const double x1[5]
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:91
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:2524
#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:242
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:206
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.
Class for user-defined tracks attached to a geometry.
Definition TGeoTrack.h:27
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:70
virtual void Draw(Option_t *option="")
Draw this track over-imposed on a geometry, according to option.
TGeoTrack()
Default constructor.
Definition TGeoTrack.cxx:36
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:79
@ 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:44
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition TList.cxx:822
An array of TObjects.
Definition TObjArray.h:37
Int_t GetEntriesFast() const
Definition TObjArray.h:64
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:37
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:187
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:130
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:107
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:696
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:893
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:1145
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:2331
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()