Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Ellipsoid.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// class G4Ellipsoid
27//
28// Implementation of G4Ellipsoid class
29//
30// 10.11.99 G.Horton-Smith: first writing, based on G4Sphere class
31// 25.02.05 G.Guerrieri: Revised
32// 15.12.19 E.Tcherniaev: Complete revision
33// --------------------------------------------------------------------
34
35#include "G4Ellipsoid.hh"
36
37#if !(defined(G4GEOM_USE_UELLIPSOID) && defined(G4GEOM_USE_SYS_USOLIDS))
38
39#include "globals.hh"
40
41#include "G4VoxelLimits.hh"
42#include "G4AffineTransform.hh"
44#include "G4BoundingEnvelope.hh"
45#include "G4QuickRand.hh"
46
48
49#include "G4VGraphicsScene.hh"
50#include "G4VisExtent.hh"
51
52#include "G4AutoLock.hh"
53
54namespace
55{
56 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
57 G4Mutex lateralareaMutex = G4MUTEX_INITIALIZER;
58 G4Mutex ellipseMutex = G4MUTEX_INITIALIZER;
59}
60
61using namespace CLHEP;
62
63//////////////////////////////////////////////////////////////////////////
64//
65// Constructor
66
68 G4double xSemiAxis,
69 G4double ySemiAxis,
70 G4double zSemiAxis,
71 G4double zBottomCut,
72 G4double zTopCut)
73 : G4VSolid(name), fDx(xSemiAxis), fDy(ySemiAxis), fDz(zSemiAxis),
74 fZBottomCut(zBottomCut), fZTopCut(zTopCut)
75{
76 CheckParameters();
77}
78
79//////////////////////////////////////////////////////////////////////////
80//
81// Fake default constructor - sets only member data and allocates memory
82// for usage restricted to object persistency
83
85 : G4VSolid(a), fDx(0.), fDy(0.), fDz(0.), fZBottomCut(0.), fZTopCut(0.)
86{
87}
88
89//////////////////////////////////////////////////////////////////////////
90//
91// Destructor
92
94{
95 delete fpPolyhedron; fpPolyhedron = nullptr;
96}
97
98//////////////////////////////////////////////////////////////////////////
99//
100// Copy constructor
101
103 : G4VSolid(rhs),
104 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz),
105 fZBottomCut(rhs.fZBottomCut), fZTopCut(rhs.fZTopCut),
106 halfTolerance(rhs.halfTolerance),
107 fXmax(rhs.fXmax), fYmax(rhs.fYmax),
108 fRsph(rhs.fRsph), fR(rhs.fR),
109 fSx(rhs.fSx), fSy(rhs.fSy), fSz(rhs.fSz),
110 fZMidCut(rhs.fZMidCut), fZDimCut(rhs.fZDimCut),
111 fQ1(rhs.fQ1), fQ2(rhs.fQ2),
112 fCubicVolume(rhs.fCubicVolume),
113 fSurfaceArea(rhs.fSurfaceArea),
114 fLateralArea(rhs.fLateralArea)
115{
116}
117
118//////////////////////////////////////////////////////////////////////////
119//
120// Assignment operator
121
123{
124 // Check assignment to self
125 //
126 if (this == &rhs) { return *this; }
127
128 // Copy base class data
129 //
131
132 // Copy data
133 //
134 fDx = rhs.fDx;
135 fDy = rhs.fDy;
136 fDz = rhs.fDz;
137 fZBottomCut = rhs.fZBottomCut;
138 fZTopCut = rhs.fZTopCut;
139
140 halfTolerance = rhs.halfTolerance;
141 fXmax = rhs.fXmax;
142 fYmax = rhs.fYmax;
143 fRsph = rhs.fRsph;
144 fR = rhs.fR;
145 fSx = rhs.fSx;
146 fSy = rhs.fSy;
147 fSz = rhs.fSz;
148 fZMidCut = rhs.fZMidCut;
149 fZDimCut = rhs.fZDimCut;
150 fQ1 = rhs.fQ1;
151 fQ2 = rhs.fQ2;
152
153 fCubicVolume = rhs.fCubicVolume;
154 fSurfaceArea = rhs.fSurfaceArea;
155 fLateralArea = rhs.fLateralArea;
156
157 fRebuildPolyhedron = false;
158 delete fpPolyhedron; fpPolyhedron = nullptr;
159
160 return *this;
161}
162
163//////////////////////////////////////////////////////////////////////////
164//
165// Check parameters and make precalculation
166
167void G4Ellipsoid::CheckParameters()
168{
169 halfTolerance = 0.5 * kCarTolerance; // half tolerance
170 G4double dmin = 2. * kCarTolerance;
171
172 // Check dimensions
173 //
174 if (fDx < dmin || fDy < dmin || fDz < dmin)
175 {
176 std::ostringstream message;
177 message << "Invalid (too small or negative) dimensions for Solid: "
178 << GetName() << "\n"
179 << " semi-axis x: " << fDx << "\n"
180 << " semi-axis y: " << fDy << "\n"
181 << " semi-axis z: " << fDz;
182 G4Exception("G4Ellipsoid::CheckParameters()", "GeomSolids0002",
183 FatalException, message);
184 }
185 G4double A = fDx;
186 G4double B = fDy;
187 G4double C = fDz;
188
189 // Check cuts
190 //
191 if (fZBottomCut == 0. && fZTopCut == 0.)
192 {
193 fZBottomCut = -C;
194 fZTopCut = C;
195 }
196 if (fZBottomCut >= C || fZTopCut <= -C || fZBottomCut >= fZTopCut)
197 {
198 std::ostringstream message;
199 message << "Invalid Z cuts for Solid: "
200 << GetName() << "\n"
201 << " bottom cut: " << fZBottomCut << "\n"
202 << " top cut: " << fZTopCut;
203 G4Exception("G4Ellipsoid::CheckParameters()", "GeomSolids0002",
204 FatalException, message);
205
206 }
207 fZBottomCut = std::max(fZBottomCut, -C);
208 fZTopCut = std::min(fZTopCut, C);
209
210 // Set extent in x and y
211 fXmax = A;
212 fYmax = B;
213 if (fZBottomCut > 0.)
214 {
215 G4double ratio = fZBottomCut / C;
216 G4double scale = std::sqrt((1. - ratio) * (1 + ratio));
217 fXmax *= scale;
218 fYmax *= scale;
219 }
220 if (fZTopCut < 0.)
221 {
222 G4double ratio = fZTopCut / C;
223 G4double scale = std::sqrt((1. - ratio) * (1 + ratio));
224 fXmax *= scale;
225 fYmax *= scale;
226 }
227
228 // Set scale factors
229 fRsph = std::max(std::max(A, B), C); // bounding sphere
230 fR = std::min(std::min(A, B), C); // radius of sphere after scaling
231 fSx = fR / A; // X scale factor
232 fSy = fR / B; // Y scale factor
233 fSz = fR / C; // Z scale factor
234
235 // Scaled cuts
236 fZMidCut = 0.5 * (fZTopCut + fZBottomCut) * fSz; // middle position
237 fZDimCut = 0.5 * (fZTopCut - fZBottomCut) * fSz; // half distance
238
239 // Coefficients for approximation of distance: Q1 * (x^2 + y^2 + z^2) - Q2
240 fQ1 = 0.5 / fR;
241 fQ2 = 0.5 * fR + halfTolerance * halfTolerance * fQ1;
242
243 fCubicVolume = 0.; // volume
244 fSurfaceArea = 0.; // surface area
245 fLateralArea = 0.; // lateral surface area
246}
247
248//////////////////////////////////////////////////////////////////////////
249//
250// Dispatch to parameterisation for replication mechanism dimension
251// computation & modification
252
254 const G4int n,
255 const G4VPhysicalVolume* pRep)
256{
257 p->ComputeDimensions(*this,n,pRep);
258}
259
260//////////////////////////////////////////////////////////////////////////
261//
262// Get bounding box
263
265 G4ThreeVector& pMax) const
266{
267 pMin.set(-fXmax,-fYmax, fZBottomCut);
268 pMax.set( fXmax, fYmax, fZTopCut);
269}
270
271//////////////////////////////////////////////////////////////////////////
272//
273// Calculate extent under transform and specified limits
274
275G4bool
277 const G4VoxelLimits& pVoxelLimit,
278 const G4AffineTransform& pTransform,
279 G4double& pMin, G4double& pMax) const
280{
281 G4ThreeVector bmin, bmax;
282
283 // Get bounding box
284 BoundingLimits(bmin,bmax);
285
286 // Find extent
287 G4BoundingEnvelope bbox(bmin,bmax);
288 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
289}
290
291//////////////////////////////////////////////////////////////////////////
292//
293// Return position of point: inside/outside/on surface
294
296{
297 G4double x = p.x() * fSx;
298 G4double y = p.y() * fSy;
299 G4double z = p.z() * fSz;
300 G4double rr = x * x + y * y + z * z;
301 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
302 G4double distR = fQ1 * rr - fQ2;
303 G4double dist = std::max(distZ, distR);
304
305 if (dist > halfTolerance) { return kOutside; }
306 return (dist > -halfTolerance) ? kSurface : kInside;
307}
308
309//////////////////////////////////////////////////////////////////////////
310//
311// Return unit normal to surface at p
312
314{
315 G4ThreeVector norm(0., 0., 0.);
316 G4int nsurf = 0;
317
318 // Check cuts
319 G4double x = p.x() * fSx;
320 G4double y = p.y() * fSy;
321 G4double z = p.z() * fSz;
322 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
323 if (std::abs(distZ) <= halfTolerance)
324 {
325 norm.setZ(std::copysign(1., z - fZMidCut));
326 ++nsurf;
327 }
328
329 // Check lateral surface
330 G4double distR = fQ1*(x*x + y*y + z*z) - fQ2;
331 if (std::abs(distR) <= halfTolerance)
332 {
333 // normal = (p.x/a^2, p.y/b^2, p.z/c^2)
334 norm += G4ThreeVector(x*fSx, y*fSy, z*fSz).unit();
335 ++nsurf;
336 }
337
338 // Return normal
339 if (nsurf == 1)
340 {
341 return norm;
342 }
343 if (nsurf > 1)
344 {
345 return norm.unit(); // edge
346 }
347
348#ifdef G4SPECSDEBUG
349 std::ostringstream message;
350 G4long oldprc = message.precision(16);
351 message << "Point p is not on surface (!?) of solid: "
352 << GetName() << "\n";
353 message << "Position:\n";
354 message << " p.x() = " << p.x()/mm << " mm\n";
355 message << " p.y() = " << p.y()/mm << " mm\n";
356 message << " p.z() = " << p.z()/mm << " mm";
357 G4cout.precision(oldprc);
358 G4Exception("G4Ellipsoid::SurfaceNormal(p)", "GeomSolids1002",
359 JustWarning, message );
360 DumpInfo();
361#endif
362 return ApproxSurfaceNormal(p);
363}
364
365//////////////////////////////////////////////////////////////////////////
366//
367// Find surface nearest to point and return corresponding normal.
368// This method normally should not be called.
369
370G4ThreeVector G4Ellipsoid::ApproxSurfaceNormal(const G4ThreeVector& p) const
371{
372 G4double x = p.x() * fSx;
373 G4double y = p.y() * fSy;
374 G4double z = p.z() * fSz;
375 G4double rr = x * x + y * y + z * z;
376 G4double distZ = std::abs(z - fZMidCut) - fZDimCut;
377 G4double distR = std::sqrt(rr) - fR;
378 if (distR > distZ && rr > 0.)
379 { // distR > distZ is correct!
380 return G4ThreeVector(x*fSx, y*fSy, z*fSz).unit();
381 }
382 return { 0., 0., std::copysign(1., z - fZMidCut) };
383}
384
385//////////////////////////////////////////////////////////////////////////
386//
387// Calculate distance to shape from outside along normalised vector
388
390 const G4ThreeVector& v) const
391{
392 G4double offset = 0.;
393 G4ThreeVector pcur = p;
394
395 // Check if point is flying away, relative to bounding box
396 //
397 G4double safex = std::abs(p.x()) - fXmax;
398 G4double safey = std::abs(p.y()) - fYmax;
399 G4double safet = p.z() - fZTopCut;
400 G4double safeb = fZBottomCut - p.z();
401
402 if (safex >= -halfTolerance && p.x() * v.x() >= 0.) { return kInfinity; }
403 if (safey >= -halfTolerance && p.y() * v.y() >= 0.) { return kInfinity; }
404 if (safet >= -halfTolerance && v.z() >= 0.) { return kInfinity; }
405 if (safeb >= -halfTolerance && v.z() <= 0.) { return kInfinity; }
406
407 // Relocate point, if required
408 //
409 G4double safe = std::max(std::max(std::max(safex, safey), safet), safeb);
410 if (safe > 32. * fRsph)
411 {
412 offset = (1. - 1.e-08) * safe - 2. * fRsph;
413 pcur += offset * v;
414 G4double dist = DistanceToIn(pcur, v);
415 return (dist == kInfinity) ? kInfinity : dist + offset;
416 }
417
418 // Scale ellipsoid to sphere
419 //
420 G4double px = pcur.x() * fSx;
421 G4double py = pcur.y() * fSy;
422 G4double pz = pcur.z() * fSz;
423 G4double vx = v.x() * fSx;
424 G4double vy = v.y() * fSy;
425 G4double vz = v.z() * fSz;
426
427 // Check if point is leaving the solid
428 //
429 G4double dzcut = fZDimCut;
430 G4double pzcut = pz - fZMidCut;
431 G4double distZ = std::abs(pzcut) - dzcut;
432 if (distZ >= -halfTolerance && pzcut * vz >= 0.) { return kInfinity; }
433
434 G4double rr = px * px + py * py + pz * pz;
435 G4double pv = px * vx + py * vy + pz * vz;
436 G4double distR = fQ1 * rr - fQ2;
437 if (distR >= -halfTolerance && pv >= 0.) { return kInfinity; }
438
439 G4double A = vx * vx + vy * vy + vz * vz;
440 G4double B = pv;
441 G4double C = rr - fR * fR;
442 G4double D = B * B - A * C;
443 // scratch^2 = R^2 - (R - halfTolerance)^2 = 2 * R * halfTolerance
444 G4double EPS = A * A * fR * kCarTolerance; // discriminant at scratching
445 if (D <= EPS) { return kInfinity; } // no intersection or scratching
446
447 // Find intersection with Z planes
448 //
449 G4double invz = (vz == 0) ? DBL_MAX : -1./vz;
450 G4double dz = std::copysign(dzcut, invz);
451 G4double tzmin = (pzcut - dz) * invz;
452 G4double tzmax = (pzcut + dz) * invz;
453
454 // Find intersection with lateral surface
455 //
456 G4double tmp = -B - std::copysign(std::sqrt(D), B);
457 G4double t1 = tmp / A;
458 G4double t2 = C / tmp;
459 G4double trmin = std::min(t1, t2);
460 G4double trmax = std::max(t1, t2);
461
462 // Return distance
463 //
464 G4double tmin = std::max(tzmin, trmin);
465 G4double tmax = std::min(tzmax, trmax);
466
467 if (tmax - tmin <= halfTolerance) { return kInfinity; } // touch or no hit
468
469 return (tmin < halfTolerance) ? offset : tmin + offset;
470}
471
472//////////////////////////////////////////////////////////////////////////
473//
474// Estimate distance to surface from outside
475
477{
478 G4double px = p.x();
479 G4double py = p.y();
480 G4double pz = p.z();
481
482 // Safety distance to bounding box
483 G4double distX = std::abs(px) - fXmax;
484 G4double distY = std::abs(py) - fYmax;
485 G4double distZ = std::max(pz - fZTopCut, fZBottomCut - pz);
486 G4double distB = std::max(std::max(distX, distY), distZ);
487
488 // Safety distance to lateral surface
489 G4double x = px * fSx;
490 G4double y = py * fSy;
491 G4double z = pz * fSz;
492 G4double distR = std::sqrt(x*x + y*y + z*z) - fR;
493
494 // Return safety to in
495 G4double dist = std::max(distB, distR);
496 return (dist < 0.) ? 0. : dist;
497}
498
499//////////////////////////////////////////////////////////////////////////
500//
501// Calculate distance to surface from inside along normalised vector
502
504 const G4ThreeVector& v,
505 const G4bool calcNorm,
506 G4bool* validNorm,
507 G4ThreeVector* n ) const
508{
509 // Check if point flying away relative to Z planes
510 //
511 G4double pz = p.z() * fSz;
512 G4double vz = v.z() * fSz;
513 G4double dzcut = fZDimCut;
514 G4double pzcut = pz - fZMidCut;
515 G4double distZ = std::abs(pzcut) - dzcut;
516 if (distZ >= -halfTolerance && pzcut * vz > 0.)
517 {
518 if (calcNorm)
519 {
520 *validNorm = true;
521 n->set(0., 0., std::copysign(1., pzcut));
522 }
523 return 0.;
524 }
525
526 // Check if point is flying away relative to lateral surface
527 //
528 G4double px = p.x() * fSx;
529 G4double py = p.y() * fSy;
530 G4double vx = v.x() * fSx;
531 G4double vy = v.y() * fSy;
532 G4double rr = px * px + py * py + pz * pz;
533 G4double pv = px * vx + py * vy + pz * vz;
534 G4double distR = fQ1 * rr - fQ2;
535 if (distR >= -halfTolerance && pv > 0.)
536 {
537 if (calcNorm)
538 {
539 *validNorm = true;
540 *n = G4ThreeVector(px*fSx, py*fSy, pz*fSz).unit();
541 }
542 return 0.;
543 }
544
545 // Just in case check if point is outside (normally it should never be)
546 //
547 if (std::max(distZ, distR) > halfTolerance)
548 {
549#ifdef G4SPECSDEBUG
550 std::ostringstream message;
551 G4long oldprc = message.precision(16);
552 message << "Point p is outside (!?) of solid: "
553 << GetName() << G4endl;
554 message << "Position: " << p << G4endl;;
555 message << "Direction: " << v;
556 G4cout.precision(oldprc);
557 G4Exception("G4Ellipsoid::DistanceToOut(p,v)", "GeomSolids1002",
558 JustWarning, message );
559 DumpInfo();
560#endif
561 if (calcNorm)
562 {
563 *validNorm = true;
564 *n = ApproxSurfaceNormal(p);
565 }
566 return 0.;
567 }
568
569 // Set coefficients of quadratic equation: A t^2 + 2B t + C = 0
570 //
571 G4double A = vx * vx + vy * vy + vz * vz;
572 G4double B = pv;
573 G4double C = rr - fR * fR;
574 G4double D = B * B - A * C;
575 // It is expected that the point is located inside the sphere, so
576 // max term in the expression for discriminant is A * R^2 and
577 // max calculation error can be derived as follows:
578 // A * (1 + 2e) * R^2 * (1 + 2e) = A * R^2 + (4 * A * R^2 * e)
579 G4double EPS = 4. * A * fR * fR * DBL_EPSILON; // calculation error
580
581 if (D <= EPS) // no intersection
582 {
583 if (calcNorm)
584 {
585 *validNorm = true;
586 *n = G4ThreeVector(px*fSx, py*fSy, pz*fSz).unit();
587 }
588 return 0.;
589 }
590
591 // Find intersection with Z cuts
592 //
593 G4double tzmax = (vz == 0.)
594 ? DBL_MAX
595 : (std::copysign(dzcut, vz) - pzcut) / vz;
596
597 // Find intersection with lateral surface
598 //
599 G4double tmp = -B - std::copysign(std::sqrt(D), B);
600 G4double trmax = (tmp < 0.) ? C/tmp : tmp/A;
601
602 // Find distance and set normal, if required
603 //
604 G4double tmax = std::min(tzmax, trmax);
605 //if (tmax < halfTolerance) tmax = 0.;
606
607 if (calcNorm)
608 {
609 *validNorm = true;
610 if (tmax == tzmax)
611 {
612 G4double pznew = pz + tmax * vz;
613 n->set(0., 0., (pznew > fZMidCut) ? 1. : -1.);
614 }
615 else
616 {
617 G4double nx = (px + tmax * vx) * fSx;
618 G4double ny = (py + tmax * vy) * fSy;
619 G4double nz = (pz + tmax * vz) * fSz;
620 *n = G4ThreeVector(nx, ny, nz).unit();
621 }
622 }
623 return tmax;
624}
625
626//////////////////////////////////////////////////////////////////////////
627//
628// Estimate distance to surface from inside
629
631{
632 // Safety distance in z direction
633 G4double distZ = std::min(fZTopCut - p.z(), p.z() - fZBottomCut);
634
635 // Safety distance to lateral surface
636 G4double x = p.x() * fSx;
637 G4double y = p.y() * fSy;
638 G4double z = p.z() * fSz;
639 G4double distR = fR - std::sqrt(x*x + y*y + z*z);
640
641 // Return safety to out
642 G4double dist = std::min(distZ, distR);
643 return (dist < 0.) ? 0. : dist;
644}
645
646//////////////////////////////////////////////////////////////////////////
647//
648// Return entity type
649
651{
652 return {"G4Ellipsoid"};
653}
654
655//////////////////////////////////////////////////////////////////////////
656//
657// Make a clone of the object
658
660{
661 return new G4Ellipsoid(*this);
662}
663
664//////////////////////////////////////////////////////////////////////////
665//
666// Stream object contents to output stream
667
668std::ostream& G4Ellipsoid::StreamInfo( std::ostream& os ) const
669{
670 G4long oldprc = os.precision(16);
671 os << "-----------------------------------------------------------\n"
672 << " *** Dump for solid - " << GetName() << " ***\n"
673 << " ===================================================\n"
674 << " Solid type: " << GetEntityType() << "\n"
675 << " Parameters: \n"
676 << " semi-axis x: " << GetDx()/mm << " mm \n"
677 << " semi-axis y: " << GetDy()/mm << " mm \n"
678 << " semi-axis z: " << GetDz()/mm << " mm \n"
679 << " lower cut in z: " << GetZBottomCut()/mm << " mm \n"
680 << " upper cut in z: " << GetZTopCut()/mm << " mm \n"
681 << "-----------------------------------------------------------\n";
682 os.precision(oldprc);
683 return os;
684}
685
686//////////////////////////////////////////////////////////////////////////
687//
688// Calculate area of lateral surface
689
690G4double G4Ellipsoid::LateralSurfaceArea() const
691{
692 constexpr G4int NPHI = 1000.;
693 constexpr G4double dPhi = CLHEP::halfpi/NPHI;
694 constexpr G4double eps = 4.*DBL_EPSILON;
695
696 G4double aa = fDx*fDx;
697 G4double bb = fDy*fDy;
698 G4double cc = fDz*fDz;
699 G4double ab = fDx*fDy;
700 G4double cc_aa = cc/aa;
701 G4double cc_bb = cc/bb;
702 G4double zmax = std::min(fZTopCut, fDz);
703 G4double zmin = std::max(fZBottomCut,-fDz);
704 G4double zmax_c = zmax/fDz;
705 G4double zmin_c = zmin/fDz;
706 G4double area = 0.;
707
708 if (aa == bb) // spheroid, use analytical expression
709 {
710 G4double k = fDz/fDx;
711 G4double kk = k*k;
712 if (kk < 1. - eps)
713 {
714 G4double invk = fDx/fDz;
715 G4double root = std::sqrt(1. - kk);
716 G4double tmax = zmax_c*root;
717 G4double tmin = zmin_c*root;
718 area = CLHEP::pi*ab*
719 ((zmax_c*std::sqrt(kk + tmax*tmax) - zmin_c*std::sqrt(kk + tmin*tmin)) +
720 (std::asinh(tmax*invk) - std::asinh(tmin*invk))*kk/root);
721 }
722 else if (kk > 1. + eps)
723 {
724 G4double invk = fDx/fDz;
725 G4double root = std::sqrt(kk - 1.);
726 G4double tmax = zmax_c*root;
727 G4double tmin = zmin_c*root;
728 area = CLHEP::pi*ab*
729 ((zmax_c*std::sqrt(kk - tmax*tmax) - zmin_c*std::sqrt(kk - tmin*tmin)) +
730 (std::asin(tmax*invk) - std::asin(tmin*invk))*kk/root);
731 }
732 else
733 {
734 area = CLHEP::twopi*fDx*(zmax - zmin);
735 }
736 return area;
737 }
738
739 // ellipsoid, integration along phi
740 for (G4int i = 0; i < NPHI; ++i)
741 {
742 G4double sinPhi = std::sin(dPhi*(i + 0.5));
743 G4double kk = cc_aa + (cc_bb - cc_aa)*sinPhi*sinPhi;
744 if (kk < 1. - eps)
745 {
746 G4double root = std::sqrt(1. - kk);
747 G4double tmax = zmax_c*root;
748 G4double tmin = zmin_c*root;
749 G4double invk = 1./std::sqrt(kk);
750 area += 2.*ab*dPhi*
751 ((zmax_c*std::sqrt(kk + tmax*tmax) - zmin_c*std::sqrt(kk + tmin*tmin)) +
752 (std::asinh(tmax*invk) - std::asinh(tmin*invk))*kk/root);
753 }
754 else if (kk > 1. + eps)
755 {
756 G4double root = std::sqrt(kk - 1.);
757 G4double tmax = zmax_c*root;
758 G4double tmin = zmin_c*root;
759 G4double invk = 1./std::sqrt(kk);
760 area += 2.*ab*dPhi*
761 ((zmax_c*std::sqrt(kk - tmax*tmax) - zmin_c*std::sqrt(kk - tmin*tmin)) +
762 (std::asin(tmax*invk) - std::asin(tmin*invk))*kk/root);
763 }
764 else
765 {
766 area += 4.*ab*dPhi*(zmax_c - zmin_c);
767 }
768 }
769 return area;
770}
771
772//////////////////////////////////////////////////////////////////////////
773//
774// Return volume
775
777{
778 if (fCubicVolume == 0)
779 {
780 G4AutoLock l(&ellipseMutex);
781 G4double piAB_3 = CLHEP::pi * fDx * fDy / 3.;
782 fCubicVolume = 4. * piAB_3 * fDz;
783 if (fZBottomCut > -fDz)
784 {
785 G4double hbot = 1. + fZBottomCut / fDz;
786 fCubicVolume -= piAB_3 * hbot * hbot * (2. * fDz - fZBottomCut);
787 }
788 if (fZTopCut < fDz)
789 {
790 G4double htop = 1. - fZTopCut / fDz;
791 fCubicVolume -= piAB_3 * htop * htop * (2. * fDz + fZTopCut);
792 }
793 l.unlock();
794 }
795 return fCubicVolume;
796}
797
798//////////////////////////////////////////////////////////////////////////
799//
800// Return surface area
801
803{
804 if (fSurfaceArea == 0.)
805 {
806 G4AutoLock l(&ellipseMutex);
807 G4double piAB = CLHEP::pi * fDx * fDy;
808 fSurfaceArea = LateralSurfaceArea();
809 if (fZBottomCut > -fDz)
810 {
811 G4double hbot = 1. + fZBottomCut / fDz;
812 fSurfaceArea += piAB * hbot * (2. - hbot);
813 }
814 if (fZTopCut < fDz)
815 {
816 G4double htop = 1. - fZTopCut / fDz;
817 fSurfaceArea += piAB * htop * (2. - htop);
818 }
819 l.unlock();
820 }
821 return fSurfaceArea;
822}
823
824//////////////////////////////////////////////////////////////////////////
825//
826// Return random point on surface
827
829{
830 G4double A = GetDx();
831 G4double B = GetDy();
832 G4double C = GetDz();
833 G4double Zbot = GetZBottomCut();
834 G4double Ztop = GetZTopCut();
835
836 // Calculate cut areas
837 G4double invC = 1. / C;
838 G4double Hbot = 1. + Zbot * invC;
839 G4double Htop = 1. - Ztop * invC;
840 G4double piAB = CLHEP::pi * A * B;
841 G4double Sbot = piAB * Hbot * (2. - Hbot);
842 G4double Stop = piAB * Htop * (2. - Htop);
843
844 // Get area of lateral surface
845 if (fLateralArea == 0.)
846 {
847 G4AutoLock l(&lateralareaMutex);
848 fLateralArea = LateralSurfaceArea();
849 l.unlock();
850 }
851 G4double Slat = fLateralArea;
852
853 // Select surface (0 - bottom cut, 1 - lateral surface, 2 - top cut)
854 G4double select = (Sbot + Slat + Stop) * G4QuickRand();
855 G4int k = 0;
856 k += (G4int)(select > Sbot);
857 k += (G4int)(select > Sbot + Slat);
858
859 // Pick random point on selected surface
860 G4double phi = CLHEP::twopi * G4QuickRand();
861 G4double cosphi = std::cos(phi);
862 G4double sinphi = std::sin(phi);
864 switch (k)
865 {
866 case 0: // bootom z-cut, ellipse
867 {
868 G4double scale = std::sqrt(Hbot * (2. - Hbot));
869 G4double rho = scale*std::sqrt(G4QuickRand());
870 p.set(A * rho * cosphi, B * rho * sinphi, Zbot);
871 break;
872 }
873 case 1: // lateral surface (rejection sampling)
874 {
875 G4double x, y, z;
876 G4double s_max = std::max(std::max(A * B, A * C), B * C);
877 for (G4int i = 0; i < 10000; ++i)
878 {
879 // generate random point on unit sphere
880 z = (Zbot + (Ztop - Zbot) * G4QuickRand()) * invC;
881 G4double rho = std::sqrt((1. + z) * (1. - z));
882 x = rho * cosphi;
883 y = rho * sinphi;
884 // check acceptance
885 G4double ss = sqr(B * C * x) + sqr(A * C * y) + sqr(A * B * z);
886 if (sqr(s_max * G4QuickRand()) <= ss) { break; }
887 }
888 p.set(A * x, B * y, C * z);
889 break;
890 }
891 case 2: // top z-cut, ellipse
892 {
893 G4double scale = std::sqrt(Htop * (2. - Htop));
894 G4double rho = scale*std::sqrt(G4QuickRand());
895 p.set(A * rho * cosphi, B * rho * sinphi, Ztop);
896 break;
897 }
898 }
899 return p;
900}
901
902//////////////////////////////////////////////////////////////////////////
903//
904// Methods for visualisation
905
907{
908 scene.AddSolid(*this);
909}
910
911//////////////////////////////////////////////////////////////////////////
912//
913// Return vis extent
914
916{
917 return { -fXmax, fXmax, -fYmax, fYmax, fZBottomCut, fZTopCut };
918}
919
920//////////////////////////////////////////////////////////////////////////
921//
922// Create polyhedron
923
925{
926 return new G4PolyhedronEllipsoid(fDx, fDy, fDz, fZBottomCut, fZTopCut);
927}
928
929//////////////////////////////////////////////////////////////////////////
930//
931// Return pointer to polyhedron
932
934{
935 if (fpPolyhedron == nullptr ||
936 fRebuildPolyhedron ||
937 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
938 fpPolyhedron->GetNumberOfRotationSteps())
939 {
940 G4AutoLock l(&polyhedronMutex);
941 delete fpPolyhedron;
942 fpPolyhedron = CreatePolyhedron();
943 fRebuildPolyhedron = false;
944 l.unlock();
945 }
946 return fpPolyhedron;
947}
948
949#endif // !defined(G4GEOM_USE_UELLIPSOID) || !defined(G4GEOM_USE_SYS_USOLIDS)
G4TemplateAutoLock< G4Mutex > G4AutoLock
const G4double kCarTolerance
G4double C(G4double temp)
G4double B(G4double temperature)
G4double D(G4double temp)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ThreadLocal T * G4GeomSplitter< T >::offset
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
long G4long
Definition G4Types.hh:87
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
void setZ(double)
void set(double x, double y, double z)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4double GetDx() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
G4Ellipsoid(const G4String &name, G4double xSemiAxis, G4double ySemiAxis, G4double zSemiAxis, G4double zBottomCut=0., G4double zTopCut=0.)
EInside Inside(const G4ThreeVector &p) const override
G4Polyhedron * CreatePolyhedron() const override
~G4Ellipsoid() override
G4double GetDy() const
G4double GetSurfaceArea() override
G4VSolid * Clone() const override
G4ThreeVector GetPointOnSurface() const override
G4double GetZTopCut() const
G4double GetZBottomCut() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
G4double GetCubicVolume() override
G4GeometryType GetEntityType() const override
G4double GetDz() const
G4Polyhedron * GetPolyhedron() const override
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4VisExtent GetExtent() const override
G4Ellipsoid & operator=(const G4Ellipsoid &rhs)
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
std::ostream & StreamInfo(std::ostream &os) const override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
virtual void AddSolid(const G4Box &)=0
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
T sqr(const T &x)
Definition templates.hh:128
#define DBL_EPSILON
Definition templates.hh:66
#define DBL_MAX
Definition templates.hh:62