Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoBBox.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$// Author: Andrei Gheata 24/10/01
2
3// Contains() and DistFromOutside/Out() implemented by Mihaela Gheata
4
5/*************************************************************************
6 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
13/** \class TGeoBBox
14\ingroup Shapes_classes
15\brief Box class.
16
17 - [Building boxes](\ref GEOB00)
18 - [Creation of boxes](\ref GEOB01)
19 - [Divisions of boxes](\ref GEOB02)
20
21All shape primitives inherit from this, their
22 constructor filling automatically the parameters of the box that bounds
23 the given shape. Defined by 6 parameters :
24
25```
26TGeoBBox(Double_t dx,Double_t dy,Double_t dz,Double_t *origin=0);
27```
28
29 - `fDX`, `fDY`, `fDZ` : half lengths on X, Y and Z axis
30 - `fOrigin[3]` : position of box origin
31
32\anchor GEOB00
33### Building boxes
34
35Normally a box has to be built only with 3 parameters: `DX,DY,DZ`
36representing the half-lengths on X, Y and Z-axes. In this case, the
37origin of the box will match the one of its reference frame and the box
38will range from: `-DX` to `DX` on X-axis, from `-DY` to `DY` on Y and
39from `-DZ` to `DZ` on Z. On the other hand, any other shape needs to
40compute and store the parameters of their minimal bounding box. The
41bounding boxes are essential to optimize navigation algorithms.
42Therefore all other primitives derive from **`TGeoBBox`**. Since the
43minimal bounding box is not necessary centered in the origin, any box
44allows an origin translation `(Ox`,`Oy`,`Oz)`. All primitive
45constructors automatically compute the bounding box parameters. Users
46should be aware that building a translated box that will represent a
47primitive shape by itself would affect any further positioning of other
48shapes inside. Therefore it is highly recommendable to build
49non-translated boxes as primitives and translate/rotate their
50corresponding volumes only during positioning stage.
51
52\anchor GEOB01
53#### Creation of boxes
54
55```
56 TGeoBBox *box = new TGeoBBox("BOX", 20, 30, 40);
57```
58
59Begin_Macro
60{
61 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
62 new TGeoManager("box", "poza1");
63 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
64 TGeoMedium *med = new TGeoMedium("MED",1,mat);
65 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
66 gGeoManager->SetTopVolume(top);
67 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
68 vol->SetLineWidth(2);
69 top->AddNode(vol,1);
70 gGeoManager->CloseGeometry();
71 gGeoManager->SetNsegments(80);
72 top->Draw();
73 TView *view = gPad->GetView();
74 view->ShowAxis();
75}
76End_Macro
77
78A volume having a box shape can be built in one step:
79
80```
81 TGeoVolume *vbox = gGeoManager->MakeBox("vbox", ptrMed, 20,30,40);
82```
83
84\anchor GEOB02
85#### Divisions of boxes
86
87 Volumes having box shape can be divided with equal-length slices on
88X, Y or Z axis. The following options are supported:
89
90 - Dividing the full range of one axis in N slices
91```
92 TGeoVolume *divx = vbox->Divide("SLICEX", 1, N);
93```
94here `1` stands for the division axis (1-X, 2-Y, 3-Z)
95
96Begin_Macro
97{
98 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
99 new TGeoManager("box", "poza1");
100 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
101 TGeoMedium *med = new TGeoMedium("MED",1,mat);
102 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
103 gGeoManager->SetTopVolume(top);
104 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
105 vol->SetLineWidth(2);
106 top->AddNode(vol,1);
107 TGeoVolume *divx = vol->Divide("SLICE",1,8,0,0);
108 gGeoManager->CloseGeometry();
109 gGeoManager->SetNsegments(80);
110 top->Draw();
111 TView *view = gPad->GetView();
112 view->ShowAxis();
113}
114End_Macro
115
116 - Dividing in a limited range - general case.
117```
118 TGeoVolume *divy = vbox->Divide("SLICEY",2,N,start,step);
119```
120 - start = starting offset within (-fDY, fDY)
121 - step = slicing step
122
123Begin_Macro
124{
125 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
126 new TGeoManager("box", "poza1");
127 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
128 TGeoMedium *med = new TGeoMedium("MED",1,mat);
129 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
130 gGeoManager->SetTopVolume(top);
131 TGeoVolume *vol = gGeoManager->MakeBox("BOX",med, 20,30,40);
132 vol->SetLineWidth(2);
133 top->AddNode(vol,1);
134 TGeoVolume *divx = vol->Divide("SLICE",2,8,2,3);
135 gGeoManager->CloseGeometry();
136 gGeoManager->SetNsegments(80);
137 top->Draw();
138 TView *view = gPad->GetView();
139 view->ShowAxis();
140}
141End_Macro
142
143Both cases are supported by all shapes.
144See also class TGeoShape for utility methods provided by any particular shape.
145*/
146
147#include <iostream>
148
149#include "TGeoManager.h"
150#include "TGeoMatrix.h"
151#include "TGeoVolume.h"
152#include "TVirtualGeoPainter.h"
153#include "TGeoBBox.h"
154#include "TBuffer3D.h"
155#include "TBuffer3DTypes.h"
156#include "TMath.h"
157#include "TRandom.h"
158
160
161////////////////////////////////////////////////////////////////////////////////
162/// Default constructor
163
165{
167 fDX = fDY = fDZ = 0;
168 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
169}
170
171////////////////////////////////////////////////////////////////////////////////
172/// Constructor where half-lengths are provided.
173
175{
176 SetShapeBit(TGeoShape::kGeoBox);
177 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
178 SetBoxDimensions(dx, dy, dz, origin);
179}
180
181////////////////////////////////////////////////////////////////////////////////
182/// Constructor with shape name.
183
184TGeoBBox::TGeoBBox(const char *name, Double_t dx, Double_t dy, Double_t dz, Double_t *origin) : TGeoShape(name)
185{
186 SetShapeBit(TGeoShape::kGeoBox);
187 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
188 SetBoxDimensions(dx, dy, dz, origin);
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// Constructor based on the array of parameters.
193/// - param[0] - half-length in x
194/// - param[1] - half-length in y
195/// - param[2] - half-length in z
196
198{
199 SetShapeBit(TGeoShape::kGeoBox);
200 fOrigin[0] = fOrigin[1] = fOrigin[2] = 0.0;
201 SetDimensions(param);
202}
203
204////////////////////////////////////////////////////////////////////////////////
205/// Destructor
206
208
209////////////////////////////////////////////////////////////////////////////////
210/// Computes capacity of the shape in [length^3].
211
213{
214 return (8. * fDX * fDY * fDZ);
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// Computes normal to closest surface from POINT.
219
220void TGeoBBox::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
221{
222 memset(norm, 0, 3 * sizeof(Double_t));
223 Double_t saf[3];
224 Int_t i;
225 saf[0] = TMath::Abs(TMath::Abs(point[0] - fOrigin[0]) - fDX);
226 saf[1] = TMath::Abs(TMath::Abs(point[1] - fOrigin[1]) - fDY);
227 saf[2] = TMath::Abs(TMath::Abs(point[2] - fOrigin[2]) - fDZ);
228 i = TMath::LocMin(3, saf);
229 norm[i] = (dir[i] > 0) ? 1 : (-1);
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// Decides fast if the bounding box could be crossed by a vector.
234
235Bool_t TGeoBBox::CouldBeCrossed(const Double_t *point, const Double_t *dir) const
236{
237 Double_t mind = fDX;
238 if (fDY < mind)
239 mind = fDY;
240 if (fDZ < mind)
241 mind = fDZ;
242 Double_t dx = fOrigin[0] - point[0];
243 Double_t dy = fOrigin[1] - point[1];
244 Double_t dz = fOrigin[2] - point[2];
245 Double_t do2 = dx * dx + dy * dy + dz * dz;
246 if (do2 <= (mind * mind))
247 return kTRUE;
248 Double_t rmax2 = fDX * fDX + fDY * fDY + fDZ * fDZ;
249 if (do2 <= rmax2)
250 return kTRUE;
251 // inside bounding sphere
252 Double_t doct = dx * dir[0] + dy * dir[1] + dz * dir[2];
253 // leaving ray
254 if (doct <= 0)
255 return kFALSE;
256 Double_t dirnorm = dir[0] * dir[0] + dir[1] * dir[1] + dir[2] * dir[2];
257 if ((doct * doct) >= (do2 - rmax2) * dirnorm)
258 return kTRUE;
259 return kFALSE;
260}
261
262////////////////////////////////////////////////////////////////////////////////
263/// Compute closest distance from point px,py to each corner.
264
266{
267 const Int_t numPoints = 8;
268 return ShapeDistancetoPrimitive(numPoints, px, py);
269}
270
271////////////////////////////////////////////////////////////////////////////////
272/// Divide this box shape belonging to volume "voldiv" into ndiv equal volumes
273/// called divname, from start position with the given step. Returns pointer
274/// to created division cell volume. In case a wrong division axis is supplied,
275/// returns pointer to volume to be divided.
276
278TGeoBBox::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
279{
280 TGeoShape *shape; //--- shape to be created
281 TGeoVolume *vol; //--- division volume to be created
282 TGeoVolumeMulti *vmulti; //--- generic divided volume
283 TGeoPatternFinder *finder; //--- finder to be attached
284 TString opt = ""; //--- option to be attached
285 Double_t end = start + ndiv * step;
286 switch (iaxis) {
287 case 1: //--- divide on X
288 shape = new TGeoBBox(step / 2., fDY, fDZ);
289 finder = new TGeoPatternX(voldiv, ndiv, start, end);
290 opt = "X";
291 break;
292 case 2: //--- divide on Y
293 shape = new TGeoBBox(fDX, step / 2., fDZ);
294 finder = new TGeoPatternY(voldiv, ndiv, start, end);
295 opt = "Y";
296 break;
297 case 3: //--- divide on Z
298 shape = new TGeoBBox(fDX, fDY, step / 2.);
299 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
300 opt = "Z";
301 break;
302 default: Error("Divide", "Wrong axis type for division"); return nullptr;
303 }
304 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
305 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
306 vmulti->AddVolume(vol);
307 voldiv->SetFinder(finder);
308 finder->SetDivIndex(voldiv->GetNdaughters());
309 for (Int_t ic = 0; ic < ndiv; ic++) {
310 voldiv->AddNodeOffset(vol, ic, start + step / 2. + ic * step, opt.Data());
311 ((TGeoNodeOffset *)voldiv->GetNodes()->At(voldiv->GetNdaughters() - 1))->SetFinder(finder);
312 }
313 return vmulti;
314}
315
316////////////////////////////////////////////////////////////////////////////////
317/// Compute bounding box - nothing to do in this case.
318
319void TGeoBBox::ComputeBBox() {}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Test if point is inside this shape.
323
324Bool_t TGeoBBox::Contains(const Double_t *point) const
325{
326 if (TMath::Abs(point[2] - fOrigin[2]) > fDZ)
327 return kFALSE;
328 if (TMath::Abs(point[0] - fOrigin[0]) > fDX)
329 return kFALSE;
330 if (TMath::Abs(point[1] - fOrigin[1]) > fDY)
331 return kFALSE;
332 return kTRUE;
333}
334
335////////////////////////////////////////////////////////////////////////////////
336/// Static method to check if point[3] is located inside a box of having dx, dy, dz
337/// as half-lengths.
338
339Bool_t TGeoBBox::Contains(const Double_t *point, Double_t dx, Double_t dy, Double_t dz, const Double_t *origin)
340{
341 if (TMath::Abs(point[2] - origin[2]) > dz)
342 return kFALSE;
343 if (TMath::Abs(point[0] - origin[0]) > dx)
344 return kFALSE;
345 if (TMath::Abs(point[1] - origin[1]) > dy)
346 return kFALSE;
347 return kTRUE;
348}
349
350////////////////////////////////////////////////////////////////////////////////
351/// Compute distance from inside point to surface of the box.
352/// Boundary safe algorithm.
353
355TGeoBBox::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
356{
357 Double_t s, smin, saf[6];
358 Double_t newpt[3];
359 Int_t i;
360 for (i = 0; i < 3; i++)
361 newpt[i] = point[i] - fOrigin[i];
362 saf[0] = fDX + newpt[0];
363 saf[1] = fDX - newpt[0];
364 saf[2] = fDY + newpt[1];
365 saf[3] = fDY - newpt[1];
366 saf[4] = fDZ + newpt[2];
367 saf[5] = fDZ - newpt[2];
368 if (iact < 3 && safe) {
369 smin = saf[0];
370 // compute safe distance
371 for (i = 1; i < 6; i++)
372 if (saf[i] < smin)
373 smin = saf[i];
374 *safe = smin;
375 if (smin < 0)
376 *safe = 0.0;
377 if (iact == 0)
378 return TGeoShape::Big();
379 if (iact == 1 && step < *safe)
380 return TGeoShape::Big();
381 }
382 // compute distance to surface
383 smin = TGeoShape::Big();
384 for (i = 0; i < 3; i++) {
385 if (dir[i] != 0) {
386 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
387 if (s < 0)
388 return 0.0;
389 if (s < smin)
390 smin = s;
391 }
392 }
393 return smin;
394}
395
396////////////////////////////////////////////////////////////////////////////////
397/// Compute distance from inside point to surface of the box.
398/// Boundary safe algorithm.
399
400Double_t TGeoBBox::DistFromInside(const Double_t *point, const Double_t *dir, Double_t dx, Double_t dy, Double_t dz,
401 const Double_t *origin, Double_t /*stepmax*/)
402{
403 Double_t s, smin, saf[6];
404 Double_t newpt[3];
405 Int_t i;
406 for (i = 0; i < 3; i++)
407 newpt[i] = point[i] - origin[i];
408 saf[0] = dx + newpt[0];
409 saf[1] = dx - newpt[0];
410 saf[2] = dy + newpt[1];
411 saf[3] = dy - newpt[1];
412 saf[4] = dz + newpt[2];
413 saf[5] = dz - newpt[2];
414 // compute distance to surface
415 smin = TGeoShape::Big();
416 for (i = 0; i < 3; i++) {
417 if (dir[i] != 0) {
418 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
419 if (s < 0)
420 return 0.0;
421 if (s < smin)
422 smin = s;
423 }
424 }
425 return smin;
426}
427
428////////////////////////////////////////////////////////////////////////////////
429/// Compute distance from outside point to surface of the box.
430/// Boundary safe algorithm.
431
433TGeoBBox::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
434{
435 Bool_t in = kTRUE;
436 Double_t saf[3];
437 Double_t par[3];
438 Double_t newpt[3];
439 Int_t i, j;
440 for (i = 0; i < 3; i++)
441 newpt[i] = point[i] - fOrigin[i];
442 par[0] = fDX;
443 par[1] = fDY;
444 par[2] = fDZ;
445 for (i = 0; i < 3; i++) {
446 saf[i] = TMath::Abs(newpt[i]) - par[i];
447 if (saf[i] >= step)
448 return TGeoShape::Big();
449 if (in && saf[i] > 0)
450 in = kFALSE;
451 }
452 if (iact < 3 && safe) {
453 // compute safe distance
454 if (in) {
455 *safe = 0.0;
456 } else {
457 *safe = saf[0];
458 if (saf[1] > *safe)
459 *safe = saf[1];
460 if (saf[2] > *safe)
461 *safe = saf[2];
462 }
463 if (iact == 0)
464 return TGeoShape::Big();
465 if (iact == 1 && step < *safe)
466 return TGeoShape::Big();
467 }
468 // compute distance from point to box
469 // protection in case point is actually inside box
470 if (in) {
471 j = 0;
472 Double_t ss = saf[0];
473 if (saf[1] > ss) {
474 ss = saf[1];
475 j = 1;
476 }
477 if (saf[2] > ss)
478 j = 2;
479 if (newpt[j] * dir[j] > 0)
480 return TGeoShape::Big(); // in fact exiting
481 return 0.0;
482 }
483 for (i = 0; i < 3; i++) {
484 if (saf[i] < 0)
485 continue;
486 if (newpt[i] * dir[i] >= 0)
487 continue;
488 Double_t snxt = saf[i] / TMath::Abs(dir[i]);
489 Int_t ibreak = 0;
490 for (j = 0; j < 3; j++) {
491 if (j == i)
492 continue;
493 Double_t coord = newpt[j] + snxt * dir[j];
494 if (TMath::Abs(coord) > par[j]) {
495 ibreak = 1;
496 break;
497 }
498 }
499 if (!ibreak)
500 return snxt;
501 }
502 return TGeoShape::Big();
503}
504
505////////////////////////////////////////////////////////////////////////////////
506/// Compute distance from outside point to surface of the box.
507/// Boundary safe algorithm.
508
510 const Double_t *origin, Double_t stepmax)
511{
512 Bool_t in = kTRUE;
513 Double_t saf[3];
514 Double_t par[3];
515 Double_t newpt[3];
516 Int_t i, j;
517 for (i = 0; i < 3; i++)
518 newpt[i] = point[i] - origin[i];
519 par[0] = dx;
520 par[1] = dy;
521 par[2] = dz;
522 for (i = 0; i < 3; i++) {
523 saf[i] = TMath::Abs(newpt[i]) - par[i];
524 if (saf[i] >= stepmax)
525 return TGeoShape::Big();
526 if (in && saf[i] > 0)
527 in = kFALSE;
528 }
529 // In case point is inside return ZERO
530 if (in)
531 return 0.0;
532 Double_t coord, snxt = TGeoShape::Big();
533 Int_t ibreak = 0;
534 for (i = 0; i < 3; i++) {
535 if (saf[i] < 0)
536 continue;
537 if (newpt[i] * dir[i] >= 0)
538 continue;
539 snxt = saf[i] / TMath::Abs(dir[i]);
540 ibreak = 0;
541 for (j = 0; j < 3; j++) {
542 if (j == i)
543 continue;
544 coord = newpt[j] + snxt * dir[j];
545 if (TMath::Abs(coord) > par[j]) {
546 ibreak = 1;
547 break;
548 }
549 }
550 if (!ibreak)
551 return snxt;
552 }
553 return TGeoShape::Big();
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// Returns name of axis IAXIS.
558
559const char *TGeoBBox::GetAxisName(Int_t iaxis) const
560{
561 switch (iaxis) {
562 case 1: return "X";
563 case 2: return "Y";
564 case 3: return "Z";
565 default: return "UNDEFINED";
566 }
567}
568
569////////////////////////////////////////////////////////////////////////////////
570/// Get range of shape for a given axis.
571
573{
574 xlo = 0;
575 xhi = 0;
576 Double_t dx = 0;
577 switch (iaxis) {
578 case 1:
579 xlo = fOrigin[0] - fDX;
580 xhi = fOrigin[0] + fDX;
581 dx = 2 * fDX;
582 return dx;
583 case 2:
584 xlo = fOrigin[1] - fDY;
585 xhi = fOrigin[1] + fDY;
586 dx = 2 * fDY;
587 return dx;
588 case 3:
589 xlo = fOrigin[2] - fDZ;
590 xhi = fOrigin[2] + fDZ;
591 dx = 2 * fDZ;
592 return dx;
593 }
594 return dx;
595}
596
597////////////////////////////////////////////////////////////////////////////////
598/// Fill vector param[4] with the bounding cylinder parameters. The order
599/// is the following : Rmin, Rmax, Phi1, Phi2
600
602{
603 param[0] = 0.; // Rmin
604 param[1] = fDX * fDX + fDY * fDY; // Rmax
605 param[2] = 0.; // Phi1
606 param[3] = 360.; // Phi2
607}
608
609////////////////////////////////////////////////////////////////////////////////
610/// Get area in internal units of the facet with a given index.
611/// Possible index values:
612/// - 0 - all facets together
613/// - 1 to 6 - facet index from bottom to top Z
614
616{
617 Double_t area = 0.;
618 switch (index) {
619 case 0: area = 8. * (fDX * fDY + fDX * fDZ + fDY * fDZ); return area;
620 case 1:
621 case 6: area = 4. * fDX * fDY; return area;
622 case 2:
623 case 4: area = 4. * fDX * fDZ; return area;
624 case 3:
625 case 5: area = 4. * fDY * fDZ; return area;
626 }
627 return area;
628}
629
630////////////////////////////////////////////////////////////////////////////////
631/// Fills array with n random points located on the surface of indexed facet.
632/// The output array must be provided with a length of minimum 3*npoints. Returns
633/// true if operation succeeded.
634/// Possible index values:
635/// - 0 - all facets together
636/// - 1 to 6 - facet index from bottom to top Z
637
639{
640 if (index < 0 || index > 6)
641 return kFALSE;
642 Double_t surf[6];
643 Double_t area = 0.;
644 if (index == 0) {
645 for (Int_t isurf = 0; isurf < 6; isurf++) {
646 surf[isurf] = TGeoBBox::GetFacetArea(isurf + 1);
647 if (isurf > 0)
648 surf[isurf] += surf[isurf - 1];
649 }
650 area = surf[5];
651 }
652
653 for (Int_t i = 0; i < npoints; i++) {
654 // Generate randomly a surface index if needed.
655 Double_t *point = &array[3 * i];
656 Int_t surfindex = index;
657 if (surfindex == 0) {
658 Double_t val = area * gRandom->Rndm();
659 surfindex = 2 + TMath::BinarySearch(6, surf, val);
660 if (surfindex > 6)
661 surfindex = 6;
662 }
663 switch (surfindex) {
664 case 1:
665 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
666 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
667 point[2] = -fDZ;
668 break;
669 case 2:
670 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
671 point[1] = -fDY;
672 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
673 break;
674 case 3:
675 point[0] = -fDX;
676 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
677 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
678 break;
679 case 4:
680 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
681 point[1] = fDY;
682 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
683 break;
684 case 5:
685 point[0] = fDX;
686 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
687 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
688 break;
689 case 6:
690 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
691 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
692 point[2] = fDZ;
693 break;
694 }
695 }
696 return kTRUE;
697}
698
699////////////////////////////////////////////////////////////////////////////////
700/// Fills array with n random points located on the line segments of the shape mesh.
701/// The output array must be provided with a length of minimum 3*npoints. Returns
702/// true if operation is implemented.
703
705{
706 if (npoints < GetNmeshVertices()) {
707 Error("GetPointsOnSegments", "You should require at least %d points", GetNmeshVertices());
708 return kFALSE;
709 }
711 Int_t npnts = buff.NbPnts();
712 Int_t nsegs = buff.NbSegs();
713 // Copy buffered points in the array
714 memcpy(array, buff.fPnts, 3 * npnts * sizeof(Double_t));
715 Int_t ipoints = npoints - npnts;
716 Int_t icrt = 3 * npnts;
717 Int_t nperseg = (Int_t)(Double_t(ipoints) / nsegs);
718 Double_t *p0, *p1;
719 Double_t x, y, z, dx, dy, dz;
720 for (Int_t i = 0; i < nsegs; i++) {
721 p0 = &array[3 * buff.fSegs[3 * i + 1]];
722 p1 = &array[3 * buff.fSegs[3 * i + 2]];
723 if (i == (nsegs - 1))
724 nperseg = ipoints;
725 dx = (p1[0] - p0[0]) / (nperseg + 1);
726 dy = (p1[1] - p0[1]) / (nperseg + 1);
727 dz = (p1[2] - p0[2]) / (nperseg + 1);
728 for (Int_t j = 0; j < nperseg; j++) {
729 x = p0[0] + (j + 1) * dx;
730 y = p0[1] + (j + 1) * dy;
731 z = p0[2] + (j + 1) * dz;
732 array[icrt++] = x;
733 array[icrt++] = y;
734 array[icrt++] = z;
735 ipoints--;
736 }
737 }
738 return kTRUE;
739}
740
741////////////////////////////////////////////////////////////////////////////////
742/// Fills real parameters of a positioned box inside this one. Returns 0 if successful.
743
744Int_t TGeoBBox::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
745{
746 dx = dy = dz = 0;
747 if (mat->IsRotation()) {
748 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
749 return 1; // ### rotation not accepted ###
750 }
751 //--> translate the origin of the parametrized box to the frame of this box.
752 Double_t origin[3];
753 mat->LocalToMaster(parambox->GetOrigin(), origin);
754 if (!Contains(origin)) {
755 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
756 return 1; // ### wrong matrix ###
757 }
758 //--> now we have to get the valid range for all parametrized axis
759 Double_t xlo = 0, xhi = 0;
760 Double_t dd[3];
761 dd[0] = parambox->GetDX();
762 dd[1] = parambox->GetDY();
763 dd[2] = parambox->GetDZ();
764 for (Int_t iaxis = 0; iaxis < 3; iaxis++) {
765 if (dd[iaxis] >= 0)
766 continue;
767 TGeoBBox::GetAxisRange(iaxis + 1, xlo, xhi);
768 //-> compute best fitting parameter
769 dd[iaxis] = TMath::Min(origin[iaxis] - xlo, xhi - origin[iaxis]);
770 if (dd[iaxis] < 0) {
771 Error("GetFittingBox", "wrong matrix");
772 return 1;
773 }
774 }
775 dx = dd[0];
776 dy = dd[1];
777 dz = dd[2];
778 return 0;
779}
780
781////////////////////////////////////////////////////////////////////////////////
782/// In case shape has some negative parameters, these has to be computed
783/// in order to fit the mother
784
786{
788 return nullptr;
789 Double_t dx, dy, dz;
790 Int_t ierr = mother->GetFittingBox(this, mat, dx, dy, dz);
791 if (ierr) {
792 Error("GetMakeRuntimeShape", "cannot fit this to mother");
793 return nullptr;
794 }
795 return (new TGeoBBox(dx, dy, dz));
796}
797
798////////////////////////////////////////////////////////////////////////////////
799/// Returns numbers of vertices, segments and polygons composing the shape mesh.
800
801void TGeoBBox::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
802{
803 nvert = 8;
804 nsegs = 12;
805 npols = 6;
806}
807
808////////////////////////////////////////////////////////////////////////////////
809/// Prints shape parameters
810
811void TGeoBBox::InspectShape() const
812{
813 printf("*** Shape %s: TGeoBBox ***\n", GetName());
814 printf(" dX = %11.5f\n", fDX);
815 printf(" dY = %11.5f\n", fDY);
816 printf(" dZ = %11.5f\n", fDZ);
817 printf(" origin: x=%11.5f y=%11.5f z=%11.5f\n", fOrigin[0], fOrigin[1], fOrigin[2]);
818}
819
820////////////////////////////////////////////////////////////////////////////////
821/// Creates a TBuffer3D describing *this* shape.
822/// Coordinates are in local reference frame.
823
825{
826 TBuffer3D *buff = new TBuffer3D(TBuffer3DTypes::kGeneric, 8, 24, 12, 36, 6, 36);
827 if (buff) {
828 SetPoints(buff->fPnts);
829 SetSegsAndPols(*buff);
830 }
831
832 return buff;
833}
834
835////////////////////////////////////////////////////////////////////////////////
836/// Fills TBuffer3D structure for segments and polygons.
837
838void TGeoBBox::SetSegsAndPols(TBuffer3D &buff) const
839{
841
842 buff.fSegs[0] = c;
843 buff.fSegs[1] = 0;
844 buff.fSegs[2] = 1;
845 buff.fSegs[3] = c + 1;
846 buff.fSegs[4] = 1;
847 buff.fSegs[5] = 2;
848 buff.fSegs[6] = c + 1;
849 buff.fSegs[7] = 2;
850 buff.fSegs[8] = 3;
851 buff.fSegs[9] = c;
852 buff.fSegs[10] = 3;
853 buff.fSegs[11] = 0;
854 buff.fSegs[12] = c + 2;
855 buff.fSegs[13] = 4;
856 buff.fSegs[14] = 5;
857 buff.fSegs[15] = c + 2;
858 buff.fSegs[16] = 5;
859 buff.fSegs[17] = 6;
860 buff.fSegs[18] = c + 3;
861 buff.fSegs[19] = 6;
862 buff.fSegs[20] = 7;
863 buff.fSegs[21] = c + 3;
864 buff.fSegs[22] = 7;
865 buff.fSegs[23] = 4;
866 buff.fSegs[24] = c;
867 buff.fSegs[25] = 0;
868 buff.fSegs[26] = 4;
869 buff.fSegs[27] = c + 2;
870 buff.fSegs[28] = 1;
871 buff.fSegs[29] = 5;
872 buff.fSegs[30] = c + 1;
873 buff.fSegs[31] = 2;
874 buff.fSegs[32] = 6;
875 buff.fSegs[33] = c + 3;
876 buff.fSegs[34] = 3;
877 buff.fSegs[35] = 7;
878
879 buff.fPols[0] = c;
880 buff.fPols[1] = 4;
881 buff.fPols[2] = 0;
882 buff.fPols[3] = 9;
883 buff.fPols[4] = 4;
884 buff.fPols[5] = 8;
885 buff.fPols[6] = c + 1;
886 buff.fPols[7] = 4;
887 buff.fPols[8] = 1;
888 buff.fPols[9] = 10;
889 buff.fPols[10] = 5;
890 buff.fPols[11] = 9;
891 buff.fPols[12] = c;
892 buff.fPols[13] = 4;
893 buff.fPols[14] = 2;
894 buff.fPols[15] = 11;
895 buff.fPols[16] = 6;
896 buff.fPols[17] = 10;
897 buff.fPols[18] = c + 1;
898 buff.fPols[19] = 4;
899 buff.fPols[20] = 3;
900 buff.fPols[21] = 8;
901 buff.fPols[22] = 7;
902 buff.fPols[23] = 11;
903 buff.fPols[24] = c + 2;
904 buff.fPols[25] = 4;
905 buff.fPols[26] = 0;
906 buff.fPols[27] = 3;
907 buff.fPols[28] = 2;
908 buff.fPols[29] = 1;
909 buff.fPols[30] = c + 3;
910 buff.fPols[31] = 4;
911 buff.fPols[32] = 4;
912 buff.fPols[33] = 5;
913 buff.fPols[34] = 6;
914 buff.fPols[35] = 7;
915}
916
917////////////////////////////////////////////////////////////////////////////////
918/// Computes the closest distance from given point to this shape.
919
920Double_t TGeoBBox::Safety(const Double_t *point, Bool_t in) const
921{
922 Double_t safe, safy, safz;
923 if (in) {
924 safe = fDX - TMath::Abs(point[0] - fOrigin[0]);
925 safy = fDY - TMath::Abs(point[1] - fOrigin[1]);
926 safz = fDZ - TMath::Abs(point[2] - fOrigin[2]);
927 if (safy < safe)
928 safe = safy;
929 if (safz < safe)
930 safe = safz;
931 } else {
932 safe = -fDX + TMath::Abs(point[0] - fOrigin[0]);
933 safy = -fDY + TMath::Abs(point[1] - fOrigin[1]);
934 safz = -fDZ + TMath::Abs(point[2] - fOrigin[2]);
935 if (safy > safe)
936 safe = safy;
937 if (safz > safe)
938 safe = safz;
939 }
940 return safe;
941}
942
943////////////////////////////////////////////////////////////////////////////////
944/// Save a primitive as a C++ statement(s) on output stream "out".
945
946void TGeoBBox::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
947{
949 return;
950 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
951 out << " dx = " << fDX << ";" << std::endl;
952 out << " dy = " << fDY << ";" << std::endl;
953 out << " dz = " << fDZ << ";" << std::endl;
956 out << " origin[0] = " << fOrigin[0] << ";" << std::endl;
957 out << " origin[1] = " << fOrigin[1] << ";" << std::endl;
958 out << " origin[2] = " << fOrigin[2] << ";" << std::endl;
959 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz,origin);"
960 << std::endl;
961 } else {
962 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz);" << std::endl;
963 }
965}
966
967////////////////////////////////////////////////////////////////////////////////
968/// Set parameters of the box.
969
971{
972 fDX = dx;
973 fDY = dy;
974 fDZ = dz;
975 if (origin) {
976 fOrigin[0] = origin[0];
977 fOrigin[1] = origin[1];
978 fOrigin[2] = origin[2];
979 }
982 return;
983 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
985}
986
987////////////////////////////////////////////////////////////////////////////////
988/// Set dimensions based on the array of parameters
989/// param[0] - half-length in x
990/// param[1] - half-length in y
991/// param[2] - half-length in z
992
994{
995 if (!param) {
996 Error("SetDimensions", "null parameters");
997 return;
998 }
999 fDX = param[0];
1000 fDY = param[1];
1001 fDZ = param[2];
1004 return;
1005 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
1007}
1008
1009////////////////////////////////////////////////////////////////////////////////
1010/// Fill box vertices to an array.
1011
1013{
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Fill box points.
1019
1021{
1022 if (!points)
1023 return;
1024 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1025 xmin = -fDX + fOrigin[0];
1026 xmax = fDX + fOrigin[0];
1027 ymin = -fDY + fOrigin[1];
1028 ymax = fDY + fOrigin[1];
1029 zmin = -fDZ + fOrigin[2];
1030 zmax = fDZ + fOrigin[2];
1031 points[0] = xmin;
1032 points[1] = ymin;
1033 points[2] = zmin;
1034 points[3] = xmin;
1035 points[4] = ymax;
1036 points[5] = zmin;
1037 points[6] = xmax;
1038 points[7] = ymax;
1039 points[8] = zmin;
1040 points[9] = xmax;
1041 points[10] = ymin;
1042 points[11] = zmin;
1043 points[12] = xmin;
1044 points[13] = ymin;
1045 points[14] = zmax;
1046 points[15] = xmin;
1047 points[16] = ymax;
1048 points[17] = zmax;
1049 points[18] = xmax;
1050 points[19] = ymax;
1051 points[20] = zmax;
1052 points[21] = xmax;
1053 points[22] = ymin;
1054 points[23] = zmax;
1055}
1056
1057////////////////////////////////////////////////////////////////////////////////
1058/// Fill box points.
1059
1061{
1062 if (!points)
1063 return;
1064 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1065 xmin = -fDX + fOrigin[0];
1066 xmax = fDX + fOrigin[0];
1067 ymin = -fDY + fOrigin[1];
1068 ymax = fDY + fOrigin[1];
1069 zmin = -fDZ + fOrigin[2];
1070 zmax = fDZ + fOrigin[2];
1071 points[0] = xmin;
1072 points[1] = ymin;
1073 points[2] = zmin;
1074 points[3] = xmin;
1075 points[4] = ymax;
1076 points[5] = zmin;
1077 points[6] = xmax;
1078 points[7] = ymax;
1079 points[8] = zmin;
1080 points[9] = xmax;
1081 points[10] = ymin;
1082 points[11] = zmin;
1083 points[12] = xmin;
1084 points[13] = ymin;
1085 points[14] = zmax;
1086 points[15] = xmin;
1087 points[16] = ymax;
1088 points[17] = zmax;
1089 points[18] = xmax;
1090 points[19] = ymax;
1091 points[20] = zmax;
1092 points[21] = xmax;
1093 points[22] = ymin;
1094 points[23] = zmax;
1095}
1096
1097////////////////////////////////////////////////////////////////////////////////
1098////// fill size of this 3-D object
1099//// TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
1100//// if (painter) painter->AddSize3D(8, 12, 6);
1101
1102void TGeoBBox::Sizeof3D() const {}
1103
1104////////////////////////////////////////////////////////////////////////////////
1105/// Fills a static 3D buffer and returns a reference.
1106
1107const TBuffer3D &TGeoBBox::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1108{
1109 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1110
1111 FillBuffer3D(buffer, reqSections, localFrame);
1112
1113 // TODO: A box itself has has nothing more as already described
1114 // by bounding box. How will viewer interpret?
1115 if (reqSections & TBuffer3D::kRawSizes) {
1116 if (buffer.SetRawSizes(8, 3 * 8, 12, 3 * 12, 6, 6 * 6)) {
1117 buffer.SetSectionsValid(TBuffer3D::kRawSizes);
1118 }
1119 }
1120 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1121 SetPoints(buffer.fPnts);
1122 if (!buffer.fLocalFrame) {
1123 TransformPoints(buffer.fPnts, buffer.NbPnts());
1124 }
1125
1126 SetSegsAndPols(buffer);
1127 buffer.SetSectionsValid(TBuffer3D::kRaw);
1128 }
1129
1130 return buffer;
1131}
1132
1133////////////////////////////////////////////////////////////////////////////////
1134/// Fills the supplied buffer, with sections in desired frame
1135/// See TBuffer3D.h for explanation of sections, frame etc.
1136
1137void TGeoBBox::FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
1138{
1139 TGeoShape::FillBuffer3D(buffer, reqSections, localFrame);
1140
1141 if (reqSections & TBuffer3D::kBoundingBox) {
1142 Double_t halfLengths[3] = {fDX, fDY, fDZ};
1143 buffer.SetAABoundingBox(fOrigin, halfLengths);
1144
1145 if (!buffer.fLocalFrame) {
1146 TransformPoints(buffer.fBBVertex[0], 8);
1147 }
1149 }
1150}
1151
1152////////////////////////////////////////////////////////////////////////////////
1153/// Check the inside status for each of the points in the array.
1154/// Input: Array of point coordinates + vector size
1155/// Output: Array of Booleans for the inside of each point
1156
1157void TGeoBBox::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1158{
1159 for (Int_t i = 0; i < vecsize; i++)
1160 inside[i] = Contains(&points[3 * i]);
1161}
1162
1163////////////////////////////////////////////////////////////////////////////////
1164/// Compute the normal for an array o points so that norm.dot.dir is positive
1165/// Input: Arrays of point coordinates and directions + vector size
1166/// Output: Array of normal directions
1167
1168void TGeoBBox::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1169{
1170 for (Int_t i = 0; i < vecsize; i++)
1171 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
1172}
1173
1174////////////////////////////////////////////////////////////////////////////////
1175/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1176
1177void TGeoBBox::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1178 Double_t *step) const
1179{
1180 for (Int_t i = 0; i < vecsize; i++)
1181 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1182}
1183
1184////////////////////////////////////////////////////////////////////////////////
1185/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1186
1187void TGeoBBox::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1188 Double_t *step) const
1189{
1190 for (Int_t i = 0; i < vecsize; i++)
1191 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1192}
1193
1194////////////////////////////////////////////////////////////////////////////////
1195/// Compute safe distance from each of the points in the input array.
1196/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1197/// Output: Safety values
1198
1199void TGeoBBox::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1200{
1201 for (Int_t i = 0; i < vecsize; i++)
1202 safe[i] = Safety(&points[3 * i], inside[i]);
1203}
#define c(i)
Definition RSha256.hxx:101
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
float xmin
float ymin
float xmax
float ymax
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
UInt_t NbPnts() const
Definition TBuffer3D.h:80
UInt_t NbSegs() const
Definition TBuffer3D.h:81
@ kBoundingBox
Definition TBuffer3D.h:51
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:114
Bool_t fLocalFrame
Definition TBuffer3D.h:90
void SetAABoundingBox(const Double_t origin[3], const Double_t halfLengths[3])
Set fBBVertex in kBoundingBox section to a axis aligned (local) BB using supplied origin and box half...
Double_t * fPnts
Definition TBuffer3D.h:113
Double_t fBBVertex[8][3]
Definition TBuffer3D.h:108
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fill the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections,...
virtual Double_t GetFacetArea(Int_t index=0) const
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
virtual const Double_t * GetOrigin() const
Definition TGeoBBox.h:79
Double_t DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Double_t fDX
Definition TGeoBBox.h:20
const char * GetAxisName(Int_t iaxis) const override
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
void GetBoundingCylinder(Double_t *param) const override
void SetPoints(Double_t *points) const override
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=nullptr)
Double_t DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1, Double_t step=TGeoShape::Big(), Double_t *safe=nullptr) const override
Double_t Capacity() const override
virtual Double_t GetDX() const
Definition TGeoBBox.h:76
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Computes distance from point (px,py) to the object.
Bool_t CouldBeCrossed(const Double_t *point, const Double_t *dir) const override
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
void SetBoxPoints(Double_t *points) const
Bool_t GetPointsOnSegments(Int_t npoints, Double_t *array) const override
virtual Double_t GetDZ() const
Definition TGeoBBox.h:78
virtual Double_t GetDY() const
Definition TGeoBBox.h:77
Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const override
Double_t fOrigin[3]
Definition TGeoBBox.h:23
Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const override
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) override
TBuffer3D * MakeBuffer3D() const override
void InspectShape() const override
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Bool_t Contains(const Double_t *point) const override
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const override
void Sizeof3D() const override
TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const override
Double_t fDY
Definition TGeoBBox.h:21
void SetSegsAndPols(TBuffer3D &buffer) const override
virtual Bool_t GetPointsOnFacet(Int_t index, Int_t npoints, Double_t *array) const
void SetDimensions(Double_t *param) override
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Stub implementation to avoid forcing implementation at this stage.
Int_t GetNmeshVertices() const override
Definition TGeoBBox.h:75
Double_t fDZ
Definition TGeoBBox.h:22
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
~TGeoBBox() override
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
void ComputeBBox() override
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition TGeoMatrix.h:38
Bool_t IsRotation() const
Definition TGeoMatrix.h:65
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
Node containing an offset.
Definition TGeoNode.h:184
Base finder class for patterns.
void SetDivIndex(Int_t index)
Base abstract class for all shapes.
Definition TGeoShape.h:25
static Double_t Big()
Definition TGeoShape.h:87
Int_t GetBasicColor() const
Get the basic color (0-7).
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.
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,...
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
virtual Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const =0
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:64
@ kGeoRunTimeShape
Definition TGeoShape.h:40
static Double_t Tolerance()
Definition TGeoShape.h:90
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:167
Volume families.
Definition TGeoVolume.h:266
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:175
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:244
Int_t GetNdaughters() const
Definition TGeoVolume.h:362
TObjArray * GetNodes()
Definition TGeoVolume.h:169
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:798
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1005
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:559
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:986
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:347
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123