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/// Check if 2 positioned boxes overlap.
211
212Bool_t
213TGeoBBox::AreOverlapping(const TGeoBBox *box1, const TGeoMatrix *mat1, const TGeoBBox *box2, const TGeoMatrix *mat2)
214{
215 Double_t master[3];
216 Double_t local[3];
217 Double_t ldir1[3], ldir2[3];
218 const Double_t *o1 = box1->GetOrigin();
219 const Double_t *o2 = box2->GetOrigin();
220 // Convert center of first box to the local frame of second
221 mat1->LocalToMaster(o1, master);
222 mat2->MasterToLocal(master, local);
223 if (TGeoBBox::Contains(local, box2->GetDX(), box2->GetDY(), box2->GetDZ(), o2))
224 return kTRUE;
225 Double_t distsq = (local[0] - o2[0]) * (local[0] - o2[0]) + (local[1] - o2[1]) * (local[1] - o2[1]) +
226 (local[2] - o2[2]) * (local[2] - o2[2]);
227 // Compute distance between box centers and compare with max value
228 Double_t rmaxsq = (box1->GetDX() + box2->GetDX()) * (box1->GetDX() + box2->GetDX()) +
229 (box1->GetDY() + box2->GetDY()) * (box1->GetDY() + box2->GetDY()) +
230 (box1->GetDZ() + box2->GetDZ()) * (box1->GetDZ() + box2->GetDZ());
231 if (distsq > rmaxsq + TGeoShape::Tolerance())
232 return kFALSE;
233 // We are still not sure: shoot a ray from the center of "1" towards the
234 // center of 2.
235 Double_t dir[3];
236 mat1->LocalToMaster(o1, ldir1);
237 mat2->LocalToMaster(o2, ldir2);
238 distsq = 1. / TMath::Sqrt(distsq);
239 dir[0] = (ldir2[0] - ldir1[0]) * distsq;
240 dir[1] = (ldir2[1] - ldir1[1]) * distsq;
241 dir[2] = (ldir2[2] - ldir1[2]) * distsq;
242 mat1->MasterToLocalVect(dir, ldir1);
243 mat2->MasterToLocalVect(dir, ldir2);
244 // Distance to exit from o1
245 Double_t dist1 = TGeoBBox::DistFromInside(o1, ldir1, box1->GetDX(), box1->GetDY(), box1->GetDZ(), o1);
246 // Distance to enter from o2
247 Double_t dist2 = TGeoBBox::DistFromOutside(local, ldir2, box2->GetDX(), box2->GetDY(), box2->GetDZ(), o2);
248 if (dist1 > dist2)
249 return kTRUE;
250 return kFALSE;
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// Computes capacity of the shape in [length^3].
255
257{
258 return (8. * fDX * fDY * fDZ);
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// Computes normal to closest surface from POINT.
263
264void TGeoBBox::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
265{
266 memset(norm, 0, 3 * sizeof(Double_t));
267 Double_t saf[3];
268 Int_t i;
269 saf[0] = TMath::Abs(TMath::Abs(point[0] - fOrigin[0]) - fDX);
270 saf[1] = TMath::Abs(TMath::Abs(point[1] - fOrigin[1]) - fDY);
271 saf[2] = TMath::Abs(TMath::Abs(point[2] - fOrigin[2]) - fDZ);
272 i = TMath::LocMin(3, saf);
273 norm[i] = (dir[i] > 0) ? 1 : (-1);
274}
275
276////////////////////////////////////////////////////////////////////////////////
277/// Decides fast if the bounding box could be crossed by a vector.
278
279Bool_t TGeoBBox::CouldBeCrossed(const Double_t *point, const Double_t *dir) const
280{
281 Double_t mind = fDX;
282 if (fDY < mind)
283 mind = fDY;
284 if (fDZ < mind)
285 mind = fDZ;
286 Double_t dx = fOrigin[0] - point[0];
287 Double_t dy = fOrigin[1] - point[1];
288 Double_t dz = fOrigin[2] - point[2];
289 Double_t do2 = dx * dx + dy * dy + dz * dz;
290 if (do2 <= (mind * mind))
291 return kTRUE;
292 Double_t rmax2 = fDX * fDX + fDY * fDY + fDZ * fDZ;
293 if (do2 <= rmax2)
294 return kTRUE;
295 // inside bounding sphere
296 Double_t doct = dx * dir[0] + dy * dir[1] + dz * dir[2];
297 // leaving ray
298 if (doct <= 0)
299 return kFALSE;
300 Double_t dirnorm = dir[0] * dir[0] + dir[1] * dir[1] + dir[2] * dir[2];
301 if ((doct * doct) >= (do2 - rmax2) * dirnorm)
302 return kTRUE;
303 return kFALSE;
304}
305
306////////////////////////////////////////////////////////////////////////////////
307/// Compute closest distance from point px,py to each corner.
308
310{
311 const Int_t numPoints = 8;
312 return ShapeDistancetoPrimitive(numPoints, px, py);
313}
314
315////////////////////////////////////////////////////////////////////////////////
316/// Divide this box shape belonging to volume "voldiv" into ndiv equal volumes
317/// called divname, from start position with the given step. Returns pointer
318/// to created division cell volume. In case a wrong division axis is supplied,
319/// returns pointer to volume to be divided.
320
322TGeoBBox::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
323{
324 TGeoShape *shape; //--- shape to be created
325 TGeoVolume *vol; //--- division volume to be created
326 TGeoVolumeMulti *vmulti; //--- generic divided volume
327 TGeoPatternFinder *finder; //--- finder to be attached
328 TString opt = ""; //--- option to be attached
329 Double_t end = start + ndiv * step;
330 switch (iaxis) {
331 case 1: //--- divide on X
332 shape = new TGeoBBox(step / 2., fDY, fDZ);
333 finder = new TGeoPatternX(voldiv, ndiv, start, end);
334 opt = "X";
335 break;
336 case 2: //--- divide on Y
337 shape = new TGeoBBox(fDX, step / 2., fDZ);
338 finder = new TGeoPatternY(voldiv, ndiv, start, end);
339 opt = "Y";
340 break;
341 case 3: //--- divide on Z
342 shape = new TGeoBBox(fDX, fDY, step / 2.);
343 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
344 opt = "Z";
345 break;
346 default: Error("Divide", "Wrong axis type for division"); return nullptr;
347 }
348 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
349 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
350 vmulti->AddVolume(vol);
351 voldiv->SetFinder(finder);
352 finder->SetDivIndex(voldiv->GetNdaughters());
353 for (Int_t ic = 0; ic < ndiv; ic++) {
354 voldiv->AddNodeOffset(vol, ic, start + step / 2. + ic * step, opt.Data());
355 ((TGeoNodeOffset *)voldiv->GetNodes()->At(voldiv->GetNdaughters() - 1))->SetFinder(finder);
356 }
357 return vmulti;
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// Compute bounding box - nothing to do in this case.
362
363void TGeoBBox::ComputeBBox() {}
364
365////////////////////////////////////////////////////////////////////////////////
366/// Test if point is inside this shape.
367
368Bool_t TGeoBBox::Contains(const Double_t *point) const
369{
370 if (TMath::Abs(point[2] - fOrigin[2]) > fDZ)
371 return kFALSE;
372 if (TMath::Abs(point[0] - fOrigin[0]) > fDX)
373 return kFALSE;
374 if (TMath::Abs(point[1] - fOrigin[1]) > fDY)
375 return kFALSE;
376 return kTRUE;
377}
378
379////////////////////////////////////////////////////////////////////////////////
380/// Static method to check if point[3] is located inside a box of having dx, dy, dz
381/// as half-lengths.
382
383Bool_t TGeoBBox::Contains(const Double_t *point, Double_t dx, Double_t dy, Double_t dz, const Double_t *origin)
384{
385 if (TMath::Abs(point[2] - origin[2]) > dz)
386 return kFALSE;
387 if (TMath::Abs(point[0] - origin[0]) > dx)
388 return kFALSE;
389 if (TMath::Abs(point[1] - origin[1]) > dy)
390 return kFALSE;
391 return kTRUE;
392}
393
394////////////////////////////////////////////////////////////////////////////////
395/// Compute distance from inside point to surface of the box.
396/// Boundary safe algorithm.
397
399TGeoBBox::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
400{
401 Double_t s, smin, saf[6];
402 Double_t newpt[3];
403 Int_t i;
404 for (i = 0; i < 3; i++)
405 newpt[i] = point[i] - fOrigin[i];
406 saf[0] = fDX + newpt[0];
407 saf[1] = fDX - newpt[0];
408 saf[2] = fDY + newpt[1];
409 saf[3] = fDY - newpt[1];
410 saf[4] = fDZ + newpt[2];
411 saf[5] = fDZ - newpt[2];
412 if (iact < 3 && safe) {
413 smin = saf[0];
414 // compute safe distance
415 for (i = 1; i < 6; i++)
416 if (saf[i] < smin)
417 smin = saf[i];
418 *safe = smin;
419 if (smin < 0)
420 *safe = 0.0;
421 if (iact == 0)
422 return TGeoShape::Big();
423 if (iact == 1 && step < *safe)
424 return TGeoShape::Big();
425 }
426 // compute distance to surface
427 smin = TGeoShape::Big();
428 for (i = 0; i < 3; i++) {
429 if (dir[i] != 0) {
430 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
431 if (s < 0)
432 return 0.0;
433 if (s < smin)
434 smin = s;
435 }
436 }
437 return smin;
438}
439
440////////////////////////////////////////////////////////////////////////////////
441/// Compute distance from inside point to surface of the box.
442/// Boundary safe algorithm.
443
444Double_t TGeoBBox::DistFromInside(const Double_t *point, const Double_t *dir, Double_t dx, Double_t dy, Double_t dz,
445 const Double_t *origin, Double_t /*stepmax*/)
446{
447 Double_t s, smin, saf[6];
448 Double_t newpt[3];
449 Int_t i;
450 for (i = 0; i < 3; i++)
451 newpt[i] = point[i] - origin[i];
452 saf[0] = dx + newpt[0];
453 saf[1] = dx - newpt[0];
454 saf[2] = dy + newpt[1];
455 saf[3] = dy - newpt[1];
456 saf[4] = dz + newpt[2];
457 saf[5] = dz - newpt[2];
458 // compute distance to surface
459 smin = TGeoShape::Big();
460 for (i = 0; i < 3; i++) {
461 if (dir[i] != 0) {
462 s = (dir[i] > 0) ? (saf[(i << 1) + 1] / dir[i]) : (-saf[i << 1] / dir[i]);
463 if (s < 0)
464 return 0.0;
465 if (s < smin)
466 smin = s;
467 }
468 }
469 return smin;
470}
471
472////////////////////////////////////////////////////////////////////////////////
473/// Compute distance from outside point to surface of the box.
474/// Boundary safe algorithm.
475
477TGeoBBox::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
478{
479 Bool_t in = kTRUE;
480 Double_t saf[3];
481 Double_t par[3];
482 Double_t newpt[3];
483 Int_t i, j;
484 for (i = 0; i < 3; i++)
485 newpt[i] = point[i] - fOrigin[i];
486 par[0] = fDX;
487 par[1] = fDY;
488 par[2] = fDZ;
489 for (i = 0; i < 3; i++) {
490 saf[i] = TMath::Abs(newpt[i]) - par[i];
491 if (saf[i] >= step)
492 return TGeoShape::Big();
493 if (in && saf[i] > 0)
494 in = kFALSE;
495 }
496 if (iact < 3 && safe) {
497 // compute safe distance
498 if (in) {
499 *safe = 0.0;
500 } else {
501 *safe = saf[0];
502 if (saf[1] > *safe)
503 *safe = saf[1];
504 if (saf[2] > *safe)
505 *safe = saf[2];
506 }
507 if (iact == 0)
508 return TGeoShape::Big();
509 if (iact == 1 && step < *safe)
510 return TGeoShape::Big();
511 }
512 // compute distance from point to box
513 // protection in case point is actually inside box
514 if (in) {
515 j = 0;
516 Double_t ss = saf[0];
517 if (saf[1] > ss) {
518 ss = saf[1];
519 j = 1;
520 }
521 if (saf[2] > ss)
522 j = 2;
523 if (newpt[j] * dir[j] > 0)
524 return TGeoShape::Big(); // in fact exiting
525 return 0.0;
526 }
527 for (i = 0; i < 3; i++) {
528 if (saf[i] < 0)
529 continue;
530 if (newpt[i] * dir[i] >= 0)
531 continue;
532 Double_t snxt = saf[i] / TMath::Abs(dir[i]);
533 Int_t ibreak = 0;
534 for (j = 0; j < 3; j++) {
535 if (j == i)
536 continue;
537 Double_t coord = newpt[j] + snxt * dir[j];
538 if (TMath::Abs(coord) > par[j]) {
539 ibreak = 1;
540 break;
541 }
542 }
543 if (!ibreak)
544 return snxt;
545 }
546 return TGeoShape::Big();
547}
548
549////////////////////////////////////////////////////////////////////////////////
550/// Compute distance from outside point to surface of the box.
551/// Boundary safe algorithm.
552
554 const Double_t *origin, Double_t stepmax)
555{
556 Bool_t in = kTRUE;
557 Double_t saf[3];
558 Double_t par[3];
559 Double_t newpt[3];
560 Int_t i, j;
561 for (i = 0; i < 3; i++)
562 newpt[i] = point[i] - origin[i];
563 par[0] = dx;
564 par[1] = dy;
565 par[2] = dz;
566 for (i = 0; i < 3; i++) {
567 saf[i] = TMath::Abs(newpt[i]) - par[i];
568 if (saf[i] >= stepmax)
569 return TGeoShape::Big();
570 if (in && saf[i] > 0)
571 in = kFALSE;
572 }
573 // In case point is inside return ZERO
574 if (in)
575 return 0.0;
576 Double_t coord, snxt = TGeoShape::Big();
577 Int_t ibreak = 0;
578 for (i = 0; i < 3; i++) {
579 if (saf[i] < 0)
580 continue;
581 if (newpt[i] * dir[i] >= 0)
582 continue;
583 snxt = saf[i] / TMath::Abs(dir[i]);
584 ibreak = 0;
585 for (j = 0; j < 3; j++) {
586 if (j == i)
587 continue;
588 coord = newpt[j] + snxt * dir[j];
589 if (TMath::Abs(coord) > par[j]) {
590 ibreak = 1;
591 break;
592 }
593 }
594 if (!ibreak)
595 return snxt;
596 }
597 return TGeoShape::Big();
598}
599
600////////////////////////////////////////////////////////////////////////////////
601/// Returns name of axis IAXIS.
602
603const char *TGeoBBox::GetAxisName(Int_t iaxis) const
604{
605 switch (iaxis) {
606 case 1: return "X";
607 case 2: return "Y";
608 case 3: return "Z";
609 default: return "UNDEFINED";
610 }
611}
612
613////////////////////////////////////////////////////////////////////////////////
614/// Get range of shape for a given axis.
615
617{
618 xlo = 0;
619 xhi = 0;
620 Double_t dx = 0;
621 switch (iaxis) {
622 case 1:
623 xlo = fOrigin[0] - fDX;
624 xhi = fOrigin[0] + fDX;
625 dx = 2 * fDX;
626 return dx;
627 case 2:
628 xlo = fOrigin[1] - fDY;
629 xhi = fOrigin[1] + fDY;
630 dx = 2 * fDY;
631 return dx;
632 case 3:
633 xlo = fOrigin[2] - fDZ;
634 xhi = fOrigin[2] + fDZ;
635 dx = 2 * fDZ;
636 return dx;
637 }
638 return dx;
639}
640
641////////////////////////////////////////////////////////////////////////////////
642/// Fill vector param[4] with the bounding cylinder parameters. The order
643/// is the following : Rmin, Rmax, Phi1, Phi2
644
646{
647 param[0] = 0.; // Rmin
648 param[1] = fDX * fDX + fDY * fDY; // Rmax
649 param[2] = 0.; // Phi1
650 param[3] = 360.; // Phi2
651}
652
653////////////////////////////////////////////////////////////////////////////////
654/// Get area in internal units of the facet with a given index.
655/// Possible index values:
656/// - 0 - all facets together
657/// - 1 to 6 - facet index from bottom to top Z
658
660{
661 Double_t area = 0.;
662 switch (index) {
663 case 0: area = 8. * (fDX * fDY + fDX * fDZ + fDY * fDZ); return area;
664 case 1:
665 case 6: area = 4. * fDX * fDY; return area;
666 case 2:
667 case 4: area = 4. * fDX * fDZ; return area;
668 case 3:
669 case 5: area = 4. * fDY * fDZ; return area;
670 }
671 return area;
672}
673
674////////////////////////////////////////////////////////////////////////////////
675/// Fills array with n random points located on the surface of indexed facet.
676/// The output array must be provided with a length of minimum 3*npoints. Returns
677/// true if operation succeeded.
678/// Possible index values:
679/// - 0 - all facets together
680/// - 1 to 6 - facet index from bottom to top Z
681
683{
684 if (index < 0 || index > 6)
685 return kFALSE;
686 Double_t surf[6];
687 Double_t area = 0.;
688 if (index == 0) {
689 for (Int_t isurf = 0; isurf < 6; isurf++) {
690 surf[isurf] = TGeoBBox::GetFacetArea(isurf + 1);
691 if (isurf > 0)
692 surf[isurf] += surf[isurf - 1];
693 }
694 area = surf[5];
695 }
696
697 for (Int_t i = 0; i < npoints; i++) {
698 // Generate randomly a surface index if needed.
699 Double_t *point = &array[3 * i];
700 Int_t surfindex = index;
701 if (surfindex == 0) {
702 Double_t val = area * gRandom->Rndm();
703 surfindex = 2 + TMath::BinarySearch(6, surf, val);
704 if (surfindex > 6)
705 surfindex = 6;
706 }
707 switch (surfindex) {
708 case 1:
709 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
710 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
711 point[2] = -fDZ;
712 break;
713 case 2:
714 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
715 point[1] = -fDY;
716 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
717 break;
718 case 3:
719 point[0] = -fDX;
720 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
721 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
722 break;
723 case 4:
724 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
725 point[1] = fDY;
726 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
727 break;
728 case 5:
729 point[0] = fDX;
730 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
731 point[2] = -fDZ + 2 * fDZ * gRandom->Rndm();
732 break;
733 case 6:
734 point[0] = -fDX + 2 * fDX * gRandom->Rndm();
735 point[1] = -fDY + 2 * fDY * gRandom->Rndm();
736 point[2] = fDZ;
737 break;
738 }
739 }
740 return kTRUE;
741}
742
743////////////////////////////////////////////////////////////////////////////////
744/// Fills array with n random points located on the line segments of the shape mesh.
745/// The output array must be provided with a length of minimum 3*npoints. Returns
746/// true if operation is implemented.
747
749{
750 if (npoints < GetNmeshVertices()) {
751 Error("GetPointsOnSegments", "You should require at least %d points", GetNmeshVertices());
752 return kFALSE;
753 }
755 Int_t npnts = buff.NbPnts();
756 Int_t nsegs = buff.NbSegs();
757 // Copy buffered points in the array
758 memcpy(array, buff.fPnts, 3 * npnts * sizeof(Double_t));
759 Int_t ipoints = npoints - npnts;
760 Int_t icrt = 3 * npnts;
761 Int_t nperseg = (Int_t)(Double_t(ipoints) / nsegs);
762 Double_t *p0, *p1;
763 Double_t x, y, z, dx, dy, dz;
764 for (Int_t i = 0; i < nsegs; i++) {
765 p0 = &array[3 * buff.fSegs[3 * i + 1]];
766 p1 = &array[3 * buff.fSegs[3 * i + 2]];
767 if (i == (nsegs - 1))
768 nperseg = ipoints;
769 dx = (p1[0] - p0[0]) / (nperseg + 1);
770 dy = (p1[1] - p0[1]) / (nperseg + 1);
771 dz = (p1[2] - p0[2]) / (nperseg + 1);
772 for (Int_t j = 0; j < nperseg; j++) {
773 x = p0[0] + (j + 1) * dx;
774 y = p0[1] + (j + 1) * dy;
775 z = p0[2] + (j + 1) * dz;
776 array[icrt++] = x;
777 array[icrt++] = y;
778 array[icrt++] = z;
779 ipoints--;
780 }
781 }
782 return kTRUE;
783}
784
785////////////////////////////////////////////////////////////////////////////////
786/// Fills real parameters of a positioned box inside this one. Returns 0 if successful.
787
788Int_t TGeoBBox::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
789{
790 dx = dy = dz = 0;
791 if (mat->IsRotation()) {
792 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
793 return 1; // ### rotation not accepted ###
794 }
795 //--> translate the origin of the parametrized box to the frame of this box.
796 Double_t origin[3];
797 mat->LocalToMaster(parambox->GetOrigin(), origin);
798 if (!Contains(origin)) {
799 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
800 return 1; // ### wrong matrix ###
801 }
802 //--> now we have to get the valid range for all parametrized axis
803 Double_t xlo = 0, xhi = 0;
804 Double_t dd[3];
805 dd[0] = parambox->GetDX();
806 dd[1] = parambox->GetDY();
807 dd[2] = parambox->GetDZ();
808 for (Int_t iaxis = 0; iaxis < 3; iaxis++) {
809 if (dd[iaxis] >= 0)
810 continue;
811 TGeoBBox::GetAxisRange(iaxis + 1, xlo, xhi);
812 //-> compute best fitting parameter
813 dd[iaxis] = TMath::Min(origin[iaxis] - xlo, xhi - origin[iaxis]);
814 if (dd[iaxis] < 0) {
815 Error("GetFittingBox", "wrong matrix");
816 return 1;
817 }
818 }
819 dx = dd[0];
820 dy = dd[1];
821 dz = dd[2];
822 return 0;
823}
824
825////////////////////////////////////////////////////////////////////////////////
826/// In case shape has some negative parameters, these has to be computed
827/// in order to fit the mother
828
830{
832 return nullptr;
833 Double_t dx, dy, dz;
834 Int_t ierr = mother->GetFittingBox(this, mat, dx, dy, dz);
835 if (ierr) {
836 Error("GetMakeRuntimeShape", "cannot fit this to mother");
837 return nullptr;
838 }
839 return (new TGeoBBox(dx, dy, dz));
840}
841
842////////////////////////////////////////////////////////////////////////////////
843/// Returns numbers of vertices, segments and polygons composing the shape mesh.
844
845void TGeoBBox::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
846{
847 nvert = 8;
848 nsegs = 12;
849 npols = 6;
850}
851
852////////////////////////////////////////////////////////////////////////////////
853/// Prints shape parameters
854
855void TGeoBBox::InspectShape() const
856{
857 printf("*** Shape %s: TGeoBBox ***\n", GetName());
858 printf(" dX = %11.5f\n", fDX);
859 printf(" dY = %11.5f\n", fDY);
860 printf(" dZ = %11.5f\n", fDZ);
861 printf(" origin: x=%11.5f y=%11.5f z=%11.5f\n", fOrigin[0], fOrigin[1], fOrigin[2]);
862}
863
864////////////////////////////////////////////////////////////////////////////////
865/// Creates a TBuffer3D describing *this* shape.
866/// Coordinates are in local reference frame.
867
869{
870 TBuffer3D *buff = new TBuffer3D(TBuffer3DTypes::kGeneric, 8, 24, 12, 36, 6, 36);
871 if (buff) {
872 SetPoints(buff->fPnts);
873 SetSegsAndPols(*buff);
874 }
875
876 return buff;
877}
878
879////////////////////////////////////////////////////////////////////////////////
880/// Fills TBuffer3D structure for segments and polygons.
881
882void TGeoBBox::SetSegsAndPols(TBuffer3D &buff) const
883{
885
886 buff.fSegs[0] = c;
887 buff.fSegs[1] = 0;
888 buff.fSegs[2] = 1;
889 buff.fSegs[3] = c + 1;
890 buff.fSegs[4] = 1;
891 buff.fSegs[5] = 2;
892 buff.fSegs[6] = c + 1;
893 buff.fSegs[7] = 2;
894 buff.fSegs[8] = 3;
895 buff.fSegs[9] = c;
896 buff.fSegs[10] = 3;
897 buff.fSegs[11] = 0;
898 buff.fSegs[12] = c + 2;
899 buff.fSegs[13] = 4;
900 buff.fSegs[14] = 5;
901 buff.fSegs[15] = c + 2;
902 buff.fSegs[16] = 5;
903 buff.fSegs[17] = 6;
904 buff.fSegs[18] = c + 3;
905 buff.fSegs[19] = 6;
906 buff.fSegs[20] = 7;
907 buff.fSegs[21] = c + 3;
908 buff.fSegs[22] = 7;
909 buff.fSegs[23] = 4;
910 buff.fSegs[24] = c;
911 buff.fSegs[25] = 0;
912 buff.fSegs[26] = 4;
913 buff.fSegs[27] = c + 2;
914 buff.fSegs[28] = 1;
915 buff.fSegs[29] = 5;
916 buff.fSegs[30] = c + 1;
917 buff.fSegs[31] = 2;
918 buff.fSegs[32] = 6;
919 buff.fSegs[33] = c + 3;
920 buff.fSegs[34] = 3;
921 buff.fSegs[35] = 7;
922
923 buff.fPols[0] = c;
924 buff.fPols[1] = 4;
925 buff.fPols[2] = 0;
926 buff.fPols[3] = 9;
927 buff.fPols[4] = 4;
928 buff.fPols[5] = 8;
929 buff.fPols[6] = c + 1;
930 buff.fPols[7] = 4;
931 buff.fPols[8] = 1;
932 buff.fPols[9] = 10;
933 buff.fPols[10] = 5;
934 buff.fPols[11] = 9;
935 buff.fPols[12] = c;
936 buff.fPols[13] = 4;
937 buff.fPols[14] = 2;
938 buff.fPols[15] = 11;
939 buff.fPols[16] = 6;
940 buff.fPols[17] = 10;
941 buff.fPols[18] = c + 1;
942 buff.fPols[19] = 4;
943 buff.fPols[20] = 3;
944 buff.fPols[21] = 8;
945 buff.fPols[22] = 7;
946 buff.fPols[23] = 11;
947 buff.fPols[24] = c + 2;
948 buff.fPols[25] = 4;
949 buff.fPols[26] = 0;
950 buff.fPols[27] = 3;
951 buff.fPols[28] = 2;
952 buff.fPols[29] = 1;
953 buff.fPols[30] = c + 3;
954 buff.fPols[31] = 4;
955 buff.fPols[32] = 4;
956 buff.fPols[33] = 5;
957 buff.fPols[34] = 6;
958 buff.fPols[35] = 7;
959}
960
961////////////////////////////////////////////////////////////////////////////////
962/// Computes the closest distance from given point to this shape.
963
964Double_t TGeoBBox::Safety(const Double_t *point, Bool_t in) const
965{
966 Double_t safe, safy, safz;
967 if (in) {
968 safe = fDX - TMath::Abs(point[0] - fOrigin[0]);
969 safy = fDY - TMath::Abs(point[1] - fOrigin[1]);
970 safz = fDZ - TMath::Abs(point[2] - fOrigin[2]);
971 if (safy < safe)
972 safe = safy;
973 if (safz < safe)
974 safe = safz;
975 } else {
976 safe = -fDX + TMath::Abs(point[0] - fOrigin[0]);
977 safy = -fDY + TMath::Abs(point[1] - fOrigin[1]);
978 safz = -fDZ + TMath::Abs(point[2] - fOrigin[2]);
979 if (safy > safe)
980 safe = safy;
981 if (safz > safe)
982 safe = safz;
983 }
984 return safe;
985}
986
987////////////////////////////////////////////////////////////////////////////////
988/// Save a primitive as a C++ statement(s) on output stream "out".
989
990void TGeoBBox::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
991{
993 return;
994 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
995 out << " dx = " << fDX << ";" << std::endl;
996 out << " dy = " << fDY << ";" << std::endl;
997 out << " dz = " << fDZ << ";" << std::endl;
1000 out << " origin[0] = " << fOrigin[0] << ";" << std::endl;
1001 out << " origin[1] = " << fOrigin[1] << ";" << std::endl;
1002 out << " origin[2] = " << fOrigin[2] << ";" << std::endl;
1003 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz,origin);"
1004 << std::endl;
1005 } else {
1006 out << " TGeoShape *" << GetPointerName() << " = new TGeoBBox(\"" << GetName() << "\", dx,dy,dz);" << std::endl;
1007 }
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Set parameters of the box.
1013
1015{
1016 fDX = dx;
1017 fDY = dy;
1018 fDZ = dz;
1019 if (origin) {
1020 fOrigin[0] = origin[0];
1021 fOrigin[1] = origin[1];
1022 fOrigin[2] = origin[2];
1023 }
1026 return;
1027 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
1029}
1030
1031////////////////////////////////////////////////////////////////////////////////
1032/// Set dimensions based on the array of parameters
1033/// param[0] - half-length in x
1034/// param[1] - half-length in y
1035/// param[2] - half-length in z
1036
1038{
1039 if (!param) {
1040 Error("SetDimensions", "null parameters");
1041 return;
1042 }
1043 fDX = param[0];
1044 fDY = param[1];
1045 fDZ = param[2];
1048 return;
1049 if ((fDX < 0) || (fDY < 0) || (fDZ < 0))
1051}
1052
1053////////////////////////////////////////////////////////////////////////////////
1054/// Fill box vertices to an array.
1055
1057{
1059}
1060
1061////////////////////////////////////////////////////////////////////////////////
1062/// Fill box points.
1063
1065{
1066 if (!points)
1067 return;
1068 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1069 xmin = -fDX + fOrigin[0];
1070 xmax = fDX + fOrigin[0];
1071 ymin = -fDY + fOrigin[1];
1072 ymax = fDY + fOrigin[1];
1073 zmin = -fDZ + fOrigin[2];
1074 zmax = fDZ + fOrigin[2];
1075 points[0] = xmin;
1076 points[1] = ymin;
1077 points[2] = zmin;
1078 points[3] = xmin;
1079 points[4] = ymax;
1080 points[5] = zmin;
1081 points[6] = xmax;
1082 points[7] = ymax;
1083 points[8] = zmin;
1084 points[9] = xmax;
1085 points[10] = ymin;
1086 points[11] = zmin;
1087 points[12] = xmin;
1088 points[13] = ymin;
1089 points[14] = zmax;
1090 points[15] = xmin;
1091 points[16] = ymax;
1092 points[17] = zmax;
1093 points[18] = xmax;
1094 points[19] = ymax;
1095 points[20] = zmax;
1096 points[21] = xmax;
1097 points[22] = ymin;
1098 points[23] = zmax;
1099}
1100
1101////////////////////////////////////////////////////////////////////////////////
1102/// Fill box points.
1103
1105{
1106 if (!points)
1107 return;
1108 Double_t xmin, xmax, ymin, ymax, zmin, zmax;
1109 xmin = -fDX + fOrigin[0];
1110 xmax = fDX + fOrigin[0];
1111 ymin = -fDY + fOrigin[1];
1112 ymax = fDY + fOrigin[1];
1113 zmin = -fDZ + fOrigin[2];
1114 zmax = fDZ + fOrigin[2];
1115 points[0] = xmin;
1116 points[1] = ymin;
1117 points[2] = zmin;
1118 points[3] = xmin;
1119 points[4] = ymax;
1120 points[5] = zmin;
1121 points[6] = xmax;
1122 points[7] = ymax;
1123 points[8] = zmin;
1124 points[9] = xmax;
1125 points[10] = ymin;
1126 points[11] = zmin;
1127 points[12] = xmin;
1128 points[13] = ymin;
1129 points[14] = zmax;
1130 points[15] = xmin;
1131 points[16] = ymax;
1132 points[17] = zmax;
1133 points[18] = xmax;
1134 points[19] = ymax;
1135 points[20] = zmax;
1136 points[21] = xmax;
1137 points[22] = ymin;
1138 points[23] = zmax;
1139}
1140
1141////////////////////////////////////////////////////////////////////////////////
1142////// fill size of this 3-D object
1143//// TVirtualGeoPainter *painter = gGeoManager->GetGeomPainter();
1144//// if (painter) painter->AddSize3D(8, 12, 6);
1145
1146void TGeoBBox::Sizeof3D() const {}
1147
1148////////////////////////////////////////////////////////////////////////////////
1149/// Fills a static 3D buffer and returns a reference.
1150
1151const TBuffer3D &TGeoBBox::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
1152{
1153 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
1154
1155 FillBuffer3D(buffer, reqSections, localFrame);
1156
1157 // TODO: A box itself has has nothing more as already described
1158 // by bounding box. How will viewer interpret?
1159 if (reqSections & TBuffer3D::kRawSizes) {
1160 if (buffer.SetRawSizes(8, 3 * 8, 12, 3 * 12, 6, 6 * 6)) {
1161 buffer.SetSectionsValid(TBuffer3D::kRawSizes);
1162 }
1163 }
1164 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
1165 SetPoints(buffer.fPnts);
1166 if (!buffer.fLocalFrame) {
1167 TransformPoints(buffer.fPnts, buffer.NbPnts());
1168 }
1169
1170 SetSegsAndPols(buffer);
1171 buffer.SetSectionsValid(TBuffer3D::kRaw);
1172 }
1173
1174 return buffer;
1175}
1176
1177////////////////////////////////////////////////////////////////////////////////
1178/// Fills the supplied buffer, with sections in desired frame
1179/// See TBuffer3D.h for explanation of sections, frame etc.
1180
1181void TGeoBBox::FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const
1182{
1183 TGeoShape::FillBuffer3D(buffer, reqSections, localFrame);
1184
1185 if (reqSections & TBuffer3D::kBoundingBox) {
1186 Double_t halfLengths[3] = {fDX, fDY, fDZ};
1187 buffer.SetAABoundingBox(fOrigin, halfLengths);
1188
1189 if (!buffer.fLocalFrame) {
1190 TransformPoints(buffer.fBBVertex[0], 8);
1191 }
1193 }
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Check the inside status for each of the points in the array.
1198/// Input: Array of point coordinates + vector size
1199/// Output: Array of Booleans for the inside of each point
1200
1201void TGeoBBox::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
1202{
1203 for (Int_t i = 0; i < vecsize; i++)
1204 inside[i] = Contains(&points[3 * i]);
1205}
1206
1207////////////////////////////////////////////////////////////////////////////////
1208/// Compute the normal for an array o points so that norm.dot.dir is positive
1209/// Input: Arrays of point coordinates and directions + vector size
1210/// Output: Array of normal directions
1211
1212void TGeoBBox::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
1213{
1214 for (Int_t i = 0; i < vecsize; i++)
1215 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
1216}
1217
1218////////////////////////////////////////////////////////////////////////////////
1219/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1220
1221void TGeoBBox::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1222 Double_t *step) const
1223{
1224 for (Int_t i = 0; i < vecsize; i++)
1225 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1226}
1227
1228////////////////////////////////////////////////////////////////////////////////
1229/// Compute distance from array of input points having directions specified by dirs. Store output in dists
1230
1231void TGeoBBox::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
1232 Double_t *step) const
1233{
1234 for (Int_t i = 0; i < vecsize; i++)
1235 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
1236}
1237
1238////////////////////////////////////////////////////////////////////////////////
1239/// Compute safe distance from each of the points in the input array.
1240/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
1241/// Output: Safety values
1242
1243void TGeoBBox::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
1244{
1245 for (Int_t i = 0; i < vecsize; i++)
1246 safe[i] = Safety(&points[3 * i], inside[i]);
1247}
#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:377
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:82
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:79
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:81
virtual Double_t GetDY() const
Definition TGeoBBox.h:80
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.
static Bool_t AreOverlapping(const TGeoBBox *box1, const TGeoMatrix *mat1, const TGeoBBox *box2, const TGeoMatrix *mat2) R__DEPRECATED(6
Int_t GetNmeshVertices() const override
Definition TGeoBBox.h:78
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
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void MasterToLocalVect(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
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:207
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
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:982
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
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