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