Logo ROOT   6.18/05
Reference Guide
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 Geometry_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 *safe)`
93
94computes the distance to exiting a shape from a given point INSIDE, along
95a given direction. The direction is given by its director cosines with respect
96to the local shape coordinate system. This method provides additional
97information according the value of IACT input parameter :
98
99 - IACT = 0 => compute only safe distance and fill it at the location
100 given by SAFE
101 - IACT = 1 => a proposed STEP is supplied. The safe distance is computed
102 first. If this is bigger than STEP than the proposed step
103 is approved and returned by the method since it does not
104 cross the shape boundaries. Otherwise, the distance to
105 exiting the shape is computed and returned.
106 - IACT = 2 => compute both safe distance and distance to exiting, ignoring
107 the proposed step.
108 - IACT > 2 => compute only the distance to exiting, ignoring anything else.
109
110#### `Double_t TGeoShape::DistFromOutside(Double_t *point[3], Double_t *dir[3], Int_t iact, Double_t step, Double_t *safe)`
111
112computes the distance to entering a shape from a given point OUTSIDE. Acts
113in the same way as B).
114
115#### `Double_t Safety(const Double_t *point[3], Bool_t inside)`
116
117compute maximum shift of a point in any direction that does not change its
118INSIDE/OUTSIDE state (does not cross shape boundaries). The state of the point
119have to be properly supplied.
120
121#### `Double_t *Normal(Double_t *point[3], Double_t *dir[3], Bool_t inside)`
122
123returns director cosines of normal to the crossed shape surface from a
124given point towards a direction. One has to specify if the point is inside
125or outside shape. According to this, the normal will be outwards or inwards
126shape respectively. Normal components are statically stored by shape class,
127so it has to be copied after retrieval in a different array.
128
129### Dividing shapes
130
131Shapes can generally be divided along a given axis. Supported axis are
132X, Y, Z, Rxy, Phi, Rxyz. A given shape cannot be divided however on any axis.
133The general rule is that that divisions are possible on whatever axis that
134produces still known shapes as slices. The division of shapes should not be
135performed by TGeoShape::Divide() calls, but rather by TGeoVolume::Divide().
136The algorithm for dividing a specific shape is known by the shape object, but
137is always invoked in a generic way from the volume level. Details on how to
138do that can be found in TGeoVolume class. One can see how all division options
139are interpreted and which is their result inside specific shape classes.
140
141\image html geom_t_shape.png
142*/
143
144#include "TObjArray.h"
145#include "TEnv.h"
146#include "TError.h"
147
148#include "TGeoMatrix.h"
149#include "TGeoManager.h"
150#include "TGeoVolume.h"
151#include "TGeoShape.h"
152#include "TVirtualGeoPainter.h"
153#include "TBuffer3D.h"
154#include "TBuffer3DTypes.h"
155#include "TMath.h"
156
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 :TNamed(name, "")
182{
183 fShapeBits = 0;
184 fShapeId = 0;
185 if (!gGeoManager) {
186 gGeoManager = new TGeoManager("Geometry", "default geometry");
187 // gROOT->AddGeoManager(gGeoManager);
188 }
190 gGeoManager->AddShape(this);
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// Destructor
195
197{
199}
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
209void TGeoShape::CheckShape(Int_t testNo, Int_t nsamples, Option_t *option)
210{
211 if (!gGeoManager) {
212 Error("CheckShape","No geometry manager");
213 return;
214 }
215 TGeoShape *shape = (TGeoShape*)this;
216 gGeoManager->CheckShape(shape, testNo, nsamples, option);
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 }
233 fgEpsMch = mchEps;
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) return 9999;
263 return painter->ShapeDistancetoPrimitive(this, numpoints, px, py);
264}
265
266////////////////////////////////////////////////////////////////////////////////
267/// True if point is closer than epsil to one of the phi planes defined by c1,s1 or c2,s2
268
270{
271 Double_t saf1 = TGeoShape::Big();
272 Double_t saf2 = TGeoShape::Big();
273 if (point[0]*c1+point[1]*s1 >= 0) saf1 = TMath::Abs(-point[0]*s1 + point[1]*c1);
274 if (point[0]*c2+point[1]*s2 >= 0) saf2 = TMath::Abs(point[0]*s2 - point[1]*c2);
275 Double_t saf = TMath::Min(saf1,saf2);
276 if (saf<epsil) return kTRUE;
277 return kFALSE;
278}
279
280////////////////////////////////////////////////////////////////////////////////
281/// Static method to check if a point is in the phi range (phi1, phi2) [degrees]
282
284{
285 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
286 while (phi<phi1) phi+=360.;
287 Double_t ddp = phi-phi1;
288 if (ddp>phi2-phi1) return kFALSE;
289 return kTRUE;
290}
291
292////////////////////////////////////////////////////////////////////////////////
293/// Compute distance from POINT to semiplane defined by PHI angle along DIR. Computes
294/// also radius at crossing point. This might be negative in case the crossing is
295/// on the other side of the semiplane.
296
298{
299 snext = rxy = TGeoShape::Big();
300 Double_t nx = -sphi;
301 Double_t ny = cphi;
302 Double_t rxy0 = point[0]*cphi+point[1]*sphi;
303 Double_t rdotn = point[0]*nx + point[1]*ny;
304 if (TMath::Abs(rdotn)<TGeoShape::Tolerance()) {
305 snext = 0.0;
306 rxy = rxy0;
307 return kTRUE;
308 }
309 if (rdotn<0) {
310 rdotn = -rdotn;
311 } else {
312 nx = -nx;
313 ny = -ny;
314 }
315 Double_t ddotn = dir[0]*nx + dir[1]*ny;
316 if (ddotn<=0) return kFALSE;
317 snext = rdotn/ddotn;
318 rxy = rxy0+snext*(dir[0]*cphi+dir[1]*sphi);
319 if (rxy<0) return kFALSE;
320 return kTRUE;
321}
322
323////////////////////////////////////////////////////////////////////////////////
324/// Check if two numbers differ with less than a tolerance.
325
327{
328 if (TMath::Abs(a-b)<1.E-10) return kTRUE;
329 return kFALSE;
330}
331
332////////////////////////////////////////////////////////////////////////////////
333/// Check if segments (A,B) and (C,D) are crossing,
334/// where: A(x1,y1), B(x2,y2), C(x3,y3), D(x4,y4)
335
337{
339 Bool_t stand1 = kFALSE;
340 Double_t dx1 = x2-x1;
341 Bool_t stand2 = kFALSE;
342 Double_t dx2 = x4-x3;
343 Double_t xm = 0.;
344 Double_t ym = 0.;
345 Double_t a1 = 0.;
346 Double_t b1 = 0.;
347 Double_t a2 = 0.;
348 Double_t b2 = 0.;
349 if (TMath::Abs(dx1) < eps) stand1 = kTRUE;
350 if (TMath::Abs(dx2) < eps) stand2 = kTRUE;
351 if (!stand1) {
352 a1 = (x2*y1-x1*y2)/dx1;
353 b1 = (y2-y1)/dx1;
354 }
355 if (!stand2) {
356 a2 = (x4*y3-x3*y4)/dx2;
357 b2 = (y4-y3)/dx2;
358 }
359 if (stand1 && stand2) {
360 // Segments parallel and vertical
361 if (TMath::Abs(x1-x3)<eps) {
362 // Check if segments are overlapping
363 if ((y3-y1)*(y3-y2)<-eps || (y4-y1)*(y4-y2)<-eps ||
364 (y1-y3)*(y1-y4)<-eps || (y2-y3)*(y2-y4)<-eps) return kTRUE;
365 return kFALSE;
366 }
367 // Different x values
368 return kFALSE;
369 }
370
371 if (stand1) {
372 // First segment vertical
373 xm = x1;
374 ym = a2+b2*xm;
375 } else {
376 if (stand2) {
377 // Second segment vertical
378 xm = x3;
379 ym = a1+b1*xm;
380 } else {
381 // Normal crossing
382 if (TMath::Abs(b1-b2)<eps) {
383 // Parallel segments, are they aligned
384 if (TMath::Abs(y3-(a1+b1*x3))>eps) return kFALSE;
385 // Aligned segments, are they overlapping
386 if ((x3-x1)*(x3-x2)<-eps || (x4-x1)*(x4-x2)<-eps ||
387 (x1-x3)*(x1-x4)<-eps || (x2-x3)*(x2-x4)<-eps) return kTRUE;
388 return kFALSE;
389 }
390 xm = (a1-a2)/(b2-b1);
391 ym = (a1*b2-a2*b1)/(b2-b1);
392 }
393 }
394 // Check if crossing point is both between A,B and C,D
395 Double_t check = (xm-x1)*(xm-x2)+(ym-y1)*(ym-y2);
396 if (check > -eps) return kFALSE;
397 check = (xm-x3)*(xm-x4)+(ym-y3)*(ym-y4);
398 if (check > -eps) return kFALSE;
399 return kTRUE;
400}
401
402////////////////////////////////////////////////////////////////////////////////
403/// compute distance from point (inside phi) to both phi planes. Return minimum.
404
407{
410 Double_t s=0;
411 Double_t un = dir[0]*s1-dir[1]*c1;
412 if (!in) un=-un;
413 if (un>0) {
414 s=-point[0]*s1+point[1]*c1;
415 if (!in) s=-s;
416 if (s>=0) {
417 s /= un;
418 if (((point[0]+s*dir[0])*sm-(point[1]+s*dir[1])*cm)>=0) sfi1=s;
419 }
420 }
421 un = -dir[0]*s2+dir[1]*c2;
422 if (!in) un=-un;
423 if (un>0) {
424 s=point[0]*s2-point[1]*c2;
425 if (!in) s=-s;
426 if (s>=0) {
427 s /= un;
428 if ((-(point[0]+s*dir[0])*sm+(point[1]+s*dir[1])*cm)>=0) sfi2=s;
429 }
430 }
431 return TMath::Min(sfi1, sfi2);
432}
433
434////////////////////////////////////////////////////////////////////////////////
435/// Static method to compute normal to phi planes.
436
438{
439 Double_t saf1 = TGeoShape::Big();
440 Double_t saf2 = TGeoShape::Big();
441 if (point[0]*c1+point[1]*s1 >= 0) saf1 = TMath::Abs(-point[0]*s1 + point[1]*c1);
442 if (point[0]*c2+point[1]*s2 >= 0) saf2 = TMath::Abs(point[0]*s2 - point[1]*c2);
443 Double_t c,s;
444 if (saf1<saf2) {
445 c=c1;
446 s=s1;
447 } else {
448 c=c2;
449 s=s2;
450 }
451 norm[2] = 0;
452 norm[0] = -s;
453 norm[1] = c;
454 if (dir[0]*norm[0]+dir[1]*norm[1] < 0) {
455 norm[0] = s;
456 norm[1] = -c;
457 }
458}
459
460////////////////////////////////////////////////////////////////////////////////
461/// Static method to compute safety w.r.t a phi corner defined by cosines/sines
462/// of the angles phi1, phi2.
463
465{
466 Bool_t inphi = TGeoShape::IsInPhiRange(point, phi1, phi2);
467 if (inphi && !in) return -TGeoShape::Big();
468 phi1 *= TMath::DegToRad();
469 phi2 *= TMath::DegToRad();
470 Double_t c1 = TMath::Cos(phi1);
471 Double_t s1 = TMath::Sin(phi1);
472 Double_t c2 = TMath::Cos(phi2);
473 Double_t s2 = TMath::Sin(phi2);
474 Double_t rsq = point[0]*point[0]+point[1]*point[1];
475 Double_t rproj = point[0]*c1+point[1]*s1;
476 Double_t safsq = rsq-rproj*rproj;
477 if (safsq<0) return 0.;
478 Double_t saf1 = (rproj<0)?TGeoShape::Big():TMath::Sqrt(safsq);
479 rproj = point[0]*c2+point[1]*s2;
480 safsq = rsq-rproj*rproj;
481 if (safsq<0) return 0.;
482 Double_t saf2 = (rproj<0)?TGeoShape::Big():TMath::Sqrt(safsq);
483 Double_t safe = TMath::Min(saf1, saf2); // >0
484 if (safe>1E10) {
485 if (in) return TGeoShape::Big();
486 return -TGeoShape::Big();
487 }
488 return safe;
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Compute distance from point of coordinates (r,z) to segment (r1,z1):(r2,z2)
493
495{
496 Double_t crossp = (z2-z1)*(r-r1)-(z-z1)*(r2-r1);
497 crossp *= (outer) ? 1. : -1.;
498 // Positive crossp means point on the requested side of the (1,2) segment
499 if (crossp < -TGeoShape::Tolerance()) {
500// if (((z-z1)*(z2-z)) > -1.E-10) return 0;
501 if (outer) return TGeoShape::Big();
502 else return 0.;
503 }
504 // Compute (1,P) dot (1,2)
505 Double_t c1 = (z-z1)*(z2-z1)+(r-r1)*(r2-r1);
506 // Negative c1 means point (1) is closest
507 if (c1<1.E-10) return TMath::Sqrt((r-r1)*(r-r1)+(z-z1)*(z-z1));
508 // Compute (2,P) dot (1,2)
509 Double_t c2 = (z-z2)*(z2-z1)+(r-r2)*(r2-r1);
510 // Positive c2 means point (2) is closest
511 if (c2>-1.E-10) return TMath::Sqrt((r-r2)*(r-r2)+(z-z2)*(z-z2));
512 // The closest point is between (1) and (2)
513 c2 = (z2-z1)*(z2-z1)+(r2-r1)*(r2-r1);
514 // projected length factor with respect to (1,2) length
515 Double_t alpha = c1/c2;
516 Double_t rp = r1 + alpha*(r2-r1);
517 Double_t zp = z1 + alpha*(z2-z1);
518 return TMath::Sqrt((r-rp)*(r-rp)+(z-zp)*(z-zp));
519}
520
521////////////////////////////////////////////////////////////////////////////////
522/// Equivalent of TObject::SetBit.
523
525{
526 if (set) {
527 SetShapeBit(f);
528 } else {
530 }
531}
532
533////////////////////////////////////////////////////////////////////////////////
534/// Returns current transformation matrix that applies to shape.
535
537{
538 return fgTransform;
539}
540
541////////////////////////////////////////////////////////////////////////////////
542/// Set current transformation matrix that applies to shape.
543
545{
546 fgTransform = matrix;
547}
548
549////////////////////////////////////////////////////////////////////////////////
550/// Tranform a set of points (LocalToMaster)
551
553{
554 UInt_t i,j;
555 Double_t dmaster[3];
556 if (fgTransform) {
557 for (j = 0; j < NbPnts; j++) {
558 i = 3*j;
559 fgTransform->LocalToMaster(&points[i], dmaster);
560 points[i] = dmaster[0];
561 points[i+1] = dmaster[1];
562 points[i+2] = dmaster[2];
563 }
564 return;
565 }
566 if (!gGeoManager) return;
568
569 for (j = 0; j < NbPnts; j++) {
570 i = 3*j;
573 if (bomb) glmat->LocalToMasterBomb(&points[i], dmaster);
574 else glmat->LocalToMaster(&points[i], dmaster);
575 } else {
576 if (bomb) gGeoManager->LocalToMasterBomb(&points[i], dmaster);
577 else gGeoManager->LocalToMaster(&points[i],dmaster);
578 }
579 points[i] = dmaster[0];
580 points[i+1] = dmaster[1];
581 points[i+2] = dmaster[2];
582 }
583}
584
585////////////////////////////////////////////////////////////////////////////////
586/// Fill the supplied buffer, with sections in desired frame
587/// See TBuffer3D.h for explanation of sections, frame etc.
588
589void TGeoShape::FillBuffer3D(TBuffer3D & buffer, Int_t reqSections, Bool_t localFrame) const
590{
591 // Catch this common potential error here
592 // We have to set kRawSize (unless already done) to allocate buffer space
593 // before kRaw can be filled
594 if (reqSections & TBuffer3D::kRaw) {
595 if (!(reqSections & TBuffer3D::kRawSizes) && !buffer.SectionsValid(TBuffer3D::kRawSizes)) {
597 }
598 }
599
600 if (reqSections & TBuffer3D::kCore) {
601 // If writing core section all others will be invalid
602 buffer.ClearSectionsValid();
603
604 // Check/grab some objects we need
605 if (!gGeoManager) {
607 return;
608 }
609 const TGeoVolume * paintVolume = gGeoManager->GetPaintVolume();
610 if (!paintVolume) paintVolume = gGeoManager->GetTopVolume();
611 if (!paintVolume) {
612 buffer.fID = const_cast<TGeoShape *>(this);
613 buffer.fColor = 0;
614 buffer.fTransparency = 0;
615// R__ASSERT(kFALSE);
616// return;
617 } else {
618 buffer.fID = const_cast<TGeoVolume *>(paintVolume);
619 buffer.fColor = paintVolume->GetLineColor();
620
621 buffer.fTransparency = paintVolume->GetTransparency();
622 Double_t visdensity = gGeoManager->GetVisDensity();
623 if (visdensity>0 && paintVolume->GetMedium()) {
624 if (paintVolume->GetMaterial()->GetDensity() < visdensity) {
625 buffer.fTransparency = 90;
626 }
627 }
628 }
629
630 buffer.fLocalFrame = localFrame;
631 Bool_t r1,r2=kFALSE;
633 if (paintVolume && paintVolume->GetShape()) {
634 if (paintVolume->GetShape()->IsReflected()) {
635 // Temporary trick to deal with reflected shapes.
636 // Still lighting gets wrong...
637 if (buffer.Type() < TBuffer3DTypes::kTube) r2 = kTRUE;
638 }
639 }
640 buffer.fReflection = ((r1&(!r2))|(r2&!(r1)));
641
642 // Set up local -> master translation matrix
643 if (localFrame) {
644 TGeoMatrix * localMasterMat = 0;
646 localMasterMat = TGeoShape::GetTransform();
647 } else {
648 localMasterMat = gGeoManager->GetCurrentMatrix();
649
650 // For overlap drawing the correct matrix needs to obtained in
651 // from GetGLMatrix() - this should not be applied in the case
652 // of composite shapes
654 localMasterMat = gGeoManager->GetGLMatrix();
655 }
656 }
657 if (!localMasterMat) {
659 return;
660 }
661 localMasterMat->GetHomogenousMatrix(buffer.fLocalMaster);
662 } else {
663 buffer.SetLocalMasterIdentity();
664 }
665
667 }
668}
669
670////////////////////////////////////////////////////////////////////////////////
671/// Get the basic color (0-7).
672
674{
675 Int_t basicColor = 0; // TODO: Check on sensible fallback
676 if (gGeoManager) {
677 const TGeoVolume * volume = gGeoManager->GetPaintVolume();
678 if (volume) {
679 basicColor = ((volume->GetLineColor() %8) -1) * 4;
680 if (basicColor < 0) basicColor = 0;
681 }
682 }
683 return basicColor;
684}
685
686////////////////////////////////////////////////////////////////////////////////
687/// Stub implementation to avoid forcing implementation at this stage
688
689const TBuffer3D &TGeoShape::GetBuffer3D(Int_t /*reqSections*/, Bool_t /*localFrame*/) const
690{
691 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
692 Warning("GetBuffer3D", "this must be implemented for shapes in a TGeoPainter hierarchy. This will be come a pure virtual fn eventually.");
693 return buffer;
694}
695
696////////////////////////////////////////////////////////////////////////////////
697/// Provide a pointer name containing uid.
698
699const char *TGeoShape::GetPointerName() const
700{
701 static TString name;
702 Int_t uid = GetUniqueID();
703 if (uid) name = TString::Format("p%s_%d", GetName(),uid);
704 else name = TString::Format("p%s", GetName());
705 return name.Data();
706}
707
708////////////////////////////////////////////////////////////////////////////////
709/// Execute mouse actions on this shape.
710
712{
713 if (!gGeoManager) return;
715 painter->ExecuteShapeEvent(this, event, px, py);
716}
717
718////////////////////////////////////////////////////////////////////////////////
719/// Draw this shape.
720
722{
724 if (option && option[0]) {
725 painter->DrawShape(this, option);
726 } else {
727 painter->DrawShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
728 }
729}
730
731////////////////////////////////////////////////////////////////////////////////
732/// Paint this shape.
733
735{
737 if (option && option[0]) {
738 painter->PaintShape(this, option);
739 } else {
740 painter->PaintShape(this, gEnv->GetValue("Viewer3D.DefaultDrawOption",""));
741 }
742}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define c(i)
Definition: RSha256.hxx:101
#define s1(x)
Definition: RSha256.hxx:91
static const double x2[5]
static const double x4[22]
static const double x1[5]
static const double x3[11]
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
R__EXTERN TEnv * gEnv
Definition: TEnv.h:171
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:601
point * points
Definition: X3DBuffer.c:22
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
Generic 3D primitive description class.
Definition: TBuffer3D.h:18
void SetLocalMasterIdentity()
Set kRaw tessellation section of buffer with supplied sizes.
Definition: TBuffer3D.cxx:296
Int_t Type() const
Definition: TBuffer3D.h:85
Bool_t SectionsValid(UInt_t mask) const
Definition: TBuffer3D.h:67
@ kRawSizes
Definition: TBuffer3D.h:53
Double_t fLocalMaster[16]
Definition: TBuffer3D.h:92
void ClearSectionsValid()
Clear any sections marked valid.
Definition: TBuffer3D.cxx:286
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.
Definition: TCollection.h:182
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition: TEnv.cxx:491
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition: TGeoMatrix.h:421
The manager class for any TGeo geometry.
Definition: TGeoManager.h:43
TGeoHMatrix * GetGLMatrix() const
Definition: TGeoManager.h:516
Bool_t IsMatrixTransform() const
Definition: TGeoManager.h:411
void LocalToMaster(const Double_t *local, Double_t *master) const
Definition: TGeoManager.h:542
TGeoVolume * GetPaintVolume() const
Definition: TGeoManager.h:218
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
Double_t GetVisDensity() const
Definition: TGeoManager.h:220
void LocalToMasterBomb(const Double_t *local, Double_t *master) const
Definition: TGeoManager.h:544
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
TGeoHMatrix * GetCurrentMatrix() const
Definition: TGeoManager.h:515
Bool_t IsMatrixReflection() const
Definition: TGeoManager.h:412
TVirtualGeoPainter * GetPainter() const
Definition: TGeoManager.h:212
Int_t GetBombMode() const
Definition: TGeoManager.h:213
TObjArray * GetListOfShapes() const
Definition: TGeoManager.h:494
TGeoVolume * GetTopVolume() const
Definition: TGeoManager.h:531
Bool_t IsCleaning() const
Definition: TGeoManager.h:483
Int_t AddShape(const TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
virtual Double_t GetDensity() const
Definition: TGeoMaterial.h:103
Geometrical transformation package.
Definition: TGeoMatrix.h:41
virtual void LocalToMasterBomb(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
Definition: TGeoMatrix.cxx:380
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
Definition: TGeoMatrix.cxx:339
void GetHomogenousMatrix(Double_t *hmat) const
The homogenous matrix associated with the transformation is used for piling up's and visualization.
Definition: TGeoMatrix.cxx:314
Base abstract class for all shapes.
Definition: TGeoShape.h:26
virtual const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
Stub implementation to avoid forcing implementation at this stage.
Definition: TGeoShape.cxx:689
UInt_t fShapeBits
Definition: TGeoShape.h:73
static Double_t Big()
Definition: TGeoShape.h:88
Int_t GetBasicColor() const
Get the basic color (0-7).
Definition: TGeoShape.cxx:673
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),...
Definition: TGeoShape.cxx:336
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
Definition: TGeoShape.cxx:552
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
Definition: TGeoShape.cxx:524
void ResetShapeBit(UInt_t f)
Definition: TGeoShape.h:162
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.
Definition: TGeoShape.cxx:405
TGeoShape()
Default constructor.
Definition: TGeoShape.cxx:165
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,...
Definition: TGeoShape.cxx:464
Int_t fShapeId
Definition: TGeoShape.h:72
virtual ~TGeoShape()
Destructor.
Definition: TGeoShape.cxx:196
static void SetTransform(TGeoMatrix *matrix)
Set current transformation matrix that applies to shape.
Definition: TGeoShape.cxx:544
virtual Bool_t IsComposite() const
Definition: TGeoShape.h:130
static Bool_t IsSameWithinTolerance(Double_t a, Double_t b)
Check if two numbers differ with less than a tolerance.
Definition: TGeoShape.cxx:326
const char * GetPointerName() const
Provide a pointer name containing uid.
Definition: TGeoShape.cxx:699
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,...
Definition: TGeoShape.cxx:589
void CheckShape(Int_t testNo, Int_t nsamples=10000, Option_t *option="")
Test for shape navigation methods.
Definition: TGeoShape.cxx:209
static Double_t EpsMch()
static function returning the machine round-off error
Definition: TGeoShape.cxx:240
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].
Definition: TGeoShape.cxx:283
static Double_t ComputeEpsMch()
Compute machine round-off double precision error as the smallest number that if added to 1....
Definition: TGeoShape.cxx:223
static TGeoMatrix * fgTransform
Definition: TGeoShape.h:28
virtual Bool_t IsReflected() const
Definition: TGeoShape.h:137
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
Definition: TGeoShape.cxx:259
static Double_t fgEpsMch
Definition: TGeoShape.h:29
virtual const char * GetName() const
Get the shape name.
Definition: TGeoShape.cxx:248
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.
Definition: TGeoShape.cxx:437
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)
Definition: TGeoShape.cxx:494
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.
Definition: TGeoShape.cxx:297
virtual void Paint(Option_t *option="")
Paint this shape.
Definition: TGeoShape.cxx:734
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute mouse actions on this shape.
Definition: TGeoShape.cxx:711
static TGeoMatrix * GetTransform()
Returns current transformation matrix that applies to shape.
Definition: TGeoShape.cxx:536
static Double_t Tolerance()
Definition: TGeoShape.h:91
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,...
Definition: TGeoShape.cxx:269
virtual void Draw(Option_t *option="")
Draw this shape.
Definition: TGeoShape.cxx:721
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition: TGeoVolume.h:53
TGeoMedium * GetMedium() const
Definition: TGeoVolume.h:176
TGeoMaterial * GetMaterial() const
Definition: TGeoVolume.h:175
TGeoShape * GetShape() const
Definition: TGeoVolume.h:191
Char_t GetTransparency() const
Definition: TGeoVolume.h:190
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
TString fName
Definition: TNamed.h:32
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual TObject * Remove(TObject *obj)
Remove object from array.
Definition: TObjArray.cxx:718
Mother of all ROOT objects.
Definition: TObject.h:37
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition: TObject.cxx:375
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
Basic string class.
Definition: TString.h:131
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:2311
Abstract class for geometry painters.
virtual void DrawShape(TGeoShape *shape, Option_t *option="")=0
virtual void PaintShape(TGeoShape *shape, Option_t *option="")=0
virtual void ExecuteShapeEvent(TGeoShape *shape, Int_t event, Int_t px, Int_t py)=0
virtual Int_t ShapeDistancetoPrimitive(const TGeoShape *shape, Int_t numpoints, Int_t px, Int_t py) const =0
return c1
Definition: legend1.C:41
return c2
Definition: legend2.C:14
static constexpr double s
static constexpr double cm
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:667
constexpr Double_t DegToRad()
Conversion from degree to radian:
Definition: TMath.h:82
Double_t Sqrt(Double_t x)
Definition: TMath.h:679
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
Double_t Cos(Double_t)
Definition: TMath.h:629
Double_t Sin(Double_t)
Definition: TMath.h:625
constexpr Double_t RadToDeg()
Conversion from radian to degree:
Definition: TMath.h:74
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * a
Definition: textangle.C:12
#define snext(osub1, osub2)
Definition: triangle.c:1167