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