Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoPara.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 31/01/02
3// TGeoPara::Contains() 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 TGeoPara
14\ingroup Shapes_classes
15\brief Parallelepiped class.
16
17A parallelepiped is a shape having 3 pairs of parallel faces out of
18which one is parallel with the XY plane (Z faces). All faces are
19parallelograms in the general case. The Z faces have 2 edges parallel
20with the X-axis.
21
22The shape has the center in the origin and it is defined by:
23
24 - `dX, dY, dZ:` half-lengths of the projections of the edges on X, Y
25 and Z. The lower Z face is positioned at `-dZ`, while the upper at `+dZ`.
26 - `alpha:` angle between the segment defined by the centers of the
27 X-parallel edges and Y axis `[-90,90]` in degrees
28 - `theta:` theta angle of the segment defined by the centers of the Z faces;
29 - `phi:` phi angle of the same segment
30
31~~~ {.cpp}
32TGeoPara(dX,dY,dZ,alpha,theta,phi);
33~~~
34
35A box is a particular parallelepiped having the parameters:
36`(dX,dY,dZ,0.,0.,0.)`.
37
38Begin_Macro
39{
40 TCanvas *c = new TCanvas("c", "c",0,0,600,600);
41 new TGeoManager("para", "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->MakePara("PARA",med, 20,30,40,30,15,30);
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*/
58
59#include <iostream>
60
61#include "TGeoManager.h"
62#include "TGeoMatrix.h"
63#include "TGeoVolume.h"
64#include "TGeoPara.h"
65#include "TMath.h"
66
68
69////////////////////////////////////////////////////////////////////////////////
70/// Default constructor
71
73{
75 fX = fY = fZ = 0;
76 fAlpha = 0;
77 fTheta = 0;
78 fPhi = 0;
79 fTxy = 0;
80 fTxz = 0;
81 fTyz = 0;
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Default constructor specifying minimum and maximum radius
86
88 : TGeoBBox(0, 0, 0)
89{
90 SetShapeBit(TGeoShape::kGeoPara);
91 fX = dx;
92 fY = dy;
93 fZ = dz;
94 fAlpha = alpha;
95 fTheta = theta;
96 fPhi = phi;
97 fTxy = TMath::Tan(alpha * TMath::DegToRad());
98 Double_t tth = TMath::Tan(theta * TMath::DegToRad());
99 Double_t ph = phi * TMath::DegToRad();
100 fTxz = tth * TMath::Cos(ph);
101 fTyz = tth * TMath::Sin(ph);
102 if ((fX < 0) || (fY < 0) || (fZ < 0)) {
103 // printf("para : %f %f %f\n", fX, fY, fZ);
104 SetShapeBit(kGeoRunTimeShape);
105 } else
106 ComputeBBox();
107}
108
109////////////////////////////////////////////////////////////////////////////////
110/// Default constructor specifying minimum and maximum radius
111
112TGeoPara::TGeoPara(const char *name, Double_t dx, Double_t dy, Double_t dz, Double_t alpha, Double_t theta,
113 Double_t phi)
114 : TGeoBBox(name, 0, 0, 0)
115{
116 SetShapeBit(TGeoShape::kGeoPara);
117 fX = dx;
118 fY = dy;
119 fZ = dz;
120 fAlpha = alpha;
121 fTheta = theta;
122 fPhi = phi;
123 fTxy = TMath::Tan(alpha * TMath::DegToRad());
124 Double_t tth = TMath::Tan(theta * TMath::DegToRad());
125 Double_t ph = phi * TMath::DegToRad();
126 fTxz = tth * TMath::Cos(ph);
127 fTyz = tth * TMath::Sin(ph);
128 if ((fX < 0) || (fY < 0) || (fZ < 0)) {
129 // printf("para : %f %f %f\n", fX, fY, fZ);
130 SetShapeBit(kGeoRunTimeShape);
131 } else
132 ComputeBBox();
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// Default constructor
137/// - param[0] = dx
138/// - param[1] = dy
139/// - param[2] = dz
140/// - param[3] = alpha
141/// - param[4] = theta
142/// - param[5] = phi
143
144TGeoPara::TGeoPara(Double_t *param) : TGeoBBox(0, 0, 0)
145{
146 SetShapeBit(TGeoShape::kGeoPara);
147 SetDimensions(param);
148 if ((fX < 0) || (fY < 0) || (fZ < 0))
149 SetShapeBit(kGeoRunTimeShape);
150 else
151 ComputeBBox();
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// destructor
156
158
159////////////////////////////////////////////////////////////////////////////////
160/// Computes capacity of the shape in [length^3]
161
163{
164 Double_t capacity = 8. * fX * fY * fZ;
165 return capacity;
166}
167
168////////////////////////////////////////////////////////////////////////////////
169/// compute bounding box
170
172{
173 Double_t dx = fX + fY * TMath::Abs(fTxy) + fZ * TMath::Abs(fTxz);
174 Double_t dy = fY + fZ * TMath::Abs(fTyz);
175 Double_t dz = fZ;
176 TGeoBBox::SetBoxDimensions(dx, dy, dz);
177 memset(fOrigin, 0, 3 * sizeof(Double_t));
178}
179
180////////////////////////////////////////////////////////////////////////////////
181/// Compute normal to closest surface from POINT.
182
183void TGeoPara::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
184{
185 Double_t saf[3];
186 // distance from point to higher Z face
187 saf[0] = TMath::Abs(fZ - TMath::Abs(point[2])); // Z
188
189 Double_t yt = point[1] - fTyz * point[2];
190 saf[1] = TMath::Abs(fY - TMath::Abs(yt)); // Y
191 // cos of angle YZ
192 Double_t cty = 1.0 / TMath::Sqrt(1.0 + fTyz * fTyz);
193
194 Double_t xt = point[0] - fTxz * point[2] - fTxy * yt;
195 saf[2] = TMath::Abs(fX - TMath::Abs(xt)); // X
196 // cos of angle XZ
197 Double_t ctx = 1.0 / TMath::Sqrt(1.0 + fTxy * fTxy + fTxz * fTxz);
198 saf[2] *= ctx;
199 saf[1] *= cty;
200 Int_t i = TMath::LocMin(3, saf);
201 switch (i) {
202 case 0:
203 norm[0] = norm[1] = 0;
204 norm[2] = TMath::Sign(1., dir[2]);
205 return;
206 case 1:
207 norm[0] = 0;
208 norm[1] = cty;
209 norm[2] = -fTyz * cty;
210 break;
211 case 2:
214 norm[2] = -TMath::Sin(fTheta * TMath::DegToRad());
215 }
216 if (norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2] < 0) {
217 norm[0] = -norm[0];
218 norm[1] = -norm[1];
219 norm[2] = -norm[2];
220 }
221}
222
223////////////////////////////////////////////////////////////////////////////////
224/// test if point is inside this sphere
225/// test Z range
226
227Bool_t TGeoPara::Contains(const Double_t *point) const
228{
229 if (TMath::Abs(point[2]) > fZ)
230 return kFALSE;
231 // check X and Y
232 Double_t yt = point[1] - fTyz * point[2];
233 if (TMath::Abs(yt) > fY)
234 return kFALSE;
235 Double_t xt = point[0] - fTxz * point[2] - fTxy * yt;
236 if (TMath::Abs(xt) > fX)
237 return kFALSE;
238 return kTRUE;
239}
240
241////////////////////////////////////////////////////////////////////////////////
242/// compute distance from inside point to surface of the para
243/// Boundary safe algorithm.
244
246TGeoPara::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
247{
248 if (iact < 3 && safe) {
249 // compute safety
250 *safe = Safety(point, kTRUE);
251 if (iact == 0)
252 return TGeoShape::Big();
253 if (iact == 1 && step < *safe)
254 return TGeoShape::Big();
255 }
256 Double_t saf[2];
257 Double_t snxt = TGeoShape::Big();
258 Double_t s;
259 saf[0] = fZ + point[2];
260 saf[1] = fZ - point[2];
261 if (!TGeoShape::IsSameWithinTolerance(dir[2], 0)) {
262 s = (dir[2] > 0) ? (saf[1] / dir[2]) : (-saf[0] / dir[2]);
263 if (s < 0)
264 return 0.0;
265 if (s < snxt)
266 snxt = s;
267 }
268 // distance from point to center axis on Y
269 Double_t yt = point[1] - fTyz * point[2];
270 saf[0] = fY + yt;
271 saf[1] = fY - yt;
272 Double_t dy = dir[1] - fTyz * dir[2];
274 s = (dy > 0) ? (saf[1] / dy) : (-saf[0] / dy);
275 if (s < 0)
276 return 0.0;
277 if (s < snxt)
278 snxt = s;
279 }
280 // distance from point to center axis on X
281 Double_t xt = point[0] - fTxz * point[2] - fTxy * yt;
282 saf[0] = fX + xt;
283 saf[1] = fX - xt;
284 Double_t dx = dir[0] - fTxz * dir[2] - fTxy * dy;
286 s = (dx > 0) ? (saf[1] / dx) : (-saf[0] / dx);
287 if (s < 0)
288 return 0.0;
289 if (s < snxt)
290 snxt = s;
291 }
292 return snxt;
293}
294
295////////////////////////////////////////////////////////////////////////////////
296/// compute distance from inside point to surface of the para
297
299TGeoPara::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step, Double_t *safe) const
300{
301 if (iact < 3 && safe) {
302 // compute safe distance
303 *safe = Safety(point, kFALSE);
304 if (iact == 0)
305 return TGeoShape::Big();
306 if (iact == 1 && step < *safe)
307 return TGeoShape::Big();
308 }
309 Bool_t in = kTRUE;
310 Double_t safz;
311 safz = TMath::Abs(point[2]) - fZ;
312 if (safz > 0) {
313 // outside Z
314 if (point[2] * dir[2] >= 0)
315 return TGeoShape::Big();
316 in = kFALSE;
317 }
318 Double_t yt = point[1] - fTyz * point[2];
319 Double_t safy = TMath::Abs(yt) - fY;
320 Double_t dy = dir[1] - fTyz * dir[2];
321 if (safy > 0) {
322 if (yt * dy >= 0)
323 return TGeoShape::Big();
324 in = kFALSE;
325 }
326 Double_t xt = point[0] - fTxy * yt - fTxz * point[2];
327 Double_t safx = TMath::Abs(xt) - fX;
328 Double_t dx = dir[0] - fTxy * dy - fTxz * dir[2];
329 if (safx > 0) {
330 if (xt * dx >= 0)
331 return TGeoShape::Big();
332 in = kFALSE;
333 }
334 // protection in case point is actually inside
335 if (in) {
336 if (safz > safx && safz > safy) {
337 if (point[2] * dir[2] > 0)
338 return TGeoShape::Big();
339 return 0.0;
340 }
341 if (safx > safy) {
342 if (xt * dx > 0)
343 return TGeoShape::Big();
344 return 0.0;
345 }
346 if (yt * dy > 0)
347 return TGeoShape::Big();
348 return 0.0;
349 }
350 Double_t xnew, ynew, znew;
351 if (safz > 0) {
352 Double_t snxt = safz / TMath::Abs(dir[2]);
353 xnew = point[0] + snxt * dir[0];
354 ynew = point[1] + snxt * dir[1];
355 znew = (point[2] > 0) ? fZ : (-fZ);
356 Double_t ytn = ynew - fTyz * znew;
357 if (TMath::Abs(ytn) <= fY) {
358 Double_t xtn = xnew - fTxy * ytn - fTxz * znew;
359 if (TMath::Abs(xtn) <= fX)
360 return snxt;
361 }
362 }
363 if (safy > 0) {
364 Double_t snxt = safy / TMath::Abs(dy);
365 znew = point[2] + snxt * dir[2];
366 if (TMath::Abs(znew) <= fZ) {
367 Double_t ytn = (yt > 0) ? fY : (-fY);
368 xnew = point[0] + snxt * dir[0];
369 Double_t xtn = xnew - fTxy * ytn - fTxz * znew;
370 if (TMath::Abs(xtn) <= fX)
371 return snxt;
372 }
373 }
374 if (safx > 0) {
375 Double_t snxt = safx / TMath::Abs(dx);
376 znew = point[2] + snxt * dir[2];
377 if (TMath::Abs(znew) <= fZ) {
378 ynew = point[1] + snxt * dir[1];
379 Double_t ytn = ynew - fTyz * znew;
380 if (TMath::Abs(ytn) <= fY)
381 return snxt;
382 }
383 }
384 return TGeoShape::Big();
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// Divide this parallelepiped shape belonging to volume "voldiv" into ndiv equal volumes
389/// called divname, from start position with the given step. Returns pointer
390/// to created division cell volume. In case a wrong division axis is supplied,
391/// returns pointer to volume to be divided.
392
394TGeoPara::Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step)
395{
396 TGeoShape *shape; //--- shape to be created
397 TGeoVolume *vol; //--- division volume to be created
398 TGeoVolumeMulti *vmulti; //--- generic divided volume
399 TGeoPatternFinder *finder; //--- finder to be attached
400 TString opt = ""; //--- option to be attached
401 Double_t end = start + ndiv * step;
402 switch (iaxis) {
403 case 1: //--- divide on X
404 shape = new TGeoPara(step / 2, fY, fZ, fAlpha, fTheta, fPhi);
405 finder = new TGeoPatternParaX(voldiv, ndiv, start, end);
406 opt = "X";
407 break;
408 case 2: //--- divide on Y
409 shape = new TGeoPara(fX, step / 2, fZ, fAlpha, fTheta, fPhi);
410 finder = new TGeoPatternParaY(voldiv, ndiv, start, end);
411 opt = "Y";
412 break;
413 case 3: //--- divide on Z
414 shape = new TGeoPara(fX, fY, step / 2, fAlpha, fTheta, fPhi);
415 finder = new TGeoPatternParaZ(voldiv, ndiv, start, end);
416 opt = "Z";
417 break;
418 default: Error("Divide", "Wrong axis type for division"); return nullptr;
419 }
420 vol = new TGeoVolume(divname, shape, voldiv->GetMedium());
421 vmulti = gGeoManager->MakeVolumeMulti(divname, voldiv->GetMedium());
422 vmulti->AddVolume(vol);
423 voldiv->SetFinder(finder);
424 finder->SetDivIndex(voldiv->GetNdaughters());
425 for (Int_t ic = 0; ic < ndiv; ic++) {
426 voldiv->AddNodeOffset(vol, ic, start + step / 2. + ic * step, opt.Data());
427 ((TGeoNodeOffset *)voldiv->GetNodes()->At(voldiv->GetNdaughters() - 1))->SetFinder(finder);
428 }
429 return vmulti;
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Get range of shape for a given axis.
434
436{
437 xlo = 0;
438 xhi = 0;
439 Double_t dx = 0;
440 switch (iaxis) {
441 case 1:
442 xlo = -fX;
443 xhi = fX;
444 dx = xhi - xlo;
445 return dx;
446 case 2:
447 xlo = -fY;
448 xhi = fY;
449 dx = xhi - xlo;
450 return dx;
451 case 3:
452 xlo = -fZ;
453 xhi = fZ;
454 dx = xhi - xlo;
455 return dx;
456 }
457 return dx;
458}
459
460////////////////////////////////////////////////////////////////////////////////
461/// Fill vector param[4] with the bounding cylinder parameters. The order
462/// is the following : Rmin, Rmax, Phi1, Phi2
463
465{
467}
468
469////////////////////////////////////////////////////////////////////////////////
470/// Fills real parameters of a positioned box inside this. Returns 0 if successful.
471
472Int_t TGeoPara::GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const
473{
474 dx = dy = dz = 0;
475 if (mat->IsRotation()) {
476 Error("GetFittingBox", "cannot handle parametrized rotated volumes");
477 return 1; // ### rotation not accepted ###
478 }
479 //--> translate the origin of the parametrized box to the frame of this box.
480 Double_t origin[3];
481 mat->LocalToMaster(parambox->GetOrigin(), origin);
482 if (!Contains(origin)) {
483 Error("GetFittingBox", "wrong matrix - parametrized box is outside this");
484 return 1; // ### wrong matrix ###
485 }
486 //--> now we have to get the valid range for all parametrized axis
487 Double_t dd[3];
488 dd[0] = parambox->GetDX();
489 dd[1] = parambox->GetDY();
490 dd[2] = parambox->GetDZ();
491 //-> check if Z range is fixed
492 if (dd[2] < 0) {
493 dd[2] = TMath::Min(origin[2] + fZ, fZ - origin[2]);
494 if (dd[2] < 0) {
495 Error("GetFittingBox", "wrong matrix");
496 return 1;
497 }
498 }
499 if (dd[0] >= 0 && dd[1] >= 0) {
500 dx = dd[0];
501 dy = dd[1];
502 dz = dd[2];
503 return 0;
504 }
505 //-> check now range at Z = origin[2] +/- dd[2]
506 Double_t upper[8];
507 Double_t lower[8];
508 Double_t z = origin[2] - dd[2];
509 lower[0] = z * fTxz - fTxy * fY - fX;
510 lower[1] = -fY + z * fTyz;
511 lower[2] = z * fTxz + fTxy * fY - fX;
512 lower[3] = fY + z * fTyz;
513 lower[4] = z * fTxz + fTxy * fY + fX;
514 lower[5] = fY + z * fTyz;
515 lower[6] = z * fTxz - fTxy * fY + fX;
516 lower[7] = -fY + z * fTyz;
517 z = origin[2] + dd[2];
518 upper[0] = z * fTxz - fTxy * fY - fX;
519 upper[1] = -fY + z * fTyz;
520 upper[2] = z * fTxz + fTxy * fY - fX;
521 upper[3] = fY + z * fTyz;
522 upper[4] = z * fTxz + fTxy * fY + fX;
523 upper[5] = fY + z * fTyz;
524 upper[6] = z * fTxz - fTxy * fY + fX;
525 upper[7] = -fY + z * fTyz;
526
527 for (Int_t iaxis = 0; iaxis < 2; iaxis++) {
528 if (dd[iaxis] >= 0)
529 continue;
530 Double_t ddmin = TGeoShape::Big();
531 for (Int_t ivert = 0; ivert < 4; ivert++) {
532 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis] - lower[2 * ivert + iaxis]));
533 ddmin = TMath::Min(ddmin, TMath::Abs(origin[iaxis] - upper[2 * ivert + iaxis]));
534 }
535 dd[iaxis] = ddmin;
536 }
537 dx = dd[0];
538 dy = dd[1];
539 dz = dd[2];
540 return 0;
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// in case shape has some negative parameters, these has to be computed
545/// in order to fit the mother
546
548{
550 return nullptr;
551 if (!mother->TestShapeBit(kGeoPara)) {
552 Error("GetMakeRuntimeShape", "invalid mother");
553 return nullptr;
554 }
555 Double_t dx, dy, dz;
556 if (fX < 0)
557 dx = ((TGeoPara *)mother)->GetX();
558 else
559 dx = fX;
560 if (fY < 0)
561 dy = ((TGeoPara *)mother)->GetY();
562 else
563 dy = fY;
564 if (fZ < 0)
565 dz = ((TGeoPara *)mother)->GetZ();
566 else
567 dz = fZ;
568 return (new TGeoPara(dx, dy, dz, fAlpha, fTheta, fPhi));
569}
570
571////////////////////////////////////////////////////////////////////////////////
572/// print shape parameters
573
574void TGeoPara::InspectShape() const
575{
576 printf("*** Shape %s: TGeoPara ***\n", GetName());
577 printf(" dX = %11.5f\n", fX);
578 printf(" dY = %11.5f\n", fY);
579 printf(" dZ = %11.5f\n", fZ);
580 printf(" alpha = %11.5f\n", fAlpha);
581 printf(" theta = %11.5f\n", fTheta);
582 printf(" phi = %11.5f\n", fPhi);
583 printf(" Bounding box:\n");
585}
586
587////////////////////////////////////////////////////////////////////////////////
588/// computes the closest distance from given point to this shape, according
589/// to option. The matching point on the shape is stored in spoint.
590
591Double_t TGeoPara::Safety(const Double_t *point, Bool_t in) const
592{
593 Double_t saf[3];
594 // distance from point to higher Z face
595 saf[0] = fZ - TMath::Abs(point[2]); // Z
596
597 Double_t yt = point[1] - fTyz * point[2];
598 saf[1] = fY - TMath::Abs(yt); // Y
599 // cos of angle YZ
600 Double_t cty = 1.0 / TMath::Sqrt(1.0 + fTyz * fTyz);
601
602 Double_t xt = point[0] - fTxz * point[2] - fTxy * yt;
603 saf[2] = fX - TMath::Abs(xt); // X
604 // cos of angle XZ
605 Double_t ctx = 1.0 / TMath::Sqrt(1.0 + fTxy * fTxy + fTxz * fTxz);
606 saf[2] *= ctx;
607 saf[1] *= cty;
608 if (in)
609 return saf[TMath::LocMin(3, saf)];
610 for (Int_t i = 0; i < 3; i++)
611 saf[i] = -saf[i];
612 return saf[TMath::LocMax(3, saf)];
613}
614
615////////////////////////////////////////////////////////////////////////////////
616/// Save a primitive as a C++ statement(s) on output stream "out".
617
618void TGeoPara::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
619{
621 return;
622 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
623 out << " dx = " << fX << ";" << std::endl;
624 out << " dy = " << fY << ";" << std::endl;
625 out << " dz = " << fZ << ";" << std::endl;
626 out << " alpha = " << fAlpha << ";" << std::endl;
627 out << " theta = " << fTheta << ";" << std::endl;
628 out << " phi = " << fPhi << ";" << std::endl;
629 out << " TGeoShape *" << GetPointerName() << " = new TGeoPara(\"" << GetName() << "\",dx,dy,dz,alpha,theta,phi);"
630 << std::endl;
632}
633
634////////////////////////////////////////////////////////////////////////////////
635/// Set dimensions starting from an array.
636
638{
639 fX = param[0];
640 fY = param[1];
641 fZ = param[2];
642 fAlpha = param[3];
643 fTheta = param[4];
644 fPhi = param[5];
645 fTxy = TMath::Tan(param[3] * TMath::DegToRad());
646 Double_t tth = TMath::Tan(param[4] * TMath::DegToRad());
647 Double_t ph = param[5] * TMath::DegToRad();
648 fTxz = tth * TMath::Cos(ph);
649 fTyz = tth * TMath::Sin(ph);
650}
651
652////////////////////////////////////////////////////////////////////////////////
653/// Create PARA mesh points
654
656{
657 if (!points)
658 return;
659 Double_t txy = fTxy;
660 Double_t txz = fTxz;
661 Double_t tyz = fTyz;
662 *points++ = -fZ * txz - txy * fY - fX;
663 *points++ = -fY - fZ * tyz;
664 *points++ = -fZ;
665 *points++ = -fZ * txz + txy * fY - fX;
666 *points++ = +fY - fZ * tyz;
667 *points++ = -fZ;
668 *points++ = -fZ * txz + txy * fY + fX;
669 *points++ = +fY - fZ * tyz;
670 *points++ = -fZ;
671 *points++ = -fZ * txz - txy * fY + fX;
672 *points++ = -fY - fZ * tyz;
673 *points++ = -fZ;
674 *points++ = +fZ * txz - txy * fY - fX;
675 *points++ = -fY + fZ * tyz;
676 *points++ = +fZ;
677 *points++ = +fZ * txz + txy * fY - fX;
678 *points++ = +fY + fZ * tyz;
679 *points++ = +fZ;
680 *points++ = +fZ * txz + txy * fY + fX;
681 *points++ = +fY + fZ * tyz;
682 *points++ = +fZ;
683 *points++ = +fZ * txz - txy * fY + fX;
684 *points++ = -fY + fZ * tyz;
685 *points++ = +fZ;
686}
687
688////////////////////////////////////////////////////////////////////////////////
689/// create sphere mesh points
690
692{
693 if (!points)
694 return;
695 Double_t txy = fTxy;
696 Double_t txz = fTxz;
697 Double_t tyz = fTyz;
698 *points++ = -fZ * txz - txy * fY - fX;
699 *points++ = -fY - fZ * tyz;
700 *points++ = -fZ;
701 *points++ = -fZ * txz + txy * fY - fX;
702 *points++ = +fY - fZ * tyz;
703 *points++ = -fZ;
704 *points++ = -fZ * txz + txy * fY + fX;
705 *points++ = +fY - fZ * tyz;
706 *points++ = -fZ;
707 *points++ = -fZ * txz - txy * fY + fX;
708 *points++ = -fY - fZ * tyz;
709 *points++ = -fZ;
710 *points++ = +fZ * txz - txy * fY - fX;
711 *points++ = -fY + fZ * tyz;
712 *points++ = +fZ;
713 *points++ = +fZ * txz + txy * fY - fX;
714 *points++ = +fY + fZ * tyz;
715 *points++ = +fZ;
716 *points++ = +fZ * txz + txy * fY + fX;
717 *points++ = +fY + fZ * tyz;
718 *points++ = +fZ;
719 *points++ = +fZ * txz - txy * fY + fX;
720 *points++ = -fY + fZ * tyz;
721 *points++ = +fZ;
722}
723
724////////////////////////////////////////////////////////////////////////////////
725/// fill size of this 3-D object
726
727void TGeoPara::Sizeof3D() const
728{
730}
731
732////////////////////////////////////////////////////////////////////////////////
733/// Check the inside status for each of the points in the array.
734/// Input: Array of point coordinates + vector size
735/// Output: Array of Booleans for the inside of each point
736
737void TGeoPara::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
738{
739 for (Int_t i = 0; i < vecsize; i++)
740 inside[i] = Contains(&points[3 * i]);
741}
742
743////////////////////////////////////////////////////////////////////////////////
744/// Compute the normal for an array o points so that norm.dot.dir is positive
745/// Input: Arrays of point coordinates and directions + vector size
746/// Output: Array of normal directions
747
748void TGeoPara::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
749{
750 for (Int_t i = 0; i < vecsize; i++)
751 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 TGeoPara::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
758 Double_t *step) const
759{
760 for (Int_t i = 0; i < vecsize; i++)
761 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
762}
763
764////////////////////////////////////////////////////////////////////////////////
765/// Compute distance from array of input points having directions specified by dirs. Store output in dists
766
767void TGeoPara::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
768 Double_t *step) const
769{
770 for (Int_t i = 0; i < vecsize; i++)
771 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
772}
773
774////////////////////////////////////////////////////////////////////////////////
775/// Compute safe distance from each of the points in the input array.
776/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
777/// Output: Safety values
778
779void TGeoPara::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
780{
781 for (Int_t i = 0; i < vecsize; i++)
782 safe[i] = Safety(&points[3 * i], inside[i]);
783}
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
virtual const Double_t * GetOrigin() const
Definition TGeoBBox.h:82
void GetBoundingCylinder(Double_t *param) const override
void SetBoxDimensions(Double_t dx, Double_t dy, Double_t dz, Double_t *origin=nullptr)
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
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
void SetDimensions(Double_t *param) override
Double_t fTxy
Definition TGeoPara.h:26
Double_t fTyz
Definition TGeoPara.h:28
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
Double_t fX
Definition TGeoPara.h:20
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) override
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
void Sizeof3D() const override
TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const override
void GetBoundingCylinder(Double_t *param) const override
Double_t Capacity() const override
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const override
void InspectShape() const override
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
Double_t fTheta
Definition TGeoPara.h:24
Int_t GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const override
Double_t fTxz
Definition TGeoPara.h:27
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
void SetPoints(Double_t *points) const override
Double_t fAlpha
Definition TGeoPara.h:23
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
Double_t fZ
Definition TGeoPara.h:22
~TGeoPara() override
Double_t fY
Definition TGeoPara.h:21
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 fPhi
Definition TGeoPara.h:25
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
Bool_t Contains(const Double_t *point) const override
Double_t 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 ComputeBBox() override
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.
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.
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:64
@ kGeoRunTimeShape
Definition TGeoShape.h:40
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:167
Volume families.
Definition TGeoVolume.h:266
void AddVolume(TGeoVolume *vol)
Add a volume with valid shape to the list of volumes.
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
void AddNodeOffset(TGeoVolume *vol, Int_t copy_no, Double_t offset=0, Option_t *option="")
Add a division node to the list of nodes.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:175
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:244
Int_t GetNdaughters() const
Definition TGeoVolume.h:362
TObjArray * GetNodes()
Definition TGeoVolume.h:169
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:213
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
Long64_t LocMin(Long64_t n, const T *a)
Returns index of array with the minimum element.
Definition TMath.h:982
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
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
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
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:600
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123