Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Para.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// Implementation for G4Para class
27//
28// 21.03.95 P.Kent: Modified for `tolerant' geom
29// 31.10.96 V.Grichine: Modifications according G4Box/Tubs before to commit
30// 28.04.05 V.Grichine: new SurfaceNormal according to J. Apostolakis proposal
31// 29.05.17 E.Tcherniaev: complete revision, speed-up
32////////////////////////////////////////////////////////////////////////////
33
34#include "G4Para.hh"
35
36#if !defined(G4GEOM_USE_UPARA)
37
38#include "G4VoxelLimits.hh"
39#include "G4AffineTransform.hh"
40#include "G4BoundingEnvelope.hh"
41#include "G4QuickRand.hh"
42
44
45#include "G4VGraphicsScene.hh"
46#include "G4AutoLock.hh"
47
48namespace
49{
50 G4Mutex paraMutex = G4MUTEX_INITIALIZER;
51}
52
53using namespace CLHEP;
54
55//////////////////////////////////////////////////////////////////////////
56//
57// Constructor - set & check half widths
58
60 G4double pDx, G4double pDy, G4double pDz,
61 G4double pAlpha, G4double pTheta, G4double pPhi)
62 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
63{
64 SetAllParameters(pDx, pDy, pDz, pAlpha, pTheta, pPhi);
65 fRebuildPolyhedron = false; // default value for G4CSGSolid
66}
67
68//////////////////////////////////////////////////////////////////////////
69//
70// Constructor - design of trapezoid based on 8 vertices
71
73 const G4ThreeVector pt[8] )
74 : G4CSGSolid(pName), halfCarTolerance(0.5*kCarTolerance)
75{
76 // Find dimensions and trigonometric values
77 //
78 fDx = (pt[3].x() - pt[2].x())*0.5;
79 fDy = (pt[2].y() - pt[1].y())*0.5;
80 fDz = pt[7].z();
81 CheckParameters(); // check dimensions
82
83 fTalpha = (pt[2].x() + pt[3].x() - pt[1].x() - pt[0].x())*0.25/fDy;
84 fTthetaCphi = (pt[4].x() + fDy*fTalpha + fDx)/fDz;
85 fTthetaSphi = (pt[4].y() + fDy)/fDz;
86 MakePlanes();
87
88 // Recompute vertices
89 //
90 G4ThreeVector v[8];
91 G4double DyTalpha = fDy*fTalpha;
92 G4double DzTthetaSphi = fDz*fTthetaSphi;
93 G4double DzTthetaCphi = fDz*fTthetaCphi;
94 v[0].set(-DzTthetaCphi-DyTalpha-fDx, -DzTthetaSphi-fDy, -fDz);
95 v[1].set(-DzTthetaCphi-DyTalpha+fDx, -DzTthetaSphi-fDy, -fDz);
96 v[2].set(-DzTthetaCphi+DyTalpha-fDx, -DzTthetaSphi+fDy, -fDz);
97 v[3].set(-DzTthetaCphi+DyTalpha+fDx, -DzTthetaSphi+fDy, -fDz);
98 v[4].set( DzTthetaCphi-DyTalpha-fDx, DzTthetaSphi-fDy, fDz);
99 v[5].set( DzTthetaCphi-DyTalpha+fDx, DzTthetaSphi-fDy, fDz);
100 v[6].set( DzTthetaCphi+DyTalpha-fDx, DzTthetaSphi+fDy, fDz);
101 v[7].set( DzTthetaCphi+DyTalpha+fDx, DzTthetaSphi+fDy, fDz);
102
103 // Compare with original vertices
104 //
105 for (G4int i=0; i<8; ++i)
106 {
107 G4double delx = std::abs(pt[i].x() - v[i].x());
108 G4double dely = std::abs(pt[i].y() - v[i].y());
109 G4double delz = std::abs(pt[i].z() - v[i].z());
110 G4double discrepancy = std::max(std::max(delx,dely),delz);
111 if (discrepancy > 0.1*kCarTolerance)
112 {
113 std::ostringstream message;
114 G4long oldprc = message.precision(16);
115 message << "Invalid vertice coordinates for Solid: " << GetName()
116 << "\nVertix #" << i << ", discrepancy = " << discrepancy
117 << "\n original : " << pt[i]
118 << "\n recomputed : " << v[i];
119 G4cout.precision(oldprc);
120 G4Exception("G4Para::G4Para()", "GeomSolids0002",
121 FatalException, message);
122
123 }
124 }
125}
126
127//////////////////////////////////////////////////////////////////////////
128//
129// Fake default constructor - sets only member data and allocates memory
130// for usage restricted to object persistency
131
132G4Para::G4Para( __void__& a )
133 : G4CSGSolid(a), halfCarTolerance(0.5*kCarTolerance)
134{
135 SetAllParameters(1., 1., 1., 0., 0., 0.);
136 fRebuildPolyhedron = false; // default value for G4CSGSolid
137}
138
139//////////////////////////////////////////////////////////////////////////
140//
141// Copy constructor
142
144 : G4CSGSolid(rhs), halfCarTolerance(rhs.halfCarTolerance),
145 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTalpha(rhs.fTalpha),
146 fTthetaCphi(rhs.fTthetaCphi),fTthetaSphi(rhs.fTthetaSphi)
147{
148 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
149}
150
151//////////////////////////////////////////////////////////////////////////
152//
153// Assignment operator
154
156{
157 // Check assignment to self
158 //
159 if (this == &rhs) { return *this; }
160
161 // Copy base class data
162 //
164
165 // Copy data
166 //
167 halfCarTolerance = rhs.halfCarTolerance;
168 fDx = rhs.fDx;
169 fDy = rhs.fDy;
170 fDz = rhs.fDz;
171 fTalpha = rhs.fTalpha;
172 fTthetaCphi = rhs.fTthetaCphi;
173 fTthetaSphi = rhs.fTthetaSphi;
174 for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
175
176 return *this;
177}
178
179//////////////////////////////////////////////////////////////////////////
180//
181// Set all parameters, as for constructor - set and check half-widths
182
184 G4double pAlpha, G4double pTheta, G4double pPhi)
185{
186 // Reset data of the base class
187 fCubicVolume = 0;
188 fSurfaceArea = 0;
189 fRebuildPolyhedron = true;
190
191 // Set parameters
192 fDx = pDx;
193 fDy = pDy;
194 fDz = pDz;
195 fTalpha = std::tan(pAlpha);
196 fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
197 fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
198
199 CheckParameters();
200 MakePlanes();
201}
202
203//////////////////////////////////////////////////////////////////////////
204//
205// Check dimensions
206
207void G4Para::CheckParameters()
208{
209 if (fDx < 2*kCarTolerance ||
210 fDy < 2*kCarTolerance ||
211 fDz < 2*kCarTolerance)
212 {
213 std::ostringstream message;
214 message << "Invalid (too small or negative) dimensions for Solid: "
215 << GetName()
216 << "\n X - " << fDx
217 << "\n Y - " << fDy
218 << "\n Z - " << fDz;
219 G4Exception("G4Para::CheckParameters()", "GeomSolids0002",
220 FatalException, message);
221 }
222}
223
224//////////////////////////////////////////////////////////////////////////
225//
226// Set side planes
227
228void G4Para::MakePlanes()
229{
230 G4ThreeVector vx(1, 0, 0);
231 G4ThreeVector vy(fTalpha, 1, 0);
232 G4ThreeVector vz(fTthetaCphi, fTthetaSphi, 1);
233
234 // Set -Y & +Y planes
235 //
236 G4ThreeVector ynorm = (vx.cross(vz)).unit();
237
238 fPlanes[0].a = 0.;
239 fPlanes[0].b = ynorm.y();
240 fPlanes[0].c = ynorm.z();
241 fPlanes[0].d = fPlanes[0].b*fDy; // point (0,fDy,0) is on plane
242
243 fPlanes[1].a = 0.;
244 fPlanes[1].b = -fPlanes[0].b;
245 fPlanes[1].c = -fPlanes[0].c;
246 fPlanes[1].d = fPlanes[0].d;
247
248 // Set -X & +X planes
249 //
250 G4ThreeVector xnorm = (vz.cross(vy)).unit();
251
252 fPlanes[2].a = xnorm.x();
253 fPlanes[2].b = xnorm.y();
254 fPlanes[2].c = xnorm.z();
255 fPlanes[2].d = fPlanes[2].a*fDx; // point (fDx,0,0) is on plane
256
257 fPlanes[3].a = -fPlanes[2].a;
258 fPlanes[3].b = -fPlanes[2].b;
259 fPlanes[3].c = -fPlanes[2].c;
260 fPlanes[3].d = fPlanes[2].d;
261}
262
263//////////////////////////////////////////////////////////////////////////
264//
265// Dispatch to parameterisation for replication mechanism dimension
266// computation & modification
267
269 const G4int n,
270 const G4VPhysicalVolume* pRep )
271{
272 p->ComputeDimensions(*this,n,pRep);
273}
274
275//////////////////////////////////////////////////////////////////////////
276//
277// Get bounding box
278
280{
284
285 G4double x0 = dz*fTthetaCphi;
286 G4double x1 = dy*GetTanAlpha();
287 G4double xmin =
288 std::min(
289 std::min(
290 std::min(-x0-x1-dx,-x0+x1-dx),x0-x1-dx),x0+x1-dx);
291 G4double xmax =
292 std::max(
293 std::max(
294 std::max(-x0-x1+dx,-x0+x1+dx),x0-x1+dx),x0+x1+dx);
295
296 G4double y0 = dz*fTthetaSphi;
297 G4double ymin = std::min(-y0-dy,y0-dy);
298 G4double ymax = std::max(-y0+dy,y0+dy);
299
300 pMin.set(xmin,ymin,-dz);
301 pMax.set(xmax,ymax, dz);
302
303 // Check correctness of the bounding box
304 //
305 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
306 {
307 std::ostringstream message;
308 message << "Bad bounding box (min >= max) for solid: "
309 << GetName() << " !"
310 << "\npMin = " << pMin
311 << "\npMax = " << pMax;
312 G4Exception("G4Para::BoundingLimits()", "GeomMgt0001",
313 JustWarning, message);
314 DumpInfo();
315 }
316}
317
318//////////////////////////////////////////////////////////////////////////
319//
320// Calculate extent under transform and specified limit
321
323 const G4VoxelLimits& pVoxelLimit,
324 const G4AffineTransform& pTransform,
325 G4double& pMin, G4double& pMax ) const
326{
327 G4ThreeVector bmin, bmax;
328 G4bool exist;
329
330 // Check bounding box (bbox)
331 //
332 BoundingLimits(bmin,bmax);
333 G4BoundingEnvelope bbox(bmin,bmax);
334#ifdef G4BBOX_EXTENT
335 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
336#endif
337 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
338 {
339 return exist = pMin < pMax;
340 }
341
342 // Set bounding envelope (benv) and calculate extent
343 //
347
348 G4double x0 = dz*fTthetaCphi;
349 G4double x1 = dy*GetTanAlpha();
350 G4double y0 = dz*fTthetaSphi;
351
352 G4ThreeVectorList baseA(4), baseB(4);
353 baseA[0].set(-x0-x1-dx,-y0-dy,-dz);
354 baseA[1].set(-x0-x1+dx,-y0-dy,-dz);
355 baseA[2].set(-x0+x1+dx,-y0+dy,-dz);
356 baseA[3].set(-x0+x1-dx,-y0+dy,-dz);
357
358 baseB[0].set(+x0-x1-dx, y0-dy, dz);
359 baseB[1].set(+x0-x1+dx, y0-dy, dz);
360 baseB[2].set(+x0+x1+dx, y0+dy, dz);
361 baseB[3].set(+x0+x1-dx, y0+dy, dz);
362
363 std::vector<const G4ThreeVectorList *> polygons(2);
364 polygons[0] = &baseA;
365 polygons[1] = &baseB;
366
367 G4BoundingEnvelope benv(bmin,bmax,polygons);
368 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
369 return exist;
370}
371
372//////////////////////////////////////////////////////////////////////////
373//
374// Determine where is point p, inside/on_surface/outside
375//
376
378{
379 G4double xx = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()+fPlanes[2].c*p.z();
380 G4double dx = std::abs(xx) + fPlanes[2].d;
381
382 G4double yy = fPlanes[0].b*p.y()+fPlanes[0].c*p.z();
383 G4double dy = std::abs(yy) + fPlanes[0].d;
384 G4double dxy = std::max(dx,dy);
385
386 G4double dz = std::abs(p.z())-fDz;
387 G4double dist = std::max(dxy,dz);
388
389 if (dist > halfCarTolerance) { return kOutside; }
390
391 return (dist > -halfCarTolerance) ? kSurface : kInside;
392}
393
394//////////////////////////////////////////////////////////////////////////
395//
396// Determine side where point is, and return corresponding normal
397
399{
400 G4int nsurf = 0; // number of surfaces where p is placed
401
402 // Check Z faces
403 //
404 G4double nz = 0;
405 G4double dz = std::abs(p.z()) - fDz;
406 if (std::abs(dz) <= halfCarTolerance)
407 {
408 nz = (p.z() < 0) ? -1 : 1;
409 ++nsurf;
410 }
411
412 // Check Y faces
413 //
414 G4double ny = 0;
415 G4double yy = fPlanes[0].b*p.y()+fPlanes[0].c*p.z();
416 if (std::abs(fPlanes[0].d + yy) <= halfCarTolerance)
417 {
418 ny = fPlanes[0].b;
419 nz += fPlanes[0].c;
420 ++nsurf;
421 }
422 else if (std::abs(fPlanes[1].d - yy) <= halfCarTolerance)
423 {
424 ny = fPlanes[1].b;
425 nz += fPlanes[1].c;
426 ++nsurf;
427 }
428
429 // Check X faces
430 //
431 G4double nx = 0;
432 G4double xx = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()+fPlanes[2].c*p.z();
433 if (std::abs(fPlanes[2].d + xx) <= halfCarTolerance)
434 {
435 nx = fPlanes[2].a;
436 ny += fPlanes[2].b;
437 nz += fPlanes[2].c;
438 ++nsurf;
439 }
440 else if (std::abs(fPlanes[3].d - xx) <= halfCarTolerance)
441 {
442 nx = fPlanes[3].a;
443 ny += fPlanes[3].b;
444 nz += fPlanes[3].c;
445 ++nsurf;
446 }
447
448 // Return normal
449 //
450 if (nsurf == 1)
451 {
452 return {nx,ny,nz};
453 }
454 if (nsurf != 0)
455 {
456 return G4ThreeVector(nx,ny,nz).unit(); // edge or corner
457 }
458
459 // Point is not on the surface
460 //
461#ifdef G4CSGDEBUG
462 std::ostringstream message;
463 G4int oldprc = message.precision(16);
464 message << "Point p is not on surface (!?) of solid: "
465 << GetName() << G4endl;
466 message << "Position:\n";
467 message << " p.x() = " << p.x()/mm << " mm\n";
468 message << " p.y() = " << p.y()/mm << " mm\n";
469 message << " p.z() = " << p.z()/mm << " mm";
470 G4cout.precision(oldprc) ;
471 G4Exception("G4Para::SurfaceNormal(p)", "GeomSolids1002",
472 JustWarning, message );
473 DumpInfo();
474#endif
475
476 return ApproxSurfaceNormal(p);
477}
478
479//////////////////////////////////////////////////////////////////////////
480//
481// Algorithm for SurfaceNormal() following the original specification
482// for points not on the surface
483
484G4ThreeVector G4Para::ApproxSurfaceNormal( const G4ThreeVector& p ) const
485{
486 G4double dist = -DBL_MAX;
487 G4int iside = 0;
488 for (G4int i=0; i<4; ++i)
489 {
490 G4double d = fPlanes[i].a*p.x() +
491 fPlanes[i].b*p.y() +
492 fPlanes[i].c*p.z() + fPlanes[i].d;
493 if (d > dist) { dist = d; iside = i; }
494 }
495
496 G4double distz = std::abs(p.z()) - fDz;
497 if (dist > distz)
498 {
499 return { fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c };
500 }
501
502 return { 0, 0, (G4double)((p.z() < 0) ? -1 : 1) };
503}
504
505//////////////////////////////////////////////////////////////////////////
506//
507// Calculate distance to shape from outside
508// - return kInfinity if no intersection
509
511 const G4ThreeVector& v ) const
512{
513 // Z intersections
514 //
515 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() >= 0)
516 {
517 return kInfinity;
518 }
519
520 G4double invz = (-v.z() == 0) ? DBL_MAX : -1./v.z();
521 G4double dz = (invz < 0) ? fDz : -fDz;
522 G4double tzmin = (p.z() + dz)*invz;
523 G4double tzmax = (p.z() - dz)*invz;
524
525 // Y intersections
526 //
527 G4double tmin0 = tzmin, tmax0 = tzmax;
528 G4double cos0 = fPlanes[0].b*v.y() + fPlanes[0].c*v.z();
529 G4double disy = fPlanes[0].b*p.y() + fPlanes[0].c*p.z();
530 G4double dis0 = fPlanes[0].d + disy;
531 if (dis0 >= -halfCarTolerance)
532 {
533 if (cos0 >= 0) { return kInfinity; }
534 G4double tmp = -dis0/cos0;
535 if (tmin0 < tmp) { tmin0 = tmp; }
536 }
537 else if (cos0 > 0)
538 {
539 G4double tmp = -dis0/cos0;
540 if (tmax0 > tmp) { tmax0 = tmp; }
541 }
542
543 G4double tmin1 = tmin0, tmax1 = tmax0;
544 G4double cos1 = -cos0;
545 G4double dis1 = fPlanes[1].d - disy;
546 if (dis1 >= -halfCarTolerance)
547 {
548 if (cos1 >= 0) { return kInfinity; }
549 G4double tmp = -dis1/cos1;
550 if (tmin1 < tmp) { tmin1 = tmp; }
551 }
552 else if (cos1 > 0)
553 {
554 G4double tmp = -dis1/cos1;
555 if (tmax1 > tmp) { tmax1 = tmp; }
556 }
557
558 // X intersections
559 //
560 G4double tmin2 = tmin1, tmax2 = tmax1;
561 G4double cos2 = fPlanes[2].a*v.x() + fPlanes[2].b*v.y() + fPlanes[2].c*v.z();
562 G4double disx = fPlanes[2].a*p.x() + fPlanes[2].b*p.y() + fPlanes[2].c*p.z();
563 G4double dis2 = fPlanes[2].d + disx;
564 if (dis2 >= -halfCarTolerance)
565 {
566 if (cos2 >= 0) { return kInfinity; }
567 G4double tmp = -dis2/cos2;
568 if (tmin2 < tmp) { tmin2 = tmp; }
569 }
570 else if (cos2 > 0)
571 {
572 G4double tmp = -dis2/cos2;
573 if (tmax2 > tmp) { tmax2 = tmp; }
574 }
575
576 G4double tmin3 = tmin2, tmax3 = tmax2;
577 G4double cos3 = -cos2;
578 G4double dis3 = fPlanes[3].d - disx;
579 if (dis3 >= -halfCarTolerance)
580 {
581 if (cos3 >= 0) { return kInfinity; }
582 G4double tmp = -dis3/cos3;
583 if (tmin3 < tmp) { tmin3 = tmp; }
584 }
585 else if (cos3 > 0)
586 {
587 G4double tmp = -dis3/cos3;
588 if (tmax3 > tmp) { tmax3 = tmp; }
589 }
590
591 // Find distance
592 //
593 G4double tmin = tmin3, tmax = tmax3;
594 if (tmax <= tmin + halfCarTolerance) { return kInfinity; } // touch or no hit
595
596 return (tmin < halfCarTolerance ) ? 0. : tmin;
597}
598
599//////////////////////////////////////////////////////////////////////////
600//
601// Calculate exact shortest distance to any boundary from outside
602// - returns 0 is point inside
603
605{
606 G4double xx = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()+fPlanes[2].c*p.z();
607 G4double dx = std::abs(xx) + fPlanes[2].d;
608
609 G4double yy = fPlanes[0].b*p.y()+fPlanes[0].c*p.z();
610 G4double dy = std::abs(yy) + fPlanes[0].d;
611 G4double dxy = std::max(dx,dy);
612
613 G4double dz = std::abs(p.z())-fDz;
614 G4double dist = std::max(dxy,dz);
615
616 return (dist > 0) ? dist : 0.;
617}
618
619//////////////////////////////////////////////////////////////////////////
620//
621// Calculate distance to surface of shape from inside and, if required,
622// find normal at exit point
623// - when leaving the surface, return 0
624
626 const G4bool calcNorm,
627 G4bool* validNorm, G4ThreeVector* n) const
628{
629 // Z intersections
630 //
631 if ((std::abs(p.z()) - fDz) >= -halfCarTolerance && p.z()*v.z() > 0)
632 {
633 if (calcNorm)
634 {
635 *validNorm = true;
636 n->set(0, 0, (p.z() < 0) ? -1 : 1);
637 }
638 return 0.;
639 }
640 G4double vz = v.z();
641 G4double tmax = (vz == 0) ? DBL_MAX : (std::copysign(fDz,vz) - p.z())/vz;
642 G4int iside = (vz < 0) ? -4 : -2; // little trick: (-4+3)=-1, (-2+3)=+1
643
644 // Y intersections
645 //
646 G4double cos0 = fPlanes[0].b*v.y() + fPlanes[0].c*v.z();
647 if (cos0 > 0)
648 {
649 G4double dis0 = fPlanes[0].b*p.y() + fPlanes[0].c*p.z() + fPlanes[0].d;
650 if (dis0 >= -halfCarTolerance)
651 {
652 if (calcNorm)
653 {
654 *validNorm = true;
655 n->set(0, fPlanes[0].b, fPlanes[0].c);
656 }
657 return 0.;
658 }
659 G4double tmp = -dis0/cos0;
660 if (tmax > tmp) { tmax = tmp; iside = 0; }
661 }
662
663 G4double cos1 = -cos0;
664 if (cos1 > 0)
665 {
666 G4double dis1 = fPlanes[1].b*p.y() + fPlanes[1].c*p.z() + fPlanes[1].d;
667 if (dis1 >= -halfCarTolerance)
668 {
669 if (calcNorm)
670 {
671 *validNorm = true;
672 n->set(0, fPlanes[1].b, fPlanes[1].c);
673 }
674 return 0.;
675 }
676 G4double tmp = -dis1/cos1;
677 if (tmax > tmp) { tmax = tmp; iside = 1; }
678 }
679
680 // X intersections
681 //
682 G4double cos2 = fPlanes[2].a*v.x() + fPlanes[2].b*v.y() + fPlanes[2].c*v.z();
683 if (cos2 > 0)
684 {
685 G4double dis2 = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()+fPlanes[2].c*p.z()+fPlanes[2].d;
686 if (dis2 >= -halfCarTolerance)
687 {
688 if (calcNorm)
689 {
690 *validNorm = true;
691 n->set(fPlanes[2].a, fPlanes[2].b, fPlanes[2].c);
692 }
693 return 0.;
694 }
695 G4double tmp = -dis2/cos2;
696 if (tmax > tmp) { tmax = tmp; iside = 2; }
697 }
698
699 G4double cos3 = -cos2;
700 if (cos3 > 0)
701 {
702 G4double dis3 = fPlanes[3].a*p.x()+fPlanes[3].b*p.y()
703 + fPlanes[3].c*p.z()+fPlanes[3].d;
704 if (dis3 >= -halfCarTolerance)
705 {
706 if (calcNorm)
707 {
708 *validNorm = true;
709 n->set(fPlanes[3].a, fPlanes[3].b, fPlanes[3].c);
710 }
711 return 0.;
712 }
713 G4double tmp = -dis3/cos3;
714 if (tmax > tmp) { tmax = tmp; iside = 3; }
715 }
716
717 // Set normal, if required, and return distance
718 //
719 if (calcNorm)
720 {
721 *validNorm = true;
722 (iside < 0) ? (n->set(0, 0, iside + 3)) // (-4+3)=-1, (-2+3)=+1
723 : (n->set(fPlanes[iside].a, fPlanes[iside].b, fPlanes[iside].c));
724 }
725 return tmax;
726}
727
728//////////////////////////////////////////////////////////////////////////
729//
730// Calculate exact shortest distance to any boundary from inside
731// - returns 0 is point outside
732
734{
735#ifdef G4CSGDEBUG
736 if( Inside(p) == kOutside )
737 {
738 std::ostringstream message;
739 G4int oldprc = message.precision(16);
740 message << "Point p is outside (!?) of solid: " << GetName() << G4endl;
741 message << "Position:\n";
742 message << " p.x() = " << p.x()/mm << " mm\n";
743 message << " p.y() = " << p.y()/mm << " mm\n";
744 message << " p.z() = " << p.z()/mm << " mm";
745 G4cout.precision(oldprc) ;
746 G4Exception("G4Para::DistanceToOut(p)", "GeomSolids1002",
747 JustWarning, message );
748 DumpInfo();
749 }
750#endif
751 G4double xx = fPlanes[2].a*p.x()+fPlanes[2].b*p.y()+fPlanes[2].c*p.z();
752 G4double dx = std::abs(xx) + fPlanes[2].d;
753
754 G4double yy = fPlanes[0].b*p.y()+fPlanes[0].c*p.z();
755 G4double dy = std::abs(yy) + fPlanes[0].d;
756 G4double dxy = std::max(dx,dy);
757
758 G4double dz = std::abs(p.z())-fDz;
759 G4double dist = std::max(dxy,dz);
760
761 return (dist < 0) ? -dist : 0.;
762}
763
764//////////////////////////////////////////////////////////////////////////
765//
766// GetEntityType
767
769{
770 return {"G4Para"};
771}
772
773//////////////////////////////////////////////////////////////////////////
774//
775// IsFaceted
776
778{
779 return true;
780}
781
782//////////////////////////////////////////////////////////////////////////
783//
784// Make a clone of the object
785//
787{
788 return new G4Para(*this);
789}
790
791//////////////////////////////////////////////////////////////////////////
792//
793// Stream object contents to an output stream
794
795std::ostream& G4Para::StreamInfo( std::ostream& os ) const
796{
797 G4double alpha = std::atan(fTalpha);
798 G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi +
799 fTthetaSphi*fTthetaSphi));
800 G4double phi = std::atan2(fTthetaSphi,fTthetaCphi);
801
802 G4long oldprc = os.precision(16);
803 os << "-----------------------------------------------------------\n"
804 << " *** Dump for solid - " << GetName() << " ***\n"
805 << " ===================================================\n"
806 << " Solid type: G4Para\n"
807 << " Parameters:\n"
808 << " half length X: " << fDx/mm << " mm\n"
809 << " half length Y: " << fDy/mm << " mm\n"
810 << " half length Z: " << fDz/mm << " mm\n"
811 << " alpha: " << alpha/degree << "degrees\n"
812 << " theta: " << theta/degree << "degrees\n"
813 << " phi: " << phi/degree << "degrees\n"
814 << "-----------------------------------------------------------\n";
815 os.precision(oldprc);
816
817 return os;
818}
819
820//////////////////////////////////////////////////////////////////////////
821//
822// Return a point randomly and uniformly selected on the solid surface
823
825{
826 G4double sxy = fDx*fDy;
827 G4double sxz = fDx*fDz*std::sqrt(1. + sqr(fTthetaSphi));
828 G4double syz = fDy*fDz*std::sqrt(1. + sqr(fTalpha) + sqr(fTalpha*fTthetaSphi - fTthetaCphi));
829
830 G4double select = (sxy + sxz + syz)*G4QuickRand();
831 G4double u = 2.*G4QuickRand() - 1.;
832 G4double v = 2.*G4QuickRand() - 1.;
833
834 G4double x, y, z;
835 if (select < sxy)
836 {
837 x = u*fDx;
838 y = v*fDy;
839 z = (select < 0.5*sxy) ? -fDz : fDz;
840 }
841 else if (select < sxy + sxz)
842 {
843 x = u*fDx;
844 y = (select < sxy + 0.5*sxz) ? -fDy : fDy;
845 z = v*fDz;
846 }
847 else
848 {
849 x = (select < sxy + sxz + 0.5*syz) ? -fDx : fDx;
850 y = u*fDy;
851 z = v*fDz;
852 }
853 return { x + y*fTalpha + z*fTthetaCphi, y + z*fTthetaSphi, z };
854}
855
856//////////////////////////////////////////////////////////////////////////
857//
858// Get volume
859
861{
862 if (fCubicVolume == 0)
863 {
864 G4AutoLock l(&paraMutex);
865 fCubicVolume = 8*fDx*fDy*fDz;
866 l.unlock();
867 }
868 return fCubicVolume;
869}
870
871//////////////////////////////////////////////////////////////////////////
872//
873// Get surface area
874
876{
877 if (fSurfaceArea == 0)
878 {
879 G4AutoLock l(&paraMutex);
880 G4double sxy = fDx*fDy;
881 G4double sxz = fDx*fDz*std::sqrt(1. + sqr(fTthetaSphi));
882 G4double syz = fDy*fDz*std::sqrt(1. + sqr(fTalpha) + sqr(fTalpha*fTthetaSphi - fTthetaCphi));
883 fSurfaceArea = 8*(sxy+sxz+syz);
884 l.unlock();
885 }
886 return fSurfaceArea;
887}
888
889//////////////////////////////////////////////////////////////////////////
890//
891// Methods for visualisation
892
894{
895 scene.AddSolid (*this);
896}
897
899{
900 G4double phi = std::atan2(fTthetaSphi, fTthetaCphi);
901 G4double alpha = std::atan(fTalpha);
902 G4double theta = std::atan(std::sqrt(fTthetaCphi*fTthetaCphi +
903 fTthetaSphi*fTthetaSphi));
904
905 return new G4PolyhedronPara(fDx, fDy, fDz, alpha, theta, phi);
906}
907#endif
G4TemplateAutoLock< G4Mutex > G4AutoLock
const G4double kCarTolerance
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
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
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
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 BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4double fSurfaceArea
Definition G4CSGSolid.hh:93
G4double fCubicVolume
Definition G4CSGSolid.hh:92
G4bool fRebuildPolyhedron
Definition G4CSGSolid.hh:94
G4CSGSolid(const G4String &pName)
Definition G4CSGSolid.cc:49
G4CSGSolid & operator=(const G4CSGSolid &rhs)
Definition G4CSGSolid.cc:89
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
Definition G4Para.cc:510
std::ostream & StreamInfo(std::ostream &os) const override
Definition G4Para.cc:795
EInside Inside(const G4ThreeVector &p) const override
Definition G4Para.cc:377
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
Definition G4Para.cc:322
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
Definition G4Para.cc:625
G4double GetSurfaceArea() override
Definition G4Para.cc:875
G4Polyhedron * CreatePolyhedron() const override
Definition G4Para.cc:898
G4double GetTanAlpha() const
G4GeometryType GetEntityType() const override
Definition G4Para.cc:768
G4ThreeVector GetPointOnSurface() const override
Definition G4Para.cc:824
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
Definition G4Para.cc:398
G4Para & operator=(const G4Para &rhs)
Definition G4Para.cc:155
G4Para(const G4String &pName, G4double pDx, G4double pDy, G4double pDz, G4double pAlpha, G4double pTheta, G4double pPhi)
Definition G4Para.cc:59
G4VSolid * Clone() const override
Definition G4Para.cc:786
void SetAllParameters(G4double pDx, G4double pDy, G4double pDz, G4double pAlpha, G4double pTheta, G4double pPhi)
Definition G4Para.cc:183
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
Definition G4Para.cc:279
G4bool IsFaceted() const override
Definition G4Para.cc:777
G4double GetYHalfLength() const
G4double d
Definition G4Para.hh:269
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
Definition G4Para.cc:268
G4double GetZHalfLength() const
G4double GetCubicVolume() override
Definition G4Para.cc:860
G4double a
Definition G4Para.hh:269
G4double GetXHalfLength() const
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
Definition G4Para.cc:893
G4double b
Definition G4Para.hh:269
G4double c
Definition G4Para.hh:269
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
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_MAX
Definition templates.hh:62