Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoShape.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
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 TGeoShape
13\ingroup Shapes_classes
14Base abstract class for all shapes.
15
16 Shapes are geometrical objects that provide the basic modelling
17functionality. They provide the definition of the LOCAL frame of coordinates,
18with respect to which they are defined. Any implementation of a shape deriving
19from the base TGeoShape class has to provide methods for :
20
21 - finding out if a point defined in their local frame is or not contained
22 inside;
23 - computing the distance from a local point to getting outside/entering the
24 shape, given a known direction;
25 - computing the maximum distance in any direction from a local point that
26 does NOT result in a boundary crossing of the shape (safe distance);
27 - computing the cosines of the normal vector to the crossed shape surface,
28 given a starting local point and an ongoing direction.
29 All the features above are globally managed by the modeller in order to
30 provide navigation functionality. In addition to those, shapes have also to
31 implement additional specific abstract methods :
32 - computation of the minimal box bounding the shape, given that this box have
33 to be aligned with the local coordinates;
34 - algorithms for dividing the shape along a given axis and producing resulting
35 divisions volumes.
36
37 The modeler currently provides a set of 16 basic shapes, which we will call
38primitives. It also provides a special class allowing the creation of shapes
39made as a result of boolean operations between primitives. These are called
40composite shapes and the composition operation can be recursive (composition
41of composites). This allows the creation of a quite large number of different
42shape topologies and combinations.
43
44 Named shapes register themselves to the manager class at creation time. The
45manager is responsible for their final deletion. Shapes can be created using their
46default constructor if their retrieval by name is not needed, but in this case
47they are owned by the user. A shape may be referenced by several volumes,
48therefore its deletion is not possible once volumes were defined based on it.
49
50### Creating shapes
51
52 Shape objects embed only the minimum set of parameters that are fully
53describing a valid physical shape. For instance, a tube is represented by
54its half length, the minimum radius and the maximum radius. Shapes are used
55together with media in order to create volumes, which in their turn
56are the main components of the geometrical tree. A specific shape can be created
57stand-alone :
58
59~~~ {.cpp}
60 TGeoBBox *box = new TGeoBBox("s_box", halfX, halfY, halfZ); // named
61 TGeoTube *tub = new TGeoTube(rmin, rmax, halfZ); // no name
62 ... (see each specific shape constructors)
63~~~
64
65 Sometimes it is much easier to create a volume having a given shape in one
66step, since shapes are not directly linked in the geometrical tree but volumes
67are :
68
69~~~ {.cpp}
70 TGeoVolume *vol_box = gGeoManager->MakeBox("BOX_VOL", "mat1", halfX, halfY, halfZ);
71 TGeoVolume *vol_tub = gGeoManager->MakeTube("TUB_VOL", "mat2", rmin, rmax, halfZ);
72 ... (see MakeXXX() utilities in TGeoManager class)
73~~~
74
75### Shape queries
76
77Note that global queries related to a geometry are handled by the manager class.
78However, shape-related queries might be sometimes useful.
79
80#### `Bool_t TGeoShape::Contains(const Double_t *point[3])`
81
82this method returns true if POINT is actually inside the shape. The point
83has to be defined in the local shape reference. For instance, for a box having
84DX, DY and DZ half-lengths a point will be considered inside if :
85
86~~~ {.cpp}
87 | -DX <= point[0] <= DX
88 | -DY <= point[1] <= DY
89 | -DZ <= point[2] <= DZ
90~~~
91
92#### `Double_t TGeoShape::DistFromInside(Double_t *point[3], Double_t *dir[3], Int_t iact, Double_t step, Double_t
93*safe)`
94
95computes the distance to exiting a shape from a given point INSIDE, along
96a given direction. The direction is given by its director cosines with respect
97to the local shape coordinate system. This method provides additional
98information according the value of IACT input parameter :
99
100 - IACT = 0 => compute only safe distance and fill it at the location
101 given by SAFE
102 - IACT = 1 => a proposed STEP is supplied. The safe distance is computed
103 first. If this is bigger than STEP than the proposed step
104 is approved and returned by the method since it does not
105 cross the shape boundaries. Otherwise, the distance to
106 exiting the shape is computed and returned.
107 - IACT = 2 => compute both safe distance and distance to exiting, ignoring
108 the proposed step.
109 - IACT > 2 => compute only the distance to exiting, ignoring anything else.
110
111#### `Double_t TGeoShape::DistFromOutside(Double_t *point[3], Double_t *dir[3], Int_t iact, Double_t step, Double_t
112*safe)`
113
114computes the distance to entering a shape from a given point OUTSIDE. Acts
115in the same way as B).
116
117#### `Double_t Safety(const Double_t *point[3], Bool_t inside)`
118
119compute maximum shift of a point in any direction that does not change its
120INSIDE/OUTSIDE state (does not cross shape boundaries). The state of the point
121have to be properly supplied.
122
123#### `Double_t *Normal(Double_t *point[3], Double_t *dir[3], Bool_t inside)`
124
125returns director cosines of normal to the crossed shape surface from a
126given point towards a direction. One has to specify if the point is inside
127or outside shape. According to this, the normal will be outwards or inwards
128shape respectively. Normal components are statically stored by shape class,
129so it has to be copied after retrieval in a different array.
130
131### Dividing shapes
132
133Shapes can generally be divided along a given axis. Supported axis are
134X, Y, Z, Rxy, Phi, Rxyz. A given shape cannot be divided however on any axis.
135The general rule is that that divisions are possible on whatever axis that
136produces still known shapes as slices. The division of shapes should not be
137performed by TGeoShape::Divide() calls, but rather by TGeoVolume::Divide().
138The algorithm for dividing a specific shape is known by the shape object, but
139is always invoked in a generic way from the volume level. Details on how to
140do that can be found in TGeoVolume class. One can see how all division options
141are interpreted and which is their result inside specific shape classes.
142
143\image html geom_t_shape.png width=600px
144*/
145
146#include "TObjArray.h"
147#include "TEnv.h"
148#include "TError.h"
149
150#include "TGeoMatrix.h"
151#include "TGeoManager.h"
152#include "TGeoVolume.h"
153#include "TGeoShape.h"
154#include "TVirtualGeoPainter.h"
155#include "TBuffer3D.h"
156#include "TBuffer3DTypes.h"
157#include "TMath.h"
158
160Double_t TGeoShape::fgEpsMch = 2.220446049250313e-16;
161
162////////////////////////////////////////////////////////////////////////////////
163/// Default constructor
164
166{
167 fShapeBits = 0;
168 fShapeId = 0;
169 if (!gGeoManager) {
170 gGeoManager = new TGeoManager("Geometry", "default geometry");
171 // gROOT->AddGeoManager(gGeoManager);
172 }
173 // fShapeId = gGeoManager->GetListOfShapes()->GetSize();
174 // gGeoManager->AddShape(this);
175}
176
177////////////////////////////////////////////////////////////////////////////////
178/// Default constructor
179
181{
182 fShapeBits = 0;
183 fShapeId = 0;
184 if (!gGeoManager) {
185 gGeoManager = new TGeoManager("Geometry", "default geometry");
186 // gROOT->AddGeoManager(gGeoManager);
187 }
189 gGeoManager->AddShape(this);
190}
191
192////////////////////////////////////////////////////////////////////////////////
193/// Destructor
194
200
201////////////////////////////////////////////////////////////////////////////////
202/// Test for shape navigation methods. Summary for test numbers:
203///
204/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
205/// directions randomly in cos(theta). Compute DistFromInside and move the
206/// point with bigger distance. Compute DistFromOutside back from new point.
207/// Plot d-(d1+d2)
208
210{
211 if (!gGeoManager) {
212 Error("CheckShape", "No geometry manager");
213 return;
214 }
215 TGeoShape *shape = (TGeoShape *)this;
217}
218
219////////////////////////////////////////////////////////////////////////////////
220/// Compute machine round-off double precision error as the smallest number that
221/// if added to 1.0 is different than 1.0.
222
224{
225 Double_t temp1 = 1.0;
226 Double_t temp2 = 1.0 + temp1;
227 Double_t mchEps = 0.;
228 while (temp2 > 1.0) {
229 mchEps = temp1;
230 temp1 /= 2;
231 temp2 = 1.0 + temp1;
232 }
234 return fgEpsMch;
235}
236
237////////////////////////////////////////////////////////////////////////////////
238/// static function returning the machine round-off error
239
241{
242 return fgEpsMch;
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Get the shape name.
247
248const char *TGeoShape::GetName() const
249{
250 if (!fName[0]) {
251 return ((TObject *)this)->ClassName();
252 }
253 return TNamed::GetName();
254}
255
256////////////////////////////////////////////////////////////////////////////////
257/// Returns distance to shape primitive mesh.
258
260{
262 if (!painter)
263 return 9999;
264 return painter->ShapeDistancetoPrimitive(this, numpoints, px, py);
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// Implementation of the inside function using just Contains and GetNormal
269
271{
272 return tgeo_impl::Inside(point, this);
273}
274
275////////////////////////////////////////////////////////////////////////////////
276/// True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,s2
277
278Bool_t
280{
283 if (point[0] * c1 + point[1] * s1 >= 0)
284 saf1 = TMath::Abs(-point[0] * s1 + point[1] * c1);
285 if (point[0] * c2 + point[1] * s2 >= 0)
286 saf2 = TMath::Abs(point[0] * s2 - point[1] * c2);
288 if (saf < epsil)
289 return kTRUE;
290 return kFALSE;
291}
292
293////////////////////////////////////////////////////////////////////////////////
294/// Static method to check if a point is in the phi range (phi1, phi2) [degrees]
295
297{
298 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
299 while (phi < phi1)
300 phi += 360.;
301 Double_t ddp = phi - phi1;
302 if (ddp > phi2 - phi1)
303 return kFALSE;
304 return kTRUE;
305}
306
307////////////////////////////////////////////////////////////////////////////////
308/// Compute distance from POINT to semiplane defined by PHI angle along DIR. Computes
309/// also radius at crossing point. This might be negative in case the crossing is
310/// on the other side of the semiplane.
311
314{
316 Double_t nx = -sphi;
317 Double_t ny = cphi;
318 Double_t rxy0 = point[0] * cphi + point[1] * sphi;
319 Double_t rdotn = point[0] * nx + point[1] * ny;
321 snext = 0.0;
322 rxy = rxy0;
323 return kTRUE;
324 }
325 if (rdotn < 0) {
326 rdotn = -rdotn;
327 } else {
328 nx = -nx;
329 ny = -ny;
330 }
331 Double_t ddotn = dir[0] * nx + dir[1] * ny;
332 if (ddotn <= 0)
333 return kFALSE;
334 snext = rdotn / ddotn;
335 rxy = rxy0 + snext * (dir[0] * cphi + dir[1] * sphi);
336 if (rxy < 0)
337 return kFALSE;
338 return kTRUE;
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Check if two numbers differ with less than a tolerance.
343
345{
346 if (TMath::Abs(a - b) < 1.E-10)
347 return kTRUE;
348 return kFALSE;
349}
350
351////////////////////////////////////////////////////////////////////////////////
352/// Check if segments (A,B) and (C,D) are crossing,
353/// where: A(x1,y1), B(x2,y2), C(x3,y3), D(x4,y4)
354
357{
360 Double_t dx1 = x2 - x1;
362 Double_t dx2 = x4 - x3;
363 Double_t xm = 0.;
364 Double_t ym = 0.;
365 Double_t a1 = 0.;
366 Double_t b1 = 0.;
367 Double_t a2 = 0.;
368 Double_t b2 = 0.;
369 if (TMath::Abs(dx1) < eps)
370 stand1 = kTRUE;
371 if (TMath::Abs(dx2) < eps)
372 stand2 = kTRUE;
373 if (!stand1) {
374 a1 = (x2 * y1 - x1 * y2) / dx1;
375 b1 = (y2 - y1) / dx1;
376 }
377 if (!stand2) {
378 a2 = (x4 * y3 - x3 * y4) / dx2;
379 b2 = (y4 - y3) / dx2;
380 }
381 if (stand1 && stand2) {
382 // Segments parallel and vertical
383 if (TMath::Abs(x1 - x3) < eps) {
384 // Check if segments are overlapping
385 if ((y3 - y1) * (y3 - y2) < -eps || (y4 - y1) * (y4 - y2) < -eps || (y1 - y3) * (y1 - y4) < -eps ||
386 (y2 - y3) * (y2 - y4) < -eps)
387 return kTRUE;
388 return kFALSE;
389 }
390 // Different x values
391 return kFALSE;
392 }
393
394 if (stand1) {
395 // First segment vertical
396 xm = x1;
397 ym = a2 + b2 * xm;
398 } else {
399 if (stand2) {
400 // Second segment vertical
401 xm = x3;
402 ym = a1 + b1 * xm;
403 } else {
404 // Normal crossing
405 if (TMath::Abs(b1 - b2) < eps) {
406 // Parallel segments, are they aligned
407 if (TMath::Abs(y3 - (a1 + b1 * x3)) > eps)
408 return kFALSE;
409 // Aligned segments, are they overlapping
410 if ((x3 - x1) * (x3 - x2) < -eps || (x4 - x1) * (x4 - x2) < -eps || (x1 - x3) * (x1 - x4) < -eps ||
411 (x2 - x3) * (x2 - x4) < -eps)
412 return kTRUE;
413 return kFALSE;
414 }
415 xm = (a1 - a2) / (b2 - b1);
416 ym = (a1 * b2 - a2 * b1) / (b2 - b1);
417 }
418 }
419 // Check if crossing point is both between A,B and C,D
420 Double_t check = (xm - x1) * (xm - x2) + (ym - y1) * (ym - y2);
421 if (check > -eps)
422 return kFALSE;
423 check = (xm - x3) * (xm - x4) + (ym - y3) * (ym - y4);
424 if (check > -eps)
425 return kFALSE;
426 return kTRUE;
427}
428
429////////////////////////////////////////////////////////////////////////////////
430/// compute distance from point (inside phi) to both phi planes. Return minimum.
431
434{
437 Double_t s = 0;
438 Double_t un = dir[0] * s1 - dir[1] * c1;
439 if (!in)
440 un = -un;
441 if (un > 0) {
442 s = -point[0] * s1 + point[1] * c1;
443 if (!in)
444 s = -s;
445 if (s >= 0) {
446 s /= un;
447 if (((point[0] + s * dir[0]) * sm - (point[1] + s * dir[1]) * cm) >= 0)
448 sfi1 = s;
449 }
450 }
451 un = -dir[0] * s2 + dir[1] * c2;
452 if (!in)
453 un = -un;
454 if (un > 0) {
455 s = point[0] * s2 - point[1] * c2;
456 if (!in)
457 s = -s;
458 if (s >= 0) {
459 s /= un;
460 if ((-(point[0] + s * dir[0]) * sm + (point[1] + s * dir[1]) * cm) >= 0)
461 sfi2 = s;
462 }
463 }
464 return TMath::Min(sfi1, sfi2);
465}
466
467////////////////////////////////////////////////////////////////////////////////
468/// Static method to compute normal to phi planes.
469
472{
475 if (point[0] * c1 + point[1] * s1 >= 0)
476 saf1 = TMath::Abs(-point[0] * s1 + point[1] * c1);
477 if (point[0] * c2 + point[1] * s2 >= 0)
478 saf2 = TMath::Abs(point[0] * s2 - point[1] * c2);
479 Double_t c, s;
480 if (saf1 < saf2) {
481 c = c1;
482 s = s1;
483 } else {
484 c = c2;
485 s = s2;
486 }
487 norm[2] = 0;
488 norm[0] = -s;
489 norm[1] = c;
490 if (dir[0] * norm[0] + dir[1] * norm[1] < 0) {
491 norm[0] = s;
492 norm[1] = -c;
493 }
494}
495
496////////////////////////////////////////////////////////////////////////////////
497/// Static method to compute safety w.r.t a phi corner defined by cosines/sines
498/// of the angles phi1, phi2.
499
501{
503 if (inphi && !in)
504 return -TGeoShape::Big();
511 Double_t rsq = point[0] * point[0] + point[1] * point[1];
512 Double_t rproj = point[0] * c1 + point[1] * s1;
514 if (safsq < 0)
515 return 0.;
517 rproj = point[0] * c2 + point[1] * s2;
518 safsq = rsq - rproj * rproj;
519 if (safsq < 0)
520 return 0.;
522 Double_t safe = TMath::Min(saf1, saf2); // >0
523 if (safe > 1E10) {
524 if (in)
525 return TGeoShape::Big();
526 return -TGeoShape::Big();
527 }
528 return safe;
529}
530
531////////////////////////////////////////////////////////////////////////////////
532/// Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
533
535{
536 Double_t crossp = (z2 - z1) * (r - r1) - (z - z1) * (r2 - r1);
537 crossp *= (outer) ? 1. : -1.;
538 // Positive crossp means point on the requested side of the (1,2) segment
539 if (crossp < -TGeoShape::Tolerance()) {
540 // if (((z-z1)*(z2-z)) > -1.E-10) return 0;
541 if (outer)
542 return TGeoShape::Big();
543 else
544 return 0.;
545 }
546 // Compute (1,P) dot (1,2)
547 Double_t c1 = (z - z1) * (z2 - z1) + (r - r1) * (r2 - r1);
548 // Negative c1 means point (1) is closest
549 if (c1 < 1.E-10)
550 return TMath::Sqrt((r - r1) * (r - r1) + (z - z1) * (z - z1));
551 // Compute (2,P) dot (1,2)
552 Double_t c2 = (z - z2) * (z2 - z1) + (r - r2) * (r2 - r1);
553 // Positive c2 means point (2) is closest
554 if (c2 > -1.E-10)
555 return TMath::Sqrt((r - r2) * (r - r2) + (z - z2) * (z - z2));
556 // The closest point is between (1) and (2)
557 c2 = (z2 - z1) * (z2 - z1) + (r2 - r1) * (r2 - r1);
558 // projected length factor with respect to (1,2) length
559 Double_t alpha = c1 / c2;
560 Double_t rp = r1 + alpha * (r2 - r1);
561 Double_t zp = z1 + alpha * (z2 - z1);
562 return TMath::Sqrt((r - rp) * (r - rp) + (z - zp) * (z - zp));
563}
564
565////////////////////////////////////////////////////////////////////////////////
566/// Equivalent of TObject::SetBit.
567
569{
570 if (set) {
571 SetShapeBit(f);
572 } else {
574 }
575}
576
577////////////////////////////////////////////////////////////////////////////////
578/// Returns current transformation matrix that applies to shape.
579
584
585////////////////////////////////////////////////////////////////////////////////
586/// Set current transformation matrix that applies to shape.
587
592
593////////////////////////////////////////////////////////////////////////////////
594/// Tranform a set of points (LocalToMaster)
595
597{
598 UInt_t i, j;
599 Double_t dmaster[3];
600 if (fgTransform) {
601 for (j = 0; j < NbPnts; j++) {
602 i = 3 * j;
603 fgTransform->LocalToMaster(&points[i], dmaster);
604 points[i] = dmaster[0];
605 points[i + 1] = dmaster[1];
606 points[i + 2] = dmaster[2];
607 }
608 return;
609 }
610 if (!gGeoManager)
611 return;
613
614 for (j = 0; j < NbPnts; j++) {
615 i = 3 * j;
618 if (bomb)
619 glmat->LocalToMasterBomb(&points[i], dmaster);
620 else
621 glmat->LocalToMaster(&points[i], dmaster);
622 } else {
623 if (bomb)
625 else
627 }
628 points[i] = dmaster[0];
629 points[i + 1] = dmaster[1];
630 points[i + 2] = dmaster[2];
631 }
632}
633
634////////////////////////////////////////////////////////////////////////////////
635/// Fill the supplied buffer, with sections in desired frame
636/// See TBuffer3D.h for explanation of sections, frame etc.
637
639{
640 // Catch this common potential error here
641 // We have to set kRawSize (unless already done) to allocate buffer space
642 // before kRaw can be filled
646 }
647 }
648
650 // If writing core section all others will be invalid
651 buffer.ClearSectionsValid();
652
653 // Check/grab some objects we need
654 if (!gGeoManager) {
656 return;
657 }
659 if (!paintVolume)
661 if (!paintVolume) {
662 buffer.fID = const_cast<TGeoShape *>(this);
663 buffer.fColor = 0;
664 buffer.fTransparency = 0;
665 // R__ASSERT(kFALSE);
666 // return;
667 } else {
668 buffer.fID = const_cast<TGeoVolume *>(paintVolume);
669 buffer.fColor = paintVolume->GetLineColor();
670
671 buffer.fTransparency = paintVolume->GetTransparency();
673 if (visdensity > 0 && paintVolume->GetMedium()) {
674 if (paintVolume->GetMaterial()->GetDensity() < visdensity) {
675 buffer.fTransparency = 90;
676 }
677 }
678 }
679
680 buffer.fLocalFrame = localFrame;
681 Bool_t r1, r2 = kFALSE;
683 if (paintVolume && paintVolume->GetShape()) {
684 if (paintVolume->GetShape()->IsReflected()) {
685 // Temporary trick to deal with reflected shapes.
686 // Still lighting gets wrong...
687 if (buffer.Type() < TBuffer3DTypes::kTube)
688 r2 = kTRUE;
689 }
690 }
691 buffer.fReflection = ((r1 & (!r2)) | (r2 & !(r1)));
692
693 // Set up local -> master translation matrix
694 if (localFrame) {
695 TGeoMatrix *localMasterMat = nullptr;
698 } else {
700
701 // For overlap drawing the correct matrix needs to obtained in
702 // from GetGLMatrix() - this should not be applied in the case
703 // of composite shapes
706 }
707 }
708 if (!localMasterMat) {
710 return;
711 }
712 localMasterMat->GetHomogenousMatrix(buffer.fLocalMaster);
713 } else {
714 buffer.SetLocalMasterIdentity();
715 }
716
718 }
719}
720
721////////////////////////////////////////////////////////////////////////////////
722/// Get the basic color (0-7).
723
725{
726 Int_t basicColor = 0; // TODO: Check on sensible fallback
727 if (gGeoManager) {
728 const TGeoVolume *volume = gGeoManager->GetPaintVolume();
729 if (volume) {
730 basicColor = ((volume->GetLineColor() % 8) - 1) * 4;
731 if (basicColor < 0)
732 basicColor = 0;
733 }
734 }
735 return basicColor;
736}
737
738////////////////////////////////////////////////////////////////////////////////
739/// Stub implementation to avoid forcing implementation at this stage
740
741const TBuffer3D &TGeoShape::GetBuffer3D(Int_t /*reqSections*/, Bool_t /*localFrame*/) const
742{
743 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
744 Warning("GetBuffer3D", "this must be implemented for shapes in a TGeoPainter hierarchy. This will be come a pure "
745 "virtual fn eventually.");
746 return buffer;
747}
748
749////////////////////////////////////////////////////////////////////////////////
750/// Provide a pointer name containing uid.
751
752const char *TGeoShape::GetPointerName() const
753{
754 static TString name;
755 Int_t uid = GetUniqueID();
756 if (uid)
757 name = TString::Format("p%s_%d", GetName(), uid);
758 else
759 name = TString::Format("p%s", GetName());
760 return name.Data();
761}
762
763////////////////////////////////////////////////////////////////////////////////
764/// Execute mouse actions on this shape.
765
767{
768 if (!gGeoManager)
769 return;
771 painter->ExecuteShapeEvent(this, event, px, py);
772}
773
774////////////////////////////////////////////////////////////////////////////////
775/// Draw this shape.
776
778{
780 if (option && option[0]) {
781 painter->DrawShape(this, option);
782 } else {
783 painter->DrawShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption", ""));
784 }
785}
786
787////////////////////////////////////////////////////////////////////////////////
788/// Paint this shape.
789
791{
793 if (option && option[0]) {
794 painter->PaintShape(this, option);
795 } else {
796 painter->PaintShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption", ""));
797 }
798}
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
Generic 3D primitive description class.
Definition TBuffer3D.h:18
void SetLocalMasterIdentity()
Set kRaw tessellation section of buffer with supplied sizes.
Int_t Type() const
Definition TBuffer3D.h:85
Bool_t SectionsValid(UInt_t mask) const
Definition TBuffer3D.h:67
Double_t fLocalMaster[16]
Definition TBuffer3D.h:93
void ClearSectionsValid()
Clear any sections marked valid.
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Int_t fColor
Definition TBuffer3D.h:88
Short_t fTransparency
Definition TBuffer3D.h:89
Bool_t fReflection
Definition TBuffer3D.h:91
TObject * fID
Definition TBuffer3D.h:87
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:503
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:459
The manager class for any TGeo geometry.
Definition TGeoManager.h:46
TGeoHMatrix * GetGLMatrix() const
Bool_t IsMatrixTransform() const
void LocalToMaster(const Double_t *local, Double_t *master) const
TGeoVolume * GetPaintVolume() const
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
Double_t GetVisDensity() const
void LocalToMasterBomb(const Double_t *local, Double_t *master) const
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
TGeoHMatrix * GetCurrentMatrix() const
Bool_t IsMatrixReflection() const
TVirtualGeoPainter * GetPainter() const
Int_t GetBombMode() const
TObjArray * GetListOfShapes() const
TGeoVolume * GetTopVolume() const
Bool_t IsCleaning() const
Int_t AddShape(const TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
Geometrical transformation package.
Definition TGeoMatrix.h:39
Base abstract class for all shapes.
Definition TGeoShape.h:25
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Stub implementation to avoid forcing implementation at this stage.
UInt_t fShapeBits
Definition TGeoShape.h:80
static Double_t Big()
Definition TGeoShape.h:95
Int_t GetBasicColor() const
Get the basic color (0-7).
static Bool_t IsSegCrossing(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x3, Double_t y3, Double_t x4, Double_t y4)
Check if segments (A,B) and (C,D) are crossing, where: A(x1,y1), B(x2,y2), C(x3,y3),...
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
void ResetShapeBit(UInt_t f)
Definition TGeoShape.h:176
static Double_t DistToPhiMin(const Double_t *point, const Double_t *dir, Double_t s1, Double_t c1, Double_t s2, Double_t c2, Double_t sm, Double_t cm, Bool_t in=kTRUE)
compute distance from point (inside phi) to both phi planes. Return minimum.
TGeoShape()
Default constructor.
static Double_t SafetyPhi(const Double_t *point, Bool_t in, Double_t phi1, Double_t phi2)
Static method to compute safety w.r.t a phi corner defined by cosines/sines of the angles phi1,...
Int_t fShapeId
Definition TGeoShape.h:79
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
virtual Bool_t IsComposite() const
Definition TGeoShape.h:139
void Draw(Option_t *option="") override
Draw this shape.
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
const char * GetPointerName() const
Provide a pointer name containing uid.
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
Fill the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections,...
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Test for shape navigation methods.
static Double_t EpsMch()
static function returning the machine round-off error
virtual EInside Inside(const Double_t *point) const
Implementation of the inside function using just Contains and GetNormal.
static Bool_t IsInPhiRange(const Double_t *point, Double_t phi1, Double_t phi2)
Static method to check if a point is in the phi range (phi1, phi2) [degrees].
static Double_t ComputeEpsMch()
Compute machine round-off double precision error as the smallest number that if added to 1....
void Paint(Option_t *option="") override
Paint this shape.
~TGeoShape() override
Destructor.
static TGeoMatrix * fgTransform
Definition TGeoShape.h:27
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
static Double_t fgEpsMch
Definition TGeoShape.h:28
static void NormalPhi(const Double_t *point, const Double_t *dir, Double_t *norm, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
Static method to compute normal to phi planes.
static Double_t SafetySeg(Double_t r, Double_t z, Double_t r1, Double_t z1, Double_t r2, Double_t z2, Bool_t outer)
Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
static Bool_t IsCrossingSemiplane(const Double_t *point, const Double_t *dir, Double_t cphi, Double_t sphi, Double_t &snext, Double_t &rxy)
Compute distance from POINT to semiplane defined by PHI angle along DIR.
const char * GetName() const override
Get the shape name.
static TGeoMatrix * GetTransform()
Returns current transformation matrix that applies to shape.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute mouse actions on this shape.
static Double_t Tolerance()
Definition TGeoShape.h:98
static Bool_t IsCloseToPhi(Double_t epsil, const Double_t *point, Double_t c1, Double_t s1, Double_t c2, Double_t s2)
True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,...
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
TString fName
Definition TNamed.h:32
TObject * Remove(TObject *obj) override
Remove object from array.
Mother of all ROOT objects.
Definition TObject.h:41
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition TObject.cxx:475
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1074
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1088
Basic string class.
Definition TString.h:138
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:2384
Abstract class for geometry painters.
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:657
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:82
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:75
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:122
TGeoShape::EInside Inside(const Double_t *point, Solid const *solid)
Generic implementation of the inside function using just Contains and GetNormal.
Definition TGeoShape.h:187