Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoTrd1.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 24/10/01
3// TGeoTrd1::Contains() and DistFromInside() 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 TGeoTrd1
14\ingroup Trapezoids
15
16A trapezoid with only X varying with Z. It is defined by the
17half-length in Z, the half-length in X at the lowest and highest Z
18planes and the half-length in Y:
19
20~~~ {.cpp}
21TGeoTrd1(Double_t dx1,Double_t dx2,Double_t dy,Double_t dz);
22~~~
23
24Begin_Macro
25{
26 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
27 new TGeoManager("trd1", "poza8");
28 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
29 TGeoMedium *med = new TGeoMedium("MED",1,mat);
30 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
31 gGeoManager->SetTopVolume(top);
32 TGeoVolume *vol = gGeoManager->MakeTrd1("Trd1",med, 10,20,30,40);
33 vol->SetLineWidth(2);
34 top->AddNode(vol,1);
35 gGeoManager->CloseGeometry();
36 gGeoManager->SetNsegments(80);
37 top->Draw();
38 TView *view = gPad->GetView();
39 view->ShowAxis();
40}
41End_Macro
42*/
43
44#include <iostream>
45
46#include "TGeoManager.h"
47#include "TGeoMatrix.h"
48#include "TGeoVolume.h"
49#include "TGeoTrd1.h"
50#include "TMath.h"
51
53
54////////////////////////////////////////////////////////////////////////////////
55/// dummy ctor
56
58{
59 fDz = fDx1 = fDx2 = fDy = 0;
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// constructor.
65
67 :TGeoBBox(0,0,0)
68{
70 fDx1 = dx1;
71 fDx2 = dx2;
72 fDy = dy;
73 fDz = dz;
74 if ((dx1<0) || (dx2<0) || (dy<0) || (dz<0)) {
76 printf("trd1 : dx1=%f, dx2=%f, dy=%f, dz=%f\n",
77 dx1,dx2,dy,dz);
78 }
79 else ComputeBBox();
80}
81
82////////////////////////////////////////////////////////////////////////////////
83/// constructor.
84
86 :TGeoBBox(name, 0,0,0)
87{
89 fDx1 = dx1;
90 fDx2 = dx2;
91 fDy = dy;
92 fDz = dz;
93 if ((dx1<0) || (dx2<0) || (dy<0) || (dz<0)) {
95 printf("trd1 : dx1=%f, dx2=%f, dy=%f, dz=%f\n",
96 dx1,dx2,dy,dz);
97 }
98 else ComputeBBox();
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// ctor with an array of parameters
103/// - param[0] = dx1
104/// - param[1] = dx2
105/// - param[2] = dy
106/// - param[3] = dz
107
109 :TGeoBBox(0,0,0)
110{
112 SetDimensions(param);
113 if ((fDx1<0) || (fDx2<0) || (fDy<=0) || (fDz<=0)) SetShapeBit(kGeoRunTimeShape);
114 else ComputeBBox();
115}
116
117////////////////////////////////////////////////////////////////////////////////
118/// destructor
119
121{
122}
123
124////////////////////////////////////////////////////////////////////////////////
125/// Computes capacity of the shape in [length^3]
126
128{
129 Double_t capacity = 4.*(fDx1+fDx2)*fDy*fDz;
130 return capacity;
131}
132
133////////////////////////////////////////////////////////////////////////////////
134/// compute bounding box for a trd1
135
137{
139 fDY = fDy;
140 fDZ = fDz;
141 memset(fOrigin, 0, 3*sizeof(Double_t));
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// Compute normal to closest surface from POINT.
146
147void TGeoTrd1::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
148{
149 Double_t safe, safemin;
150 //--- Compute safety first
151 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
152 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
153 // check Z facettes
154 safe = safemin = TMath::Abs(fDz-TMath::Abs(point[2]));
155 norm[0] = norm[1] = 0;
156 norm[2] = (dir[2]>=0)?1:-1;
157 if (safe<1E-6) return;
158 // check X facettes
159 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
160 if (distx>=0) {
161 safe=TMath::Abs(distx-TMath::Abs(point[0]))*calf;
162 if (safe<safemin) {
163 safemin = safe;
164 norm[0] = (point[0]>0)?calf:(-calf);
165 norm[1] = 0;
166 norm[2] = calf*fx;
167 Double_t dot = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
168 if (dot<0) {
169 norm[0] = -norm[0];
170 norm[2] = -norm[2];
171 }
172 if (safe<1E-6) return;
173 }
174 }
175 // check Y facettes
176 safe = TMath::Abs(fDy-TMath::Abs(point[1]));
177 if (safe<safemin) {
178 norm[0] = norm[2] = 0;
179 norm[1] = (dir[1]>=0)?1:-1;
180 }
181}
182
183////////////////////////////////////////////////////////////////////////////////
184/// test if point is inside this shape
185/// check Z range
186
188{
189 if (TMath::Abs(point[2]) > fDz) return kFALSE;
190 // then y
191 if (TMath::Abs(point[1]) > fDy) return kFALSE;
192 // then x
193 Double_t dx = 0.5*(fDx2*(point[2]+fDz)+fDx1*(fDz-point[2]))/fDz;
194 if (TMath::Abs(point[0]) > dx) return kFALSE;
195 return kTRUE;
196}
197
198////////////////////////////////////////////////////////////////////////////////
199/// Compute distance from inside point to surface of the trd1
200/// Boundary safe algorithm.
201
202Double_t TGeoTrd1::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
203{
204 if (iact<3 && safe) {
205 // compute safe distance
206 *safe = Safety(point, kTRUE);
207 if (iact==0) return TGeoShape::Big();
208 if (iact==1 && step<*safe) return TGeoShape::Big();
209 }
210
211 //--- Compute safety first
212 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
213 Double_t cn;
214 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
215 //--- Compute distance to this shape
216 // first check if Z facettes are crossed
217 Double_t dist[3];
218 for (Int_t i=0; i<3; i++) dist[i]=TGeoShape::Big();
219 if (dir[2]<0) {
220 dist[0]=-(point[2]+fDz)/dir[2];
221 } else if (dir[2]>0) {
222 dist[0]=(fDz-point[2])/dir[2];
223 }
224 if (dist[0]<=0) return 0.0;
225 // now check X facettes
226 cn = -dir[0]+fx*dir[2];
227 if (cn>0) {
228 dist[1] = point[0]+distx;
229 if (dist[1]<=0) return 0.0;
230 dist[1] /= cn;
231 }
232 cn = dir[0]+fx*dir[2];
233 if (cn>0) {
234 Double_t s = distx-point[0];
235 if (s<=0) return 0.0;
236 s /= cn;
237 if (s<dist[1]) dist[1] = s;
238 }
239 // now check Y facettes
240 if (dir[1]<0) {
241 dist[2]=-(point[1]+fDy)/dir[1];
242 } else if (dir[1]>0) {
243 dist[2]=(fDy-point[1])/dir[1];
244 }
245 if (dist[2]<=0) return 0.0;
246 return dist[TMath::LocMin(3,dist)];
247}
248
249////////////////////////////////////////////////////////////////////////////////
250/// get the most visible corner from outside point and the normals
251
252void TGeoTrd1::GetVisibleCorner(const Double_t *point, Double_t *vertex, Double_t *normals) const
253{
254 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
255 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
256 Double_t salf = calf*fx;
257 // check visibility of X faces
258 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
259 memset(normals, 0, 9*sizeof(Double_t));
260 TGeoTrd1 *trd1 = (TGeoTrd1*)this;
261 if (point[0]>distx) {
262 // hi x face visible
263 trd1->SetShapeBit(kGeoVisX);
264 normals[0]=calf;
265 normals[2]=salf;
266 } else {
268 normals[0]=-calf;
269 normals[2]=salf;
270 }
271 if (point[1]>fDy) {
272 // hi y face visible
273 trd1->SetShapeBit(kGeoVisY);
274 normals[4]=1;
275 } else {
277 normals[4]=-1;
278 }
279 if (point[2]>fDz) {
280 // hi z face visible
281 trd1->SetShapeBit(kGeoVisZ);
282 normals[8]=1;
283 } else {
285 normals[8]=-1;
286 }
288}
289
290////////////////////////////////////////////////////////////////////////////////
291/// get the opposite corner of the intersected face
292
293void TGeoTrd1::GetOppositeCorner(const Double_t * /*point*/, Int_t inorm, Double_t *vertex, Double_t *normals) const
294{
295 TGeoTrd1 *trd1 = (TGeoTrd1*)this;
296 if (inorm != 0) {
297 // change x face
299 normals[0]=-normals[0];
300 }
301 if (inorm != 1) {
302 // change y face
304 normals[4]=-normals[4];
305 }
306 if (inorm != 2) {
307 // hi z face visible
309 normals[8]=-normals[8];
310 }
312}
313
314////////////////////////////////////////////////////////////////////////////////
315/// Compute distance from outside point to surface of the trd1
316/// Boundary safe algorithm
317
318Double_t TGeoTrd1::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
319{
320 if (iact<3 && safe) {
321 // compute safe distance
322 *safe = Safety(point, kFALSE);
323 if (iact==0) return TGeoShape::Big();
324 if (iact==1 && step<*safe) return TGeoShape::Big();
325 }
326 // find a visible face
327 Double_t xnew,ynew,znew;
328 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
329 Double_t cn;
330 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
331 Bool_t in = kTRUE;
332 Double_t safx = distx-TMath::Abs(point[0]);
333 Double_t safy = fDy-TMath::Abs(point[1]);
334 Double_t safz = fDz-TMath::Abs(point[2]);
335
336 //--- Compute distance to this shape
337 // first check if Z facettes are crossed
338 if (point[2]<=-fDz) {
339 if (dir[2]<=0) return TGeoShape::Big();
340 in = kFALSE;
341 Double_t snxt = -(fDz+point[2])/dir[2];
342 // find extrapolated X and Y
343 xnew = point[0]+snxt*dir[0];
344 if (TMath::Abs(xnew) <= fDx1) {
345 ynew = point[1]+snxt*dir[1];
346 if (TMath::Abs(ynew) <= fDy) return snxt;
347 }
348 } else if (point[2]>=fDz) {
349 if (dir[2]>=0) return TGeoShape::Big();
350 in = kFALSE;
351 Double_t snxt = (fDz-point[2])/dir[2];
352 // find extrapolated X and Y
353 xnew = point[0]+snxt*dir[0];
354 if (TMath::Abs(xnew) <= fDx2) {
355 ynew = point[1]+snxt*dir[1];
356 if (TMath::Abs(ynew) <= fDy) return snxt;
357 }
358 }
359 // check if X facettes are crossed
360 if (point[0]<=-distx) {
361 cn = -dir[0]+fx*dir[2];
362 if (cn>=0) return TGeoShape::Big();
363 in = kFALSE;
364 Double_t snxt = (point[0]+distx)/cn;
365 // find extrapolated Y and Z
366 ynew = point[1]+snxt*dir[1];
367 if (TMath::Abs(ynew) <= fDy) {
368 znew = point[2]+snxt*dir[2];
369 if (TMath::Abs(znew) <= fDz) return snxt;
370 }
371 }
372 if (point[0]>=distx) {
373 cn = dir[0]+fx*dir[2];
374 if (cn>=0) return TGeoShape::Big();
375 in = kFALSE;
376 Double_t snxt = (distx-point[0])/cn;
377 // find extrapolated Y and Z
378 ynew = point[1]+snxt*dir[1];
379 if (TMath::Abs(ynew) < fDy) {
380 znew = point[2]+snxt*dir[2];
381 if (TMath::Abs(znew) < fDz) return snxt;
382 }
383 }
384 // finally check Y facettes
385 if (point[1]<=-fDy) {
386 cn = -dir[1];
387 if (cn>=0) return TGeoShape::Big();
388 in = kFALSE;
389 Double_t snxt = (point[1]+fDy)/cn;
390 // find extrapolated X and Z
391 znew = point[2]+snxt*dir[2];
392 if (TMath::Abs(znew) < fDz) {
393 xnew = point[0]+snxt*dir[0];
394 Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
395 if (TMath::Abs(xnew) < dx) return snxt;
396 }
397 } else if (point[1]>=fDy) {
398 cn = dir[1];
399 if (cn>=0) return TGeoShape::Big();
400 in = kFALSE;
401 Double_t snxt = (fDy-point[1])/cn;
402 // find extrapolated X and Z
403 znew = point[2]+snxt*dir[2];
404 if (TMath::Abs(znew) < fDz) {
405 xnew = point[0]+snxt*dir[0];
406 Double_t dx = 0.5*(fDx1+fDx2)-fx*znew;
407 if (TMath::Abs(xnew) < dx) return snxt;
408 }
409 }
410 if (!in) return TGeoShape::Big();
411 // Point actually inside
412 if (safz<safx && safz<safy) {
413 if (point[2]*dir[2]>=0) return TGeoShape::Big();
414 return 0.0;
415 }
416 if (safy<safx) {
417 if (point[1]*dir[1]>=0) return TGeoShape::Big();
418 return 0.0;
419 }
420 cn = TMath::Sign(1.0,point[0])*dir[0]+fx*dir[2];
421 if (cn>=0) return TGeoShape::Big();
422 return 0.0;
423}
424
425////////////////////////////////////////////////////////////////////////////////
426/// Divide this trd1 shape belonging to volume "voldiv" into ndiv volumes
427/// called divname, from start position with the given step. Returns pointer
428/// to created division cell volume in case of Y divisions. For Z divisions just
429/// return the pointer to the volume to be divided. In case a wrong
430/// division axis is supplied, returns pointer to volume that was divided.
431
432TGeoVolume *TGeoTrd1::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
433 Double_t start, Double_t step)
434{
435 TGeoShape *shape; //--- shape to be created
436 TGeoVolume *vol; //--- division volume to be created
437 TGeoVolumeMulti *vmulti; //--- generic divided volume
438 TGeoPatternFinder *finder; //--- finder to be attached
439 TString opt = ""; //--- option to be attached
440 Double_t zmin, zmax, dx1n, dx2n;
441 Int_t id;
442 Double_t end = start+ndiv*step;
443 switch (iaxis) {
444 case 1:
445 Warning("Divide", "dividing a Trd1 on X not implemented");
446 return 0;
447 case 2:
448 finder = new TGeoPatternY(voldiv, ndiv, start, end);
449 voldiv->SetFinder(finder);
450 finder->SetDivIndex(voldiv->GetNdaughters());
451 shape = new TGeoTrd1(fDx1, fDx2, step/2, fDz);
452 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
453 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
454 vmulti->AddVolume(vol);
455 opt = "Y";
456 for (id=0; id<ndiv; id++) {
457 voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
458 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
459 }
460 return vmulti;
461 case 3:
462 finder = new TGeoPatternZ(voldiv, ndiv, start, end);
463 voldiv->SetFinder(finder);
464 finder->SetDivIndex(voldiv->GetNdaughters());
465 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
466 for (id=0; id<ndiv; id++) {
467 zmin = start+id*step;
468 zmax = start+(id+1)*step;
469 dx1n = 0.5*(fDx1*(fDz-zmin)+fDx2*(fDz+zmin))/fDz;
470 dx2n = 0.5*(fDx1*(fDz-zmax)+fDx2*(fDz+zmax))/fDz;
471 shape = new TGeoTrd1(dx1n, dx2n, fDy, step/2.);
472 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
473 vmulti->AddVolume(vol);
474 opt = "Z";
475 voldiv->AddNodeOffset(vol, id, start+step/2+id*step, opt.Data());
476 ((TGeoNodeOffset*)voldiv->GetNodes()->At(voldiv->GetNdaughters()-1))->SetFinder(finder);
477 }
478 return vmulti;
479 default:
480 Error("Divide", "Wrong axis type for division");
481 return 0;
482 }
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Get range of shape for a given axis.
487
489{
490 xlo = 0;
491 xhi = 0;
492 Double_t dx = 0;
493 switch (iaxis) {
494 case 2:
495 xlo = -fDy;
496 xhi = fDy;
497 dx = xhi-xlo;
498 return dx;
499 case 3:
500 xlo = -fDz;
501 xhi = fDz;
502 dx = xhi-xlo;
503 return dx;
504 }
505 return dx;
506}
507
508////////////////////////////////////////////////////////////////////////////////
509/// Fill vector param[4] with the bounding cylinder parameters. The order
510/// is the following : Rmin, Rmax, Phi1, Phi2
511
513{
515}
516
517////////////////////////////////////////////////////////////////////////////////
518/// Fills real parameters of a positioned box inside this. Returns 0 if successful.
519
520Int_t TGeoTrd1::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
521{
522 dx=dy=dz=0;
523 if (mat->IsRotation()) {
524 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
525 return 1; // ### rotation not accepted ###
526 }
527 //--> translate the origin of the parametrized box to the frame of this box.
528 Double_t origin[3];
529 mat->LocalToMaster(parambox->GetOrigin(), origin);
530 if (!Contains(origin)) {
531 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
532 return 1; // ### wrong matrix ###
533 }
534 //--> now we have to get the valid range for all parametrized axis
535 Double_t dd[3];
536 dd[0] = parambox->GetDX();
537 dd[1] = parambox->GetDY();
538 dd[2] = parambox->GetDZ();
539 //-> check if Z range is fixed
540 if (dd[2]<0) {
541 dd[2] = TMath::Min(origin[2]+fDz, fDz-origin[2]);
542 if (dd[2]<0) {
543 Error("GetFittingBox", "wrong matrix");
544 return 1;
545 }
546 }
547 //-> check if Y range is fixed
548 if (dd[1]<0) {
549 dd[1] = TMath::Min(origin[1]+fDy, fDy-origin[1]);
550 if (dd[1]<0) {
551 Error("GetFittingBox", "wrong matrix");
552 return 1;
553 }
554 }
555 if (dd[0]>=0) {
556 dx = dd[0];
557 dy = dd[1];
558 dz = dd[2];
559 return 0;
560 }
561 //-> check now range at Z = origin[2] +/- dd[2]
562 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
563 Double_t dx0 = 0.5*(fDx1+fDx2);
564 Double_t z=origin[2]-dd[2];
565 dd[0] = dx0-fx*z-origin[0];
566 z=origin[2]+dd[2];
567 dd[0] = TMath::Min(dd[0], dx0-fx*z-origin[0]);
568 if (dd[0]<0) {
569 Error("GetFittingBox", "wrong matrix");
570 return 1;
571 }
572 dx = dd[0];
573 dy = dd[1];
574 dz = dd[2];
575 return 0;
576}
577
578////////////////////////////////////////////////////////////////////////////////
579/// in case shape has some negative parameters, these has to be computed
580/// in order to fit the mother
581
583{
584 if (!TestShapeBit(kGeoRunTimeShape)) return 0;
585 if (!mother->TestShapeBit(kGeoTrd1)) {
586 Error("GetMakeRuntimeShape", "invalid mother");
587 return 0;
588 }
589 Double_t dx1, dx2, dy, dz;
590 if (fDx1<0) dx1=((TGeoTrd1*)mother)->GetDx1();
591 else dx1=fDx1;
592 if (fDx2<0) dx2=((TGeoTrd1*)mother)->GetDx2();
593 else dx2=fDx2;
594 if (fDy<0) dy=((TGeoTrd1*)mother)->GetDy();
595 else dy=fDy;
596 if (fDz<0) dz=((TGeoTrd1*)mother)->GetDz();
597 else dz=fDz;
598
599 return (new TGeoTrd1(dx1, dx2, dy, dz));
600}
601
602////////////////////////////////////////////////////////////////////////////////
603/// print shape parameters
604
606{
607 printf("*** Shape %s: TGeoTrd1 ***\n", GetName());
608 printf(" dx1 = %11.5f\n", fDx1);
609 printf(" dx2 = %11.5f\n", fDx2);
610 printf(" dy = %11.5f\n", fDy);
611 printf(" dz = %11.5f\n", fDz);
612 printf(" Bounding box:\n");
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// computes the closest distance from given point to this shape, according
618/// to option. The matching point on the shape is stored in spoint.
619
621{
622 Double_t saf[3];
623 //--- Compute safety first
624 // check Z facettes
625 saf[0] = fDz-TMath::Abs(point[2]);
626 Double_t fx = 0.5*(fDx1-fDx2)/fDz;
627 Double_t calf = 1./TMath::Sqrt(1.0+fx*fx);
628 // check X facettes
629 Double_t distx = 0.5*(fDx1+fDx2)-fx*point[2];
630 if (distx<0) saf[1]=TGeoShape::Big();
631 else saf[1]=(distx-TMath::Abs(point[0]))*calf;
632 // check Y facettes
633 saf[2] = fDy-TMath::Abs(point[1]);
634 if (in) return saf[TMath::LocMin(3,saf)];
635 for (Int_t i=0; i<3; i++) saf[i]=-saf[i];
636 return saf[TMath::LocMax(3,saf)];
637}
638
639////////////////////////////////////////////////////////////////////////////////
640/// Save a primitive as a C++ statement(s) on output stream "out".
641
642void TGeoTrd1::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
643{
645 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
646 out << " dx1 = " << fDx1 << ";" << std::endl;
647 out << " dx2 = " << fDx2 << ";" << std::endl;
648 out << " dy = " << fDy << ";" << std::endl;
649 out << " dz = " << fDZ << ";" << std::endl;
650 out << " TGeoShape *" << GetPointerName() << " = new TGeoTrd1(\"" << GetName() << "\", dx1,dx2,dy,dz);" << std::endl;
652}
653
654////////////////////////////////////////////////////////////////////////////////
655/// set trd1 params in one step :
656
658{
659 fDx1 = param[0];
660 fDx2 = param[1];
661 fDy = param[2];
662 fDz = param[3];
663 ComputeBBox();
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// set vertex of a corner according to visibility flags
668
670{
671 if (TestShapeBit(kGeoVisX)) {
672 if (TestShapeBit(kGeoVisZ)) {
673 vertex[0] = fDx2;
674 vertex[2] = fDz;
676 } else {
677 vertex[0] = fDx1;
678 vertex[2] = -fDz;
680 }
681 } else {
682 if (TestShapeBit(kGeoVisZ)) {
683 vertex[0] = -fDx2;
684 vertex[2] = fDz;
686 } else {
687 vertex[0] = -fDx1;
688 vertex[2] = -fDz;
690 }
691 }
692}
693
694////////////////////////////////////////////////////////////////////////////////
695/// create arb8 mesh points
696
698{
699 if (!points) return;
700 points[ 0] = -fDx1; points[ 1] = -fDy; points[ 2] = -fDz;
701 points[ 3] = -fDx1; points[ 4] = fDy; points[ 5] = -fDz;
702 points[ 6] = fDx1; points[ 7] = fDy; points[ 8] = -fDz;
703 points[ 9] = fDx1; points[10] = -fDy; points[11] = -fDz;
704 points[12] = -fDx2; points[13] = -fDy; points[14] = fDz;
705 points[15] = -fDx2; points[16] = fDy; points[17] = fDz;
706 points[18] = fDx2; points[19] = fDy; points[20] = fDz;
707 points[21] = fDx2; points[22] = -fDy; points[23] = fDz;
708}
709
710////////////////////////////////////////////////////////////////////////////////
711/// create arb8 mesh points
712
714{
715 if (!points) return;
716 points[ 0] = -fDx1; points[ 1] = -fDy; points[ 2] = -fDz;
717 points[ 3] = -fDx1; points[ 4] = fDy; points[ 5] = -fDz;
718 points[ 6] = fDx1; points[ 7] = fDy; points[ 8] = -fDz;
719 points[ 9] = fDx1; points[10] = -fDy; points[11] = -fDz;
720 points[12] = -fDx2; points[13] = -fDy; points[14] = fDz;
721 points[15] = -fDx2; points[16] = fDy; points[17] = fDz;
722 points[18] = fDx2; points[19] = fDy; points[20] = fDz;
723 points[21] = fDx2; points[22] = -fDy; points[23] = fDz;
724}
725
726////////////////////////////////////////////////////////////////////////////////
727/// fill size of this 3-D object
728
730{
732}
733
734////////////////////////////////////////////////////////////////////////////////
735/// Check the inside status for each of the points in the array.
736/// Input: Array of point coordinates + vector size
737/// Output: Array of Booleans for the inside of each point
738
739void TGeoTrd1::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
740{
741 for (Int_t i=0; i<vecsize; i++) inside[i] = Contains(&points[3*i]);
742}
743
744////////////////////////////////////////////////////////////////////////////////
745/// Compute the normal for an array o points so that norm.dot.dir is positive
746/// Input: Arrays of point coordinates and directions + vector size
747/// Output: Array of normal directions
748
749void TGeoTrd1::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
750{
751 for (Int_t i=0; i<vecsize; i++) ComputeNormal(&points[3*i], &dirs[3*i], &norms[3*i]);
752}
753
754////////////////////////////////////////////////////////////////////////////////
755/// Compute distance from array of input points having directions specified by dirs. Store output in dists
756
757void TGeoTrd1::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
758{
759 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromInside(&points[3*i], &dirs[3*i], 3, step[i]);
760}
761
762////////////////////////////////////////////////////////////////////////////////
763/// Compute distance from array of input points having directions specified by dirs. Store output in dists
764
765void TGeoTrd1::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t* step) const
766{
767 for (Int_t i=0; i<vecsize; i++) dists[i] = DistFromOutside(&points[3*i], &dirs[3*i], 3, step[i]);
768}
769
770////////////////////////////////////////////////////////////////////////////////
771/// Compute safe distance from each of the points in the input array.
772/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
773/// Output: Safety values
774
775void TGeoTrd1::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
776{
777 for (Int_t i=0; i<vecsize; i++) safe[i] = Safety(&points[3*i], inside[i]);
778}
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
XFontStruct * id
Definition TGX11.cxx:109
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
point * points
Definition X3DBuffer.c:22
Box class.
Definition TGeoBBox.h:18
virtual const Double_t * GetOrigin() const
Definition TGeoBBox.h:77
Double_t fDX
Definition TGeoBBox.h:21
virtual void InspectShape() const
Prints shape parameters.
Definition TGeoBBox.cxx:819
virtual Double_t GetDX() const
Definition TGeoBBox.h:74
virtual Double_t GetDZ() const
Definition TGeoBBox.h:76
virtual void Sizeof3D() const
virtual Double_t GetDY() const
Definition TGeoBBox.h:75
Double_t fOrigin[3]
Definition TGeoBBox.h:24
Double_t fDY
Definition TGeoBBox.h:22
Double_t fDZ
Definition TGeoBBox.h:23
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoBBox.cxx:609
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
Geometrical transformation package.
Definition TGeoMatrix.h:41
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
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
const char * GetPointerName() const
Provide a pointer name containing uid.
virtual const char * GetName() const
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:65
@ kGeoRunTimeShape
Definition TGeoShape.h:41
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:163
A trapezoid with only X varying with Z.
Definition TGeoTrd1.h:18
virtual void Sizeof3D() const
fill size of this 3-D object
Definition TGeoTrd1.cxx:729
virtual void ComputeBBox()
compute bounding box for a trd1
Definition TGeoTrd1.cxx:136
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...
Definition TGeoTrd1.cxx:765
virtual void InspectShape() const
print shape parameters
Definition TGeoTrd1.cxx:605
virtual void GetBoundingCylinder(Double_t *param) const
Fill vector param[4] with the bounding cylinder parameters.
Definition TGeoTrd1.cxx:512
virtual Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const
computes the closest distance from given point to this shape, according to option.
Definition TGeoTrd1.cxx:620
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...
Definition TGeoTrd1.cxx:749
virtual Bool_t Contains(const Double_t *point) const
test if point is inside this shape check Z range
Definition TGeoTrd1.cxx:187
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...
Definition TGeoTrd1.cxx:757
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.
Definition TGeoTrd1.cxx:739
virtual ~TGeoTrd1()
destructor
Definition TGeoTrd1.cxx:120
void SetVertex(Double_t *vertex) const
set vertex of a corner according to visibility flags
Definition TGeoTrd1.cxx:669
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 TGeoTrd1.cxx:582
void GetOppositeCorner(const Double_t *point, Int_t inorm, Double_t *vertex, Double_t *normals) const
get the opposite corner of the intersected face
Definition TGeoTrd1.cxx:293
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 trd1 Boundary safe algorithm.
Definition TGeoTrd1.cxx:202
void GetVisibleCorner(const Double_t *point, Double_t *vertex, Double_t *normals) const
get the most visible corner from outside point and the normals
Definition TGeoTrd1.cxx:252
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.
Definition TGeoTrd1.cxx:775
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const
Get range of shape for a given axis.
Definition TGeoTrd1.cxx:488
Double_t fDz
Definition TGeoTrd1.h:24
virtual TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
Divide this trd1 shape belonging to volume "voldiv" into ndiv volumes called divname,...
Definition TGeoTrd1.cxx:432
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Definition TGeoTrd1.cxx:642
TGeoTrd1()
dummy ctor
Definition TGeoTrd1.cxx:57
Double_t fDx1
Definition TGeoTrd1.h:21
virtual void SetDimensions(Double_t *param)
set trd1 params in one step :
Definition TGeoTrd1.cxx:657
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. Returns 0 if successful.
Definition TGeoTrd1.cxx:520
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 trd1 Boundary safe algorithm.
Definition TGeoTrd1.cxx:318
virtual Double_t Capacity() const
Computes capacity of the shape in [length^3].
Definition TGeoTrd1.cxx:127
Double_t fDy
Definition TGeoTrd1.h:23
virtual void SetPoints(Double_t *points) const
create arb8 mesh points
Definition TGeoTrd1.cxx:697
virtual void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
Compute normal to closest surface from POINT.
Definition TGeoTrd1.cxx:147
Double_t fDx2
Definition TGeoTrd1.h:22
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:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:949
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
Long64_t LocMin(Long64_t n, const T *a)
Return index of array with the minimum element.
Definition TMath.h:922
Short_t Max(Short_t a, Short_t b)
Definition TMathBase.h:208
T1 Sign(T1 a, T2 b)
Definition TMathBase.h:161
Long64_t LocMax(Long64_t n, const T *a)
Return index of array with the maximum element.
Definition TMath.h:950
Double_t Sqrt(Double_t x)
Definition TMath.h:641
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:176
Short_t Abs(Short_t d)
Definition TMathBase.h:120
REAL * vertex
Definition triangle.c:513