Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoParaboloid.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Mihaela Gheata 20/06/04
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoParaboloid
13\ingroup Shapes_classes
14
15A paraboloid is defined by the revolution surface generated by a
16parabola and is bounded by two planes perpendicular to Z axis. The
17parabola equation is taken in the form: `z = a·r2 + b`, where:
18`r2 = x2 + y2`. Note the missing linear term (parabola symmetric with
19respect to Z axis).
20
21The coefficients a and b are computed from the input values which are
22the radii of the circular sections cut by the planes at `+/-dz`:
23
24 - `-dz = a*r2low + b`
25 - ` dz = a*r2high + b`
26
27~~~{.cpp}
28TGeoParaboloid(Double_t rlo,Double_t rhi,Double_t dz);
29~~~
30
31Begin_Macro
32{
33 new TGeoManager("parab", "paraboloid");
34 TGeoMaterial *mat = new TGeoMaterial("Al", 26.98,13,2.7);
35 TGeoMedium *med = new TGeoMedium("MED",1,mat);
36 TGeoVolume *top = gGeoManager->MakeBox("TOP",med,100,100,100);
37 gGeoManager->SetTopVolume(top);
38 TGeoVolume *vol = gGeoManager->MakeParaboloid("PARAB",med,0, 40, 50);
39 TGeoParaboloid *par = (TGeoParaboloid*)vol->GetShape();
40 top->AddNode(vol,1);
41 gGeoManager->CloseGeometry();
42 gGeoManager->SetNsegments(80);
43 top->Draw();
44 TView *view = gPad->GetView();
45 view->ShowAxis();
46}
47End_Macro
48*/
49
50#include <iostream>
51#include "TGeoManager.h"
52#include "TGeoVolume.h"
53#include "TVirtualGeoPainter.h"
54#include "TGeoParaboloid.h"
55#include "TBuffer3D.h"
56#include "TBuffer3DTypes.h"
57#include "TMath.h"
58
60
61////////////////////////////////////////////////////////////////////////////////
62/// Dummy constructor
63
65{
66 fRlo = 0;
67 fRhi = 0;
68 fDz = 0;
69 fA = 0;
70 fB = 0;
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// Default constructor specifying X and Y semiaxis length
76
78{
79 fRlo = 0;
80 fRhi = 0;
81 fDz = 0;
82 fA = 0;
83 fB = 0;
84 SetShapeBit(TGeoShape::kGeoParaboloid);
85 SetParaboloidDimensions(rlo, rhi, dz);
86 ComputeBBox();
87}
88
89////////////////////////////////////////////////////////////////////////////////
90/// Default constructor specifying X and Y semiaxis length
91
92TGeoParaboloid::TGeoParaboloid(const char *name, Double_t rlo, Double_t rhi, Double_t dz) : TGeoBBox(name, 0, 0, 0)
93{
94 fRlo = 0;
95 fRhi = 0;
96 fDz = 0;
97 fA = 0;
98 fB = 0;
99 SetShapeBit(TGeoShape::kGeoParaboloid);
100 SetParaboloidDimensions(rlo, rhi, dz);
101 ComputeBBox();
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// Default constructor specifying minimum and maximum radius
106/// - param[0] = rlo
107/// - param[1] = rhi
108/// - param[2] = dz
109
111{
113 SetDimensions(param);
114 ComputeBBox();
115}
116
117////////////////////////////////////////////////////////////////////////////////
118/// destructor
119
121
122////////////////////////////////////////////////////////////////////////////////
123/// Computes capacity of the shape in [length^3]
124
126{
127 Double_t capacity = TMath::Pi() * fDz * (fRlo * fRlo + fRhi * fRhi);
128 return capacity;
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// compute bounding box of the tube
133
135{
137 fDY = fDX;
138 fDZ = fDz;
139}
140
141////////////////////////////////////////////////////////////////////////////////
142/// Compute normal to closest surface from POINT.
143
144void TGeoParaboloid::ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm)
145{
146 norm[0] = norm[1] = 0.0;
147 if (TMath::Abs(point[2]) > fDz) {
148 norm[2] = TMath::Sign(1., dir[2]);
149 return;
150 }
151 Double_t safz = fDz - TMath::Abs(point[2]);
152 Double_t r = TMath::Sqrt(point[0] * point[0] + point[1] * point[1]);
153 Double_t safr = TMath::Abs(r - TMath::Sqrt((point[2] - fB) / fA));
154 if (safz < safr) {
155 norm[2] = TMath::Sign(1., dir[2]);
156 return;
157 }
158 Double_t talf = -2. * fA * r;
159 Double_t calf = 1. / TMath::Sqrt(1. + talf * talf);
160 Double_t salf = talf * calf;
161 Double_t phi = TMath::ATan2(point[1], point[0]);
162
163 norm[0] = salf * TMath::Cos(phi);
164 norm[1] = salf * TMath::Sin(phi);
165 norm[2] = calf;
166 Double_t ndotd = norm[0] * dir[0] + norm[1] * dir[1] + norm[2] * dir[2];
167 if (ndotd < 0) {
168 norm[0] = -norm[0];
169 norm[1] = -norm[1];
170 norm[2] = -norm[2];
171 }
172}
173
174////////////////////////////////////////////////////////////////////////////////
175/// test if point is inside the elliptical tube
176
177Bool_t TGeoParaboloid::Contains(const Double_t *point) const
178{
179 if (TMath::Abs(point[2]) > fDz)
180 return kFALSE;
181 Double_t aa = fA * (point[2] - fB);
182 if (aa < 0)
183 return kFALSE;
184 Double_t rsq = point[0] * point[0] + point[1] * point[1];
185 if (aa < fA * fA * rsq)
186 return kFALSE;
187 return kTRUE;
188}
189
190////////////////////////////////////////////////////////////////////////////////
191/// compute closest distance from point px,py to each vertex
192
194{
196 const Int_t numPoints = n * (n + 1) + 2;
197 return ShapeDistancetoPrimitive(numPoints, px, py);
198}
199
200////////////////////////////////////////////////////////////////////////////////
201/// Compute distance from a point to the parabola given by:
202/// `z = a*rsq + b; rsq = x*x+y*y`
203
204Double_t TGeoParaboloid::DistToParaboloid(const Double_t *point, const Double_t *dir, Bool_t in) const
205{
206 Double_t rsq = point[0] * point[0] + point[1] * point[1];
207 Double_t a = fA * (dir[0] * dir[0] + dir[1] * dir[1]);
208 Double_t b = 2. * fA * (point[0] * dir[0] + point[1] * dir[1]) - dir[2];
209 Double_t c = fA * rsq + fB - point[2];
213 return dist; // big
214 dist = -c / b;
215 if (dist < 0)
216 return TGeoShape::Big();
217 return dist; // OK
218 }
219 Double_t ainv = 1. / a;
220 Double_t sum = -b * ainv;
221 Double_t prod = c * ainv;
222 Double_t delta = sum * sum - 4. * prod;
223 if (delta < 0)
224 return dist; // big
225 delta = TMath::Sqrt(delta);
226 Double_t sone = TMath::Sign(1., ainv);
227 Int_t i = -1;
228 while (i < 2) {
229 dist = 0.5 * (sum + i * sone * delta);
230 i += 2;
231 if (dist < 0)
232 continue;
233 if (dist < 1.E-8) {
234 Double_t talf = -2. * fA * TMath::Sqrt(rsq);
235 Double_t phi = TMath::ATan2(point[1], point[0]);
236 Double_t ndotd = talf * (TMath::Cos(phi) * dir[0] + TMath::Sin(phi) * dir[1]) + dir[2];
237 if (!in)
238 ndotd *= -1;
239 if (ndotd < 0)
240 return dist;
241 } else
242 return dist;
243 }
244 return TGeoShape::Big();
245}
246
247////////////////////////////////////////////////////////////////////////////////
248/// compute distance from inside point to surface of the paraboloid
249
250Double_t TGeoParaboloid::DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step,
251 Double_t *safe) const
252{
253 if (iact < 3 && safe) {
254 // compute safe distance
255 *safe = Safety(point, kTRUE);
256 if (iact == 0)
257 return TGeoShape::Big();
258 if (iact == 1 && step < *safe)
259 return TGeoShape::Big();
260 }
261
263 if (dir[2] < 0) {
264 dz = -(point[2] + fDz) / dir[2];
265 } else if (dir[2] > 0) {
266 dz = (fDz - point[2]) / dir[2];
267 }
268 Double_t dpara = DistToParaboloid(point, dir, kTRUE);
269 return TMath::Min(dz, dpara);
270}
271
272////////////////////////////////////////////////////////////////////////////////
273/// compute distance from outside point to surface of the paraboloid and safe distance
274
275Double_t TGeoParaboloid::DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact, Double_t step,
276 Double_t *safe) const
277{
278 if (iact < 3 && safe) {
279 // compute safe distance
280 *safe = Safety(point, kFALSE);
281 if (iact == 0)
282 return TGeoShape::Big();
283 if (iact == 1 && step < *safe)
284 return TGeoShape::Big();
285 }
286 Double_t xnew, ynew, znew;
287 if (point[2] <= -fDz) {
288 if (dir[2] <= 0)
289 return TGeoShape::Big();
290 Double_t snxt = -(fDz + point[2]) / dir[2];
291 // find extrapolated X and Y
292 xnew = point[0] + snxt * dir[0];
293 ynew = point[1] + snxt * dir[1];
294 if ((xnew * xnew + ynew * ynew) <= fRlo * fRlo)
295 return snxt;
296 } else if (point[2] >= fDz) {
297 if (dir[2] >= 0)
298 return TGeoShape::Big();
299 Double_t snxt = (fDz - point[2]) / dir[2];
300 // find extrapolated X and Y
301 xnew = point[0] + snxt * dir[0];
302 ynew = point[1] + snxt * dir[1];
303 if ((xnew * xnew + ynew * ynew) <= fRhi * fRhi)
304 return snxt;
305 }
306 Double_t snxt = DistToParaboloid(point, dir, kFALSE);
307 if (snxt > 1E20)
308 return snxt;
309 znew = point[2] + snxt * dir[2];
310 if (TMath::Abs(znew) <= fDz)
311 return snxt;
312 return TGeoShape::Big();
313}
314
315////////////////////////////////////////////////////////////////////////////////
316/// Divide the paraboloid along one axis.
317
318TGeoVolume *TGeoParaboloid::Divide(TGeoVolume * /*voldiv*/, const char * /*divname*/, Int_t /*iaxis*/, Int_t /*ndiv*/,
319 Double_t /*start*/, Double_t /*step*/)
320{
321 Error("Divide", "Paraboloid divisions not implemented");
322 return nullptr;
323}
324
325////////////////////////////////////////////////////////////////////////////////
326/// Fill vector param[4] with the bounding cylinder parameters. The order
327/// is the following : Rmin, Rmax, Phi1, Phi2
328
330{
331 param[0] = 0.; // Rmin
332 param[1] = fDX; // Rmax
333 param[1] *= param[1];
334 param[2] = 0.; // Phi1
335 param[3] = 360.; // Phi2
336}
337
338////////////////////////////////////////////////////////////////////////////////
339/// in case shape has some negative parameters, these has to be computed
340/// in order to fit the mother
341
343{
344 return nullptr;
345}
346
347////////////////////////////////////////////////////////////////////////////////
348/// print shape parameters
349
351{
352 printf("*** Shape %s: TGeoParaboloid ***\n", GetName());
353 printf(" rlo = %11.5f\n", fRlo);
354 printf(" rhi = %11.5f\n", fRhi);
355 printf(" dz = %11.5f\n", fDz);
356 printf(" Bounding box:\n");
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// Creates a TBuffer3D describing *this* shape.
362/// Coordinates are in local reference frame.
363
365{
367 Int_t nbPnts = n * (n + 1) + 2;
368 Int_t nbSegs = n * (2 * n + 3);
369 Int_t nbPols = n * (n + 2);
370
371 TBuffer3D *buff =
372 new TBuffer3D(TBuffer3DTypes::kGeneric, nbPnts, 3 * nbPnts, nbSegs, 3 * nbSegs, nbPols, 2 * n * 5 + n * n * 6);
373
374 if (buff) {
375 SetPoints(buff->fPnts);
376 SetSegsAndPols(*buff);
377 }
378
379 return buff;
380}
381
382////////////////////////////////////////////////////////////////////////////////
383/// Fill TBuffer3D structure for segments and polygons.
384
386{
387 Int_t indx, i, j;
389
391
392 Int_t nn1 = (n + 1) * n + 1;
393 indx = 0;
394 // Lower end-cap (n radial segments)
395 for (j = 0; j < n; j++) {
396 buff.fSegs[indx++] = c + 2;
397 buff.fSegs[indx++] = 0;
398 buff.fSegs[indx++] = j + 1;
399 }
400 // Sectors (n)
401 for (i = 0; i < n + 1; i++) {
402 // lateral (circles) segments (n)
403 for (j = 0; j < n; j++) {
404 buff.fSegs[indx++] = c;
405 buff.fSegs[indx++] = n * i + 1 + j;
406 buff.fSegs[indx++] = n * i + 1 + ((j + 1) % n);
407 }
408 if (i == n)
409 break; // skip i=n for generators
410 // generator segments (n)
411 for (j = 0; j < n; j++) {
412 buff.fSegs[indx++] = c;
413 buff.fSegs[indx++] = n * i + 1 + j;
414 buff.fSegs[indx++] = n * (i + 1) + 1 + j;
415 }
416 }
417 // Upper end-cap
418 for (j = 0; j < n; j++) {
419 buff.fSegs[indx++] = c + 1;
420 buff.fSegs[indx++] = n * n + 1 + j;
421 buff.fSegs[indx++] = nn1;
422 }
423
424 indx = 0;
425
426 // lower end-cap (n polygons)
427 for (j = 0; j < n; j++) {
428 buff.fPols[indx++] = c + 2;
429 buff.fPols[indx++] = 3;
430 buff.fPols[indx++] = n + j;
431 buff.fPols[indx++] = (j + 1) % n;
432 buff.fPols[indx++] = j;
433 }
434 // Sectors (n)
435 for (i = 0; i < n; i++) {
436 // lateral faces (n)
437 for (j = 0; j < n; j++) {
438 buff.fPols[indx++] = c;
439 buff.fPols[indx++] = 4;
440 buff.fPols[indx++] = (2 * i + 1) * n + j;
441 buff.fPols[indx++] = 2 * (i + 1) * n + j;
442 buff.fPols[indx++] = (2 * i + 3) * n + j;
443 buff.fPols[indx++] = 2 * (i + 1) * n + ((j + 1) % n);
444 }
445 }
446 // upper end-cap (n polygons)
447 for (j = 0; j < n; j++) {
448 buff.fPols[indx++] = c + 1;
449 buff.fPols[indx++] = 3;
450 buff.fPols[indx++] = 2 * n * (n + 1) + j;
451 buff.fPols[indx++] = 2 * n * (n + 1) + ((j + 1) % n);
452 buff.fPols[indx++] = (2 * n + 1) * n + j;
453 }
454}
455
456////////////////////////////////////////////////////////////////////////////////
457/// Computes the closest distance from given point to this shape.
458
459Double_t TGeoParaboloid::Safety(const Double_t *point, Bool_t in) const
460{
461 Double_t safz = fDz - TMath::Abs(point[2]);
462 if (!in)
463 safz = -safz;
464 Double_t safr;
465 Double_t rsq = point[0] * point[0] + point[1] * point[1];
466 Double_t z0 = fA * rsq + fB;
467 Double_t r0sq = (point[2] - fB) / fA;
468 if (r0sq < 0) {
469 if (in)
470 return 0.;
471 return safz;
472 }
473 Double_t dr = TMath::Sqrt(rsq) - TMath::Sqrt(r0sq);
474 if (in) {
475 if (dr > -1.E-8)
476 return 0.;
477 Double_t dz = TMath::Abs(point[2] - z0);
478 safr = -dr * dz / TMath::Sqrt(dr * dr + dz * dz);
479 } else {
480 if (dr < 1.E-8)
481 return safz;
482 Double_t talf = -2. * fA * TMath::Sqrt(r0sq);
483 Double_t salf = talf / TMath::Sqrt(1. + talf * talf);
484 safr = TMath::Abs(dr * salf);
485 }
486 if (in)
487 return TMath::Min(safr, safz);
488 return TMath::Max(safr, safz);
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Set paraboloid dimensions.
493
495{
496 if ((rlo < 0) || (rhi < 0) || (dz <= 0) || TMath::Abs(rlo - rhi) < TGeoShape::Tolerance()) {
498 Error("SetParaboloidDimensions", "Dimensions of %s invalid: check (rlo>=0) (rhi>=0) (rlo!=rhi) dz>0", GetName());
499 return;
500 }
501 fRlo = rlo;
502 fRhi = rhi;
503 fDz = dz;
504 Double_t dd = 1. / (fRhi * fRhi - fRlo * fRlo);
505 fA = 2. * fDz * dd;
506 fB = -fDz * (fRlo * fRlo + fRhi * fRhi) * dd;
507}
508
509////////////////////////////////////////////////////////////////////////////////
510/// Set paraboloid dimensions starting from an array.
511
513{
514 Double_t rlo = param[0];
515 Double_t rhi = param[1];
516 Double_t dz = param[2];
517 SetParaboloidDimensions(rlo, rhi, dz);
518}
519
520////////////////////////////////////////////////////////////////////////////////
521/// Create paraboloid mesh points.
522/// ~~~ {.cpp}
523/// Npoints = n*(n+1) + 2
524/// ifirst = 0
525/// ipoint(i,j) = 1+i*n+j; i=[0,n] j=[0,n-1]
526/// ilast = 1+n*(n+1)
527/// Nsegments = n*(2*n+3)
528/// lower: (0, j+1); j=[0,n-1]
529/// circle(i): (n*i+1+j, n*i+1+(j+1)%n); i=[0,n] j=[0,n-1]
530/// generator(i): (n*i+1+j, n*(i+1)+1+j); i,j=[0,n-1]
531/// upper: (n*n+1+j, (n+1)*n+1) j=[0,n-1]
532/// Npolygons = n*(n+2)
533/// lower: (n+j, (j+1)%n, j) j=[0,n-1]
534/// lateral(i): ((2*i+1)*n+j, 2*(i+1)*n+j, (2*i+3)*n+j, 2*(i+1)*n+(j+1)%n)
535/// i,j = [0,n-1]
536/// upper: ((2n+1)*n+j, 2*n*(n+1)+(j+1)%n, 2*n*(n+1)+j) j=[0,n-1]
537/// ~~~
538
540{
541 if (!points)
542 return;
543 Double_t ttmin, ttmax;
544 ttmin = TMath::ATan2(-fDz, fRlo);
545 ttmax = TMath::ATan2(fDz, fRhi);
547 Double_t dtt = (ttmax - ttmin) / n;
548 Double_t dphi = 360. / n;
549 Double_t tt;
550 Double_t r, z, delta;
551 Double_t phi, sph, cph;
552 Int_t indx = 0;
553 // center of the lower endcap:
554 points[indx++] = 0; // x
555 points[indx++] = 0; // y
556 points[indx++] = -fDz;
557 for (Int_t i = 0; i < n + 1; i++) { // nz planes = n+1
558 if (i == 0) {
559 r = fRlo;
560 z = -fDz;
561 } else if (i == n) {
562 r = fRhi;
563 z = fDz;
564 } else {
565 tt = TMath::Tan(ttmin + i * dtt);
566 delta = tt * tt - 4 * fA * fB; // should be always positive (a*b<0)
567 r = 0.5 * (tt + TMath::Sqrt(delta)) / fA;
568 z = r * tt;
569 }
570 for (Int_t j = 0; j < n; j++) {
571 phi = j * dphi * TMath::DegToRad();
572 sph = TMath::Sin(phi);
573 cph = TMath::Cos(phi);
574 points[indx++] = r * cph;
575 points[indx++] = r * sph;
576 points[indx++] = z;
577 }
578 }
579 // center of the upper endcap
580 points[indx++] = 0; // x
581 points[indx++] = 0; // y
582 points[indx++] = fDz;
583}
584
585////////////////////////////////////////////////////////////////////////////////
586/// Returns numbers of vertices, segments and polygons composing the shape mesh.
587
588void TGeoParaboloid::GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const
589{
591 nvert = n * (n + 1) + 2;
592 nsegs = n * (2 * n + 3);
593 npols = n * (n + 2);
594}
595
596////////////////////////////////////////////////////////////////////////////////
597/// Returns number of vertices on the paraboloid mesh.
598
600{
602 return (n * (n + 1) + 2);
603}
604
605////////////////////////////////////////////////////////////////////////////////
606/// Save a primitive as a C++ statement(s) on output stream "out".
607
608void TGeoParaboloid::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
609{
611 return;
612 out << " // Shape: " << GetName() << " type: " << ClassName() << std::endl;
613 out << " rlo = " << fRlo << ";" << std::endl;
614 out << " rhi = " << fRhi << ";" << std::endl;
615 out << " dz = " << fDZ << ";" << std::endl;
616 out << " TGeoShape *" << GetPointerName() << " = new TGeoParaboloid(\"" << GetName() << "\", rlo,rhi,dz);"
617 << std::endl;
619}
620
621////////////////////////////////////////////////////////////////////////////////
622/// Create paraboloid mesh points.
623
625{
626 if (!points)
627 return;
628 Double_t ttmin, ttmax;
629 ttmin = TMath::ATan2(-fDz, fRlo);
630 ttmax = TMath::ATan2(fDz, fRhi);
632 Double_t dtt = (ttmax - ttmin) / n;
633 Double_t dphi = 360. / n;
634 Double_t tt;
635 Double_t r, z, delta;
636 Double_t phi, sph, cph;
637 Int_t indx = 0;
638 // center of the lower endcap:
639 points[indx++] = 0; // x
640 points[indx++] = 0; // y
641 points[indx++] = -fDz;
642 for (Int_t i = 0; i < n + 1; i++) { // nz planes = n+1
643 if (i == 0) {
644 r = fRlo;
645 z = -fDz;
646 } else if (i == n) {
647 r = fRhi;
648 z = fDz;
649 } else {
650 tt = TMath::Tan(ttmin + i * dtt);
651 delta = tt * tt - 4 * fA * fB; // should be always positive (a*b<0)
652 r = 0.5 * (tt + TMath::Sqrt(delta)) / fA;
653 z = r * tt;
654 }
655 for (Int_t j = 0; j < n; j++) {
656 phi = j * dphi * TMath::DegToRad();
657 sph = TMath::Sin(phi);
658 cph = TMath::Cos(phi);
659 points[indx++] = r * cph;
660 points[indx++] = r * sph;
661 points[indx++] = z;
662 }
663 }
664 // center of the upper endcap
665 points[indx++] = 0; // x
666 points[indx++] = 0; // y
667 points[indx++] = fDz;
668}
669
670////////////////////////////////////////////////////////////////////////////////
671void TGeoParaboloid::Sizeof3D() const {}
672
673////////////////////////////////////////////////////////////////////////////////
674/// Fills a static 3D buffer and returns a reference.
675
676const TBuffer3D &TGeoParaboloid::GetBuffer3D(Int_t reqSections, Bool_t localFrame) const
677{
678 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
679 TGeoBBox::FillBuffer3D(buffer, reqSections, localFrame);
680
681 if (reqSections & TBuffer3D::kRawSizes) {
683 Int_t nbPnts = n * (n + 1) + 2;
684 Int_t nbSegs = n * (2 * n + 3);
685 Int_t nbPols = n * (n + 2);
686 if (buffer.SetRawSizes(nbPnts, 3 * nbPnts, nbSegs, 3 * nbSegs, nbPols, 2 * n * 5 + n * n * 6)) {
687 buffer.SetSectionsValid(TBuffer3D::kRawSizes);
688 }
689 }
690 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
691 SetPoints(buffer.fPnts);
692 if (!buffer.fLocalFrame) {
693 TransformPoints(buffer.fPnts, buffer.NbPnts());
694 }
695 SetSegsAndPols(buffer);
696 buffer.SetSectionsValid(TBuffer3D::kRaw);
697 }
698
699 return buffer;
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Check the inside status for each of the points in the array.
704/// Input: Array of point coordinates + vector size
705/// Output: Array of Booleans for the inside of each point
706
707void TGeoParaboloid::Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const
708{
709 for (Int_t i = 0; i < vecsize; i++)
710 inside[i] = Contains(&points[3 * i]);
711}
712
713////////////////////////////////////////////////////////////////////////////////
714/// Compute the normal for an array o points so that norm.dot.dir is positive
715/// Input: Arrays of point coordinates and directions + vector size
716/// Output: Array of normal directions
717
718void TGeoParaboloid::ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize)
719{
720 for (Int_t i = 0; i < vecsize; i++)
721 ComputeNormal(&points[3 * i], &dirs[3 * i], &norms[3 * i]);
722}
723
724////////////////////////////////////////////////////////////////////////////////
725/// Compute distance from array of input points having directions specified by dirs. Store output in dists
726
727void TGeoParaboloid::DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
728 Double_t *step) const
729{
730 for (Int_t i = 0; i < vecsize; i++)
731 dists[i] = DistFromInside(&points[3 * i], &dirs[3 * i], 3, step[i]);
732}
733
734////////////////////////////////////////////////////////////////////////////////
735/// Compute distance from array of input points having directions specified by dirs. Store output in dists
736
737void TGeoParaboloid::DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize,
738 Double_t *step) const
739{
740 for (Int_t i = 0; i < vecsize; i++)
741 dists[i] = DistFromOutside(&points[3 * i], &dirs[3 * i], 3, step[i]);
742}
743
744////////////////////////////////////////////////////////////////////////////////
745/// Compute safe distance from each of the points in the input array.
746/// Input: Array of point coordinates, array of statuses for these points, size of the arrays
747/// Output: Safety values
748
749void TGeoParaboloid::Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const
750{
751 for (Int_t i = 0; i < vecsize; i++)
752 safe[i] = Safety(&points[3 * i], inside[i]);
753}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
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 Rectangle_t WindowAttributes_t Float_t r
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
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
Int_t * fSegs
Definition TBuffer3D.h:114
void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections, Bool_t localFrame) const override
Fill the supplied buffer, with sections in desired frame See TBuffer3D.h for explanation of sections,...
Double_t fDX
Definition TGeoBBox.h:20
void InspectShape() const override
Double_t fDY
Definition TGeoBBox.h:21
Double_t fDZ
Definition TGeoBBox.h:22
Int_t GetNsegments() const
Get number of segments approximating circles.
Geometrical transformation package.
Definition TGeoMatrix.h:38
const TBuffer3D & GetBuffer3D(Int_t reqSections, Bool_t localFrame) const override
Stub implementation to avoid forcing implementation at this stage.
Double_t Capacity() const override
Double_t DistToParaboloid(const Double_t *point, const Double_t *dir, Bool_t in) const
void InspectShape() const override
TGeoShape * GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const override
TBuffer3D * MakeBuffer3D() const override
void GetMeshNumbers(Int_t &nvert, Int_t &nsegs, Int_t &npols) const override
Int_t GetNmeshVertices() const override
void Sizeof3D() const override
void SetPoints(Double_t *points) const override
void GetBoundingCylinder(Double_t *param) const override
void ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm) override
void SetParaboloidDimensions(Double_t rlo, Double_t rhi, Double_t dz)
void SetSegsAndPols(TBuffer3D &buff) 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
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Computes distance from point (px,py) to the object.
void ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize) override
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save a primitive as a C++ statement(s) on output stream "out".
~TGeoParaboloid() override
void DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
void DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const override
void ComputeBBox() 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
Bool_t Contains(const Double_t *point) const override
void Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const override
void SetDimensions(Double_t *param) override
TGeoVolume * Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step) override
Double_t Safety(const Double_t *point, Bool_t in=kTRUE) const override
void Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const override
Base abstract class for all shapes.
Definition TGeoShape.h:25
static Double_t Big()
Definition TGeoShape.h:87
Int_t GetBasicColor() const
Get the basic color (0-7).
void TransformPoints(Double_t *points, UInt_t NbPoints) const
Tranform a set of points (LocalToMaster)
void SetShapeBit(UInt_t f, Bool_t set)
Equivalent of TObject::SetBit.
const char * GetPointerName() const
Provide a pointer name containing uid.
Int_t ShapeDistancetoPrimitive(Int_t numpoints, Int_t px, Int_t py) const
Returns distance to shape primitive mesh.
const char * GetName() const override
Get the shape name.
@ kGeoSavePrimitive
Definition TGeoShape.h:64
@ kGeoParaboloid
Definition TGeoShape.h:61
@ kGeoRunTimeShape
Definition TGeoShape.h:40
static Double_t Tolerance()
Definition TGeoShape.h:90
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
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:225
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:798
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1005
const Int_t n
Definition legend1.C:16
double dist(Rotation3D const &r1, Rotation3D const &r2)
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
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:650
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:666
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:598
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:592
Double_t Tan(Double_t)
Returns the tangent of an angle of x radians.
Definition TMath.h:604
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
auto * tt
Definition textangle.C:16
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345