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