Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MultiUnion.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:$
28 //
29 // Implementation for G4MultiUnion class
30 //
31 // History:
32 //
33 // 06.04.17 G.Cosmo - Imported implementation in Geant4 for VecGeom migration
34 // 19.10.12 Marek Gayer - Original implementation from USolids module
35 // --------------------------------------------------------------------
36 
37 #include <iostream>
38 #include <sstream>
39 
40 #include "G4MultiUnion.hh"
41 #include "Randomize.hh"
42 #include "G4GeometryTolerance.hh"
43 #include "G4BoundingEnvelope.hh"
44 #include "G4AffineTransform.hh"
45 #include "G4DisplacedSolid.hh"
46 
47 #include "G4VGraphicsScene.hh"
48 #include "G4Polyhedron.hh"
49 #include "HepPolyhedronProcessor.h"
50 
51 #include "G4AutoLock.hh"
52 
53 namespace
54 {
55  G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
56 }
57 
58 //______________________________________________________________________________
60  : G4VSolid(name), fAccurate(false),
61  fRebuildPolyhedron(false), fpPolyhedron(0)
62 {
63  SetName(name);
64  fSolids.clear();
65  fTransformObjs.clear();
66  fCubicVolume = 0;
67  fSurfaceArea = 0;
69 }
70 
71 //______________________________________________________________________________
73 {
74 }
75 
76 //______________________________________________________________________________
78 {
79  fSolids.push_back(&solid);
80  fTransformObjs.push_back(trans); // Store a local copy of transformations
81 }
82 
83 //______________________________________________________________________________
85 {
86  return new G4MultiUnion(*this);
87 }
88 
89 // Copy constructor
90 //______________________________________________________________________________
92  : G4VSolid(rhs),fCubicVolume (rhs.fCubicVolume),
93  fSurfaceArea (rhs.fSurfaceArea),
94  kRadTolerance(rhs.kRadTolerance), fAccurate(false),
95  fRebuildPolyhedron(false), fpPolyhedron(0)
96 {
97 }
98 
99 // Fake default constructor for persistency
100 //______________________________________________________________________________
102  : G4VSolid(a), fCubicVolume (0.), fSurfaceArea (0.), kRadTolerance(0.),
103  fAccurate(false), fRebuildPolyhedron(false), fpPolyhedron(0)
104 {
105 }
106 
107 // Assignment operator
108 //______________________________________________________________________________
110 {
111  // Check assignment to self
112  //
113  if (this == &rhs)
114  {
115  return *this;
116  }
117 
118  // Copy base class data
119  //
120  G4VSolid::operator=(rhs);
121 
122  return *this;
123 }
124 
125 //______________________________________________________________________________
127 {
128  // Computes the cubic volume of the "G4MultiUnion" structure using
129  // random points
130 
131  if (!fCubicVolume)
132  {
133  G4ThreeVector extentMin, extentMax, d, p, point;
134  G4int inside = 0, generated;
135  BoundingLimits(extentMin, extentMax);
136  d = (extentMax - extentMin) / 2.;
137  p = (extentMax + extentMin) / 2.;
138  G4ThreeVector left = p - d;
139  G4ThreeVector length = d * 2;
140  for (generated = 0; generated < 10000; generated++)
141  {
143  point = left + G4ThreeVector(length.x()*rvec.x(),
144  length.y()*rvec.y(),
145  length.z()*rvec.z());
146  if (Inside(point) != EInside::kOutside) inside++;
147  }
148  G4double vbox = (2 * d.x()) * (2 * d.y()) * (2 * d.z());
149  fCubicVolume = inside * vbox / generated;
150  }
151  return fCubicVolume;
152 }
153 
154 //______________________________________________________________________________
155 G4double
157  const G4ThreeVector& aDirection) const
158 {
159  G4ThreeVector direction = aDirection.unit();
160  G4ThreeVector localPoint, localDirection;
161  G4double minDistance = kInfinity;
162 
163  G4int numNodes = fSolids.size();
164  for (G4int i = 0 ; i < numNodes ; ++i)
165  {
166  G4VSolid& solid = *fSolids[i];
167  const G4Transform3D& transform = fTransformObjs[i];
168 
169  localPoint = GetLocalPoint(transform, aPoint);
170  localDirection = GetLocalVector(transform, direction);
171 
172  G4double distance = solid.DistanceToIn(localPoint, localDirection);
173  if (minDistance > distance) minDistance = distance;
174  }
175  return minDistance;
176 }
177 
178 //______________________________________________________________________________
180  const G4ThreeVector& direction,
181  std::vector<G4int>& candidates,
182  G4SurfBits& bits) const
183 {
184  G4int candidatesCount = candidates.size();
185  G4ThreeVector localPoint, localDirection;
186 
187  G4double minDistance = kInfinity;
188  for (G4int i = 0 ; i < candidatesCount; ++i)
189  {
190  G4int candidate = candidates[i];
191  G4VSolid& solid = *fSolids[candidate];
192  const G4Transform3D& transform = fTransformObjs[candidate];
193 
194  localPoint = GetLocalPoint(transform, aPoint);
195  localDirection = GetLocalVector(transform, direction);
196  G4double distance = solid.DistanceToIn(localPoint, localDirection);
197  if (minDistance > distance) minDistance = distance;
198  bits.SetBitNumber(candidate);
199  if (minDistance == 0) break;
200  }
201  return minDistance;
202 }
203 
204 // Algorithm note: we have to look also for all other objects in next voxels,
205 // if the distance is not shorter ... we have to do it because,
206 // for example for objects which starts in first voxel in which they
207 // do not collide with direction line, but in second it collides...
208 // The idea of crossing voxels would be still applicable,
209 // because this way we could exclude from the testing such solids,
210 // which were found that obviously are not good candidates, because
211 // they would return infinity
212 // But if distance is smaller than the shift to next voxel, we can return
213 // it immediately
214 //______________________________________________________________________________
216  const G4ThreeVector& aDirection) const
217 {
218  G4double minDistance = kInfinity;
219  G4ThreeVector direction = aDirection.unit();
220  G4double shift = fVoxels.DistanceToFirst(aPoint, direction);
221  if (shift == kInfinity) return shift;
222 
223  G4ThreeVector currentPoint = aPoint;
224  if (shift) currentPoint += direction * shift;
225 
226  G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
227  std::vector<G4int> candidates, curVoxel(3);
228  fVoxels.GetVoxel(curVoxel, currentPoint);
229 
230  do
231  {
232  {
233  if (fVoxels.GetCandidatesVoxelArray(curVoxel, candidates, &exclusion))
234  {
235  G4double distance = DistanceToInCandidates(aPoint, direction,
236  candidates, exclusion);
237  if (minDistance > distance) minDistance = distance;
238  if (distance < shift) break;
239  }
240  }
241  shift = fVoxels.DistanceToNext(aPoint, direction, curVoxel);
242  }
243  while (minDistance > shift);
244 
245  return minDistance;
246 }
247 
248 //______________________________________________________________________________
250  const G4ThreeVector& aDirection,
251  G4ThreeVector* aNormal) const
252 {
253  // Computes distance from a point presumably outside the solid to the solid
254  // surface. Ignores first surface if the point is actually inside.
255  // Early return infinity in case the safety to any surface is found greater
256  // than the proposed step aPstep.
257  // The normal vector to the crossed surface is filled only in case the box
258  // is crossed, otherwise aNormal->IsNull() is true.
259 
260  // algorithm:
261  G4ThreeVector direction = aDirection.unit();
262  G4ThreeVector localPoint, localDirection;
263  G4int ignoredSolid = -1;
264  G4double resultDistToOut = 0;
265  G4ThreeVector currentPoint = aPoint;
266 
267  G4int numNodes = fSolids.size();
268  for (G4int i = 0; i < numNodes; ++i)
269  {
270  if (i != ignoredSolid)
271  {
272  G4VSolid& solid = *fSolids[i];
273  const G4Transform3D& transform = fTransformObjs[i];
274  localPoint = GetLocalPoint(transform, currentPoint);
275  localDirection = GetLocalVector(transform, direction);
276  EInside location = solid.Inside(localPoint);
277  if (location != EInside::kOutside)
278  {
279  G4double distance = solid.DistanceToOut(localPoint, localDirection,
280  aNormal);
281  if (distance < kInfinity)
282  {
283  if (resultDistToOut == kInfinity) resultDistToOut = 0;
284  if (distance > 0)
285  {
286  currentPoint = GetGlobalPoint(transform, localPoint
287  + distance*localDirection);
288  resultDistToOut += distance;
289  ignoredSolid = i; // skip the solid which we have just left
290  i = -1; // force the loop to continue from 0
291  }
292  }
293  }
294  }
295  }
296  return resultDistToOut;
297 }
298 
299 //______________________________________________________________________________
301  const G4ThreeVector& aDirection,
302  const G4bool /* calcNorm */,
303  G4bool* /* validNorm */,
304  G4ThreeVector* aNormal) const
305 {
306  return DistanceToOutVoxels(aPoint, aDirection, aNormal);
307 }
308 
309 //______________________________________________________________________________
311  const G4ThreeVector& aDirection,
312  G4ThreeVector* aNormal) const
313 {
314  // Computes distance from a point presumably inside the solid to the solid
315  // surface. Ignores first surface along each axis systematically (for points
316  // inside or outside. Early returns zero in case the second surface is behind
317  // the starting point.
318  // o The proposed step is ignored.
319  // o The normal vector to the crossed surface is always filled.
320 
321  // In the case the considered point is located inside the G4MultiUnion
322  // structure, the treatments are as follows:
323  // - investigation of the candidates for the passed point
324  // - progressive moving of the point towards the surface, along the
325  // passed direction
326  // - processing of the normal
327 
328  G4ThreeVector direction = aDirection.unit();
329  std::vector<G4int> candidates;
330  G4double distance = 0;
331  G4int numNodes = 2*fSolids.size();
332  G4int count=0;
333 
334  if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
335  {
336  // For normal case for which we presume the point is inside
337  G4ThreeVector localPoint, localDirection, localNormal;
338  G4ThreeVector currentPoint = aPoint;
339  G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
340  G4bool notOutside;
341  G4ThreeVector maxNormal;
342 
343  do
344  {
345  notOutside = false;
346 
347  G4double maxDistance = -kInfinity;
348  G4int maxCandidate = 0;
349  G4ThreeVector maxLocalPoint;
350 
351  G4int limit = candidates.size();
352  for (G4int i = 0 ; i < limit ; ++i)
353  {
354  G4int candidate = candidates[i];
355  // ignore the current component (that you just got out of) since
356  // numerically the propagated point will be on its surface
357 
358  G4VSolid& solid = *fSolids[candidate];
359  const G4Transform3D& transform = fTransformObjs[candidate];
360 
361  // The coordinates of the point are modified so as to fit the
362  // intrinsic solid local frame:
363  localPoint = GetLocalPoint(transform, currentPoint);
364 
365  // DistanceToOut at least for Trd sometimes return non-zero value
366  // even from points that are outside. Therefore, this condition
367  // must currently be here, otherwise it would not work.
368  // But it means it would be slower.
369 
370  if (solid.Inside(localPoint) != EInside::kOutside)
371  {
372  notOutside = true;
373 
374  localDirection = GetLocalVector(transform, direction);
375 
376  // propagate with solid.DistanceToOut
377  G4double shift = solid.DistanceToOut(localPoint, localDirection,
378  false, 0, &localNormal);
379  if (maxDistance < shift)
380  {
381  maxDistance = shift;
382  maxCandidate = candidate;
383  maxNormal = localNormal;
384  }
385  }
386  }
387 
388  if (notOutside)
389  {
390  const G4Transform3D& transform = fTransformObjs[maxCandidate];
391 
392  // convert from local normal
393  if (aNormal) *aNormal = GetGlobalVector(transform, maxNormal);
394 
395  distance += maxDistance;
396  currentPoint += maxDistance * direction;
397  if(maxDistance == 0.) count++;
398 
399  // the current component will be ignored
400  exclusion.SetBitNumber(maxCandidate);
401  EInside location = InsideWithExclusion(currentPoint, &exclusion);
402 
403  // perform a Inside
404  // it should be excluded current solid from checking
405  // we have to collect the maximum distance from all given candidates.
406  // such "maximum" candidate should be then used for finding next
407  // candidates
408  if (location == EInside::kOutside)
409  {
410  // else return cumulated distances to outside of the traversed
411  // components
412  break;
413  }
414  // if inside another component, redo 1 to 3 but add the next
415  // DistanceToOut on top of the previous.
416 
417  // and fill the candidates for the corresponding voxel (just
418  // exiting current component along direction)
419  candidates.clear();
420 
421  fVoxels.GetCandidatesVoxelArray(currentPoint, candidates, &exclusion);
422  exclusion.ResetBitNumber(maxCandidate);
423  }
424  }
425  while ((notOutside) && (count < numNodes));
426  }
427 
428  return distance;
429 }
430 
431 //______________________________________________________________________________
433  G4SurfBits* exclusion) const
434 {
435  // Classify point location with respect to solid:
436  // o eInside - inside the solid
437  // o eSurface - close to surface within tolerance
438  // o eOutside - outside the solid
439 
440  // Hitherto, it is considered that only parallelepipedic nodes
441  // can be added to the container
442 
443  // Implementation using voxelisation techniques:
444  // ---------------------------------------------
445 
446  G4ThreeVector localPoint;
447  EInside location = EInside::kOutside;
448 
449  std::vector<G4int> candidates;
450  std::vector<G4MultiUnionSurface> surfaces;
451 
452  // TODO: test if it works well and if so measure performance
453  // TODO: getPointIndex should not be used, instead GetVoxel + GetVoxelsIndex
454  // should be used
455  // TODO: than pass result to GetVoxel further to GetCandidatesVoxelArray
456  // TODO: eventually GetVoxel should be inlined here, early exit if any
457  // binary search is -1
458 
459  G4int limit = fVoxels.GetCandidatesVoxelArray(aPoint, candidates, exclusion);
460  for (G4int i = 0 ; i < limit ; ++i)
461  {
462  G4int candidate = candidates[i];
463  G4VSolid& solid = *fSolids[candidate];
464  const G4Transform3D& transform = fTransformObjs[candidate];
465 
466  // The coordinates of the point are modified so as to fit the intrinsic
467  // solid local frame:
468  localPoint = GetLocalPoint(transform, aPoint);
469  location = solid.Inside(localPoint);
470  if (location == EInside::kInside) return EInside::kInside;
471  else if (location == EInside::kSurface)
472  {
473  G4MultiUnionSurface surface;
474  surface.point = localPoint;
475  surface.solid = &solid;
476  surfaces.push_back(surface);
477  }
478  }
479 
481  // Important comment: When two solids touch each other along a flat
482  // surface, the surface points will be considered as kSurface, while points
483  // located around will correspond to kInside (cf. G4UnionSolid)
484 
485  G4int size = surfaces.size();
486  for (G4int i = 0; i < size - 1; ++i)
487  {
488  G4MultiUnionSurface& left = surfaces[i];
489  for (G4int j = i + 1; j < size; ++j)
490  {
491  G4MultiUnionSurface& right = surfaces[j];
492  G4ThreeVector n, n2;
493  n = left.solid->SurfaceNormal(left.point);
494  n2 = right.solid->SurfaceNormal(right.point);
495  if ((n + n2).mag2() < 1000 * kRadTolerance)
496  return EInside::kInside;
497  }
498  }
499 
500  location = size ? EInside::kSurface : EInside::kOutside;
501 
502  return location;
503 }
504 
505 //______________________________________________________________________________
507 {
508  // Classify point location with respect to solid:
509  // o eInside - inside the solid
510  // o eSurface - close to surface within tolerance
511  // o eOutside - outside the solid
512 
513  // Hitherto, it is considered that only parallelepipedic nodes can be
514  // added to the container
515 
516  // Implementation using voxelisation techniques:
517  // ---------------------------------------------
518 
519  // return InsideIterator(aPoint);
520 
521  EInside location = InsideWithExclusion(aPoint);
522  return location;
523 }
524 
525 //______________________________________________________________________________
527 {
528  G4ThreeVector localPoint;
529  EInside location = EInside::kOutside;
530  G4int countSurface = 0;
531 
532  G4int numNodes = fSolids.size();
533  for (G4int i = 0 ; i < numNodes ; ++i)
534  {
535  G4VSolid& solid = *fSolids[i];
536  G4Transform3D transform = GetTransformation(i);
537 
538  // The coordinates of the point are modified so as to fit the
539  // intrinsic solid local frame:
540  localPoint = GetLocalPoint(transform, aPoint);
541 
542  location = solid.Inside(localPoint);
543 
544  if (location == EInside::kSurface)
545  countSurface++;
546 
547  if (location == EInside::kInside) return EInside::kInside;
548  }
549  if (countSurface != 0) return EInside::kSurface;
550  return EInside::kOutside;
551 }
552 
553 //______________________________________________________________________________
554 void G4MultiUnion::Extent(EAxis aAxis, G4double& aMin, G4double& aMax) const
555 {
556  // Determines the bounding box for the considered instance of "UMultipleUnion"
558 
559  G4int numNodes = fSolids.size();
560  for (G4int i = 0 ; i < numNodes ; ++i)
561  {
562  G4VSolid& solid = *fSolids[i];
563  G4Transform3D transform = GetTransformation(i);
564  solid.BoundingLimits(min, max);
565 
566  TransformLimits(min, max, transform);
567 
568  if (i == 0)
569  {
570  switch (aAxis)
571  {
572  case kXAxis:
573  aMin = min.x();
574  aMax = max.x();
575  break;
576  case kYAxis:
577  aMin = min.y();
578  aMax = max.y();
579  break;
580  case kZAxis:
581  aMin = min.z();
582  aMax = max.z();
583  break;
584  default:
585  break;
586  }
587  }
588  else
589  {
590  // Determine the min/max on the considered axis:
591  switch (aAxis)
592  {
593  case kXAxis:
594  if (min.x() < aMin)
595  aMin = min.x();
596  if (max.x() > aMax)
597  aMax = max.x();
598  break;
599  case kYAxis:
600  if (min.y() < aMin)
601  aMin = min.y();
602  if (max.y() > aMax)
603  aMax = max.y();
604  break;
605  case kZAxis:
606  if (min.z() < aMin)
607  aMin = min.z();
608  if (max.z() > aMax)
609  aMax = max.z();
610  break;
611  default:
612  break;
613  }
614  }
615  }
616 }
617 
618 //______________________________________________________________________________
620  G4ThreeVector& aMax) const
621 {
622  Extent(kXAxis, aMin[0], aMax[0]);
623  Extent(kYAxis, aMin[1], aMax[1]);
624  Extent(kZAxis, aMin[2], aMax[2]);
625 }
626 
627 //______________________________________________________________________________
628 G4bool
630  const G4VoxelLimits& pVoxelLimit,
631  const G4AffineTransform& pTransform,
632  G4double& pMin, G4double& pMax) const
633 {
634  G4ThreeVector bmin, bmax;
635 
636  // Get bounding box
637  BoundingLimits(bmin,bmax);
638 
639  // Find extent
640  G4BoundingEnvelope bbox(bmin,bmax);
641  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
642 }
643 
644 //______________________________________________________________________________
646 {
647  // Computes the localNormal on a surface and returns it as a unit vector.
648  // Must return a valid vector. (even if the point is not on the surface).
649  //
650  // On an edge or corner, provide an average localNormal of all facets within
651  // tolerance
652  // NOTE: the tolerance value used in here is not yet the global surface
653  // tolerance - we will have to revise this value - TODO
654 
655  std::vector<G4int> candidates;
656  G4ThreeVector localPoint, normal, localNormal;
657  G4double safety = kInfinity;
658  G4int node = 0;
659 
661  // Important comment: Cases for which the point is located on an edge or
662  // on a vertice remain to be treated
663 
664  // determine weather we are in voxel area
665  if (fVoxels.GetCandidatesVoxelArray(aPoint, candidates))
666  {
667  G4int limit = candidates.size();
668  for (G4int i = 0 ; i < limit ; ++i)
669  {
670  G4int candidate = candidates[i];
671  const G4Transform3D& transform = fTransformObjs[candidate];
672 
673  // The coordinates of the point are modified so as to fit the intrinsic
674  // solid local frame:
675  localPoint = GetLocalPoint(transform, aPoint);
676  G4VSolid& solid = *fSolids[candidate];
677  EInside location = solid.Inside(localPoint);
678 
679  if (location == EInside::kSurface)
680  {
681  // normal case when point is on surface, we pick first solid
682  normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
683  return normal.unit();
684  }
685  else
686  {
687  // collect the smallest safety and remember solid node
688  G4double s = (location == EInside::kInside)
689  ? solid.DistanceToOut(localPoint)
690  : solid.DistanceToIn(localPoint);
691  if (s < safety)
692  {
693  safety = s;
694  node = candidate;
695  }
696  }
697  }
698  // on none of the solids, the point was not on the surface
699  G4VSolid& solid = *fSolids[node];
700  const G4Transform3D& transform = fTransformObjs[node];
701  localPoint = GetLocalPoint(transform, aPoint);
702 
703  normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
704  return normal.unit();
705  }
706  else
707  {
708  // for the case when point is certainly outside:
709 
710  // find a solid in union with the smallest safety
711  node = SafetyFromOutsideNumberNode(aPoint, safety);
712  G4VSolid& solid = *fSolids[node];
713 
714  const G4Transform3D& transform = fTransformObjs[node];
715  localPoint = GetLocalPoint(transform, aPoint);
716 
717  // evaluate normal for point at this found solid
718  // and transform multi-union coordinates
719  normal = GetGlobalVector(transform, solid.SurfaceNormal(localPoint));
720 
721  return normal.unit();
722  }
723 }
724 
725 //______________________________________________________________________________
727 {
728  // Estimates isotropic distance to the surface of the solid. This must
729  // be either accurate or an underestimate.
730  // Two modes: - default/fast mode, sacrificing accuracy for speed
731  // - "precise" mode, requests accurate value if available.
732 
733  std::vector<G4int> candidates;
734  G4ThreeVector localPoint;
735  G4double safetyMin = kInfinity;
736 
737  // In general, the value return by DistanceToIn(p) will not be the exact
738  // but only an undervalue (cf. overlaps)
739  fVoxels.GetCandidatesVoxelArray(point, candidates);
740 
741  G4int limit = candidates.size();
742  for (G4int i = 0; i < limit; ++i)
743  {
744  G4int candidate = candidates[i];
745 
746  // The coordinates of the point are modified so as to fit the intrinsic
747  // solid local frame:
748  const G4Transform3D& transform = fTransformObjs[candidate];
749  localPoint = GetLocalPoint(transform, point);
750  G4VSolid& solid = *fSolids[candidate];
751  if (solid.Inside(localPoint) == EInside::kInside)
752  {
753  G4double safety = solid.DistanceToOut(localPoint);
754  if (safetyMin > safety) safetyMin = safety;
755  }
756  }
757  if (safetyMin == kInfinity) safetyMin = 0; // we are not inside
758 
759  return safetyMin;
760 }
761 
762 //______________________________________________________________________________
764 {
765  // Estimates the isotropic safety from a point outside the current solid to
766  // any of its surfaces. The algorithm may be accurate or should provide a fast
767  // underestimate.
768 
769  if (!fAccurate) { return fVoxels.DistanceToBoundingBox(point); }
770 
771  const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
772  G4double safetyMin = kInfinity;
773  G4ThreeVector localPoint;
774 
775  G4int numNodes = fSolids.size();
776  for (G4int j = 0; j < numNodes; ++j)
777  {
778  G4ThreeVector dxyz;
779  if (j > 0)
780  {
781  const G4ThreeVector& pos = boxes[j].pos;
782  const G4ThreeVector& hlen = boxes[j].hlen;
783  for (G4int i = 0; i <= 2; ++i)
784  // distance to middle point - hlength => distance from point to border
785  // of x,y,z
786  if ((dxyz[i] = std::abs(point[i] - pos[i]) - hlen[i]) > safetyMin)
787  continue;
788 
789  G4double d2xyz = 0.;
790  for (G4int i = 0; i <= 2; ++i)
791  if (dxyz[i] > 0) d2xyz += dxyz[i] * dxyz[i];
792 
793  // minimal distance is at least this, but could be even higher. therefore,
794  // we can stop if previous was already lower, let us check if it does any
795  // chance to be better tha previous values...
796  if (d2xyz >= safetyMin * safetyMin)
797  {
798  continue;
799  }
800  }
801  const G4Transform3D& transform = fTransformObjs[j];
802  localPoint = GetLocalPoint(transform, point);
803  G4VSolid& solid = *fSolids[j];
804 
805  G4double safety = solid.DistanceToIn(localPoint);
806  if (safety <= 0) return safety;
807  // it was detected, that the point is not located outside
808  if (safetyMin > safety) safetyMin = safety;
809  }
810  return safetyMin;
811 }
812 
813 //______________________________________________________________________________
815 {
816  if (!fSurfaceArea)
817  {
818  fSurfaceArea = EstimateSurfaceArea(1000000, 0.001);
819  }
820  return fSurfaceArea;
821 }
822 
823 //______________________________________________________________________________
825 {
827 }
828 
829 //______________________________________________________________________________
831  G4double& safetyMin) const
832 {
833  // Method returning the closest node from a point located outside a
834  // G4MultiUnion.
835  // This is used to compute the normal in the case no candidate has been found.
836 
837  const std::vector<G4VoxelBox>& boxes = fVoxels.GetBoxes();
838  safetyMin = kInfinity;
839  G4int safetyNode = 0;
840  G4ThreeVector localPoint;
841 
842  G4int numNodes = fSolids.size();
843  for (G4int i = 0; i < numNodes; ++i)
844  {
845  G4double d2xyz = 0.;
846  G4double dxyz0 = std::abs(aPoint.x() - boxes[i].pos.x()) - boxes[i].hlen.x();
847  if (dxyz0 > safetyMin) continue;
848  G4double dxyz1 = std::abs(aPoint.y() - boxes[i].pos.y()) - boxes[i].hlen.y();
849  if (dxyz1 > safetyMin) continue;
850  G4double dxyz2 = std::abs(aPoint.z() - boxes[i].pos.z()) - boxes[i].hlen.z();
851  if (dxyz2 > safetyMin) continue;
852 
853  if (dxyz0 > 0) d2xyz += dxyz0 * dxyz0;
854  if (dxyz1 > 0) d2xyz += dxyz1 * dxyz1;
855  if (dxyz2 > 0) d2xyz += dxyz2 * dxyz2;
856  if (d2xyz >= safetyMin * safetyMin) continue;
857 
858  G4VSolid& solid = *fSolids[i];
859  const G4Transform3D& transform = fTransformObjs[i];
860  localPoint = GetLocalPoint(transform, aPoint);
861  fAccurate = true;
862  G4double safety = solid.DistanceToIn(localPoint);
863  fAccurate = false;
864  if (safetyMin > safety)
865  {
866  safetyMin = safety;
867  safetyNode = i;
868  }
869  }
870  return safetyNode;
871 }
872 
873 //______________________________________________________________________________
875  const G4Transform3D& transformation) const
876 {
877  // The goal of this method is to convert the quantities min and max
878  // (representing the bounding box of a given solid in its local frame)
879  // to the main frame, using "transformation"
880 
881  G4ThreeVector vertices[8] = // Detemination of the vertices thanks to
882  { // the extension of each solid:
883  G4ThreeVector(min.x(), min.y(), min.z()), // 1st vertice:
884  G4ThreeVector(min.x(), max.y(), min.z()), // 2nd vertice:
885  G4ThreeVector(max.x(), max.y(), min.z()),
886  G4ThreeVector(max.x(), min.y(), min.z()),
887  G4ThreeVector(min.x(), min.y(), max.z()),
888  G4ThreeVector(min.x(), max.y(), max.z()),
889  G4ThreeVector(max.x(), max.y(), max.z()),
890  G4ThreeVector(max.x(), min.y(), max.z())
891  };
892 
895 
896  // Loop on th vertices
897  G4int limit = sizeof(vertices) / sizeof(G4ThreeVector);
898  for (G4int i = 0 ; i < limit; ++i)
899  {
900  // From local frame to the global one:
901  // Current positions on the three axis:
902  G4ThreeVector current = GetGlobalPoint(transformation, vertices[i]);
903 
904  // If need be, replacement of the min & max values:
905  if (current.x() > max.x()) max.setX(current.x());
906  if (current.x() < min.x()) min.setX(current.x());
907 
908  if (current.y() > max.y()) max.setY(current.y());
909  if (current.y() < min.y()) min.setY(current.y());
910 
911  if (current.z() > max.z()) max.setZ(current.z());
912  if (current.z() < min.z()) min.setZ(current.z());
913  }
914 }
915 
916 // Stream object contents to an output stream
917 //______________________________________________________________________________
918 std::ostream& G4MultiUnion::StreamInfo(std::ostream& os) const
919 {
920  G4int oldprc = os.precision(16);
921  os << "-----------------------------------------------------------\n"
922  << " *** Dump for solid - " << GetName() << " ***\n"
923  << " ===================================================\n"
924  << " Solid type: G4MultiUnion\n"
925  << " Parameters: \n";
926  G4int numNodes = fSolids.size();
927  for (G4int i = 0 ; i < numNodes ; ++i)
928  {
929  G4VSolid& solid = *fSolids[i];
930  solid.StreamInfo(os);
931  const G4Transform3D& transform = fTransformObjs[i];
932  os << " Translation is " << transform.getTranslation() << " \n";
933  os << " Rotation is :" << " \n";
934  os << " " << transform.getRotation() << "\n";
935  }
936  os << " \n"
937  << "-----------------------------------------------------------\n";
938  os.precision(oldprc);
939 
940  return os;
941 }
942 
943 //______________________________________________________________________________
945 {
946  G4ThreeVector point;
947 
948  G4long size = fSolids.size();
949 
950  do
951  {
952  G4long rnd = G4RandFlat::shootInt(G4long(0), size);
953  G4VSolid& solid = *fSolids[rnd];
954  point = solid.GetPointOnSurface();
955  const G4Transform3D& transform = fTransformObjs[rnd];
956  point = GetGlobalPoint(transform, point);
957  }
958  while (Inside(point) != EInside::kSurface);
959 
960  return point;
961 }
962 
963 //______________________________________________________________________________
964 void
966 {
967  scene.AddSolid (*this);
968 }
969 
970 //______________________________________________________________________________
972 {
975 
976  G4VSolid* solidA = GetSolid(0);
977  const G4Transform3D transform0=GetTransformation(0);
978  G4DisplacedSolid dispSolidA("placedA",solidA,transform0);
979 
980  G4Polyhedron* top = new G4Polyhedron(*dispSolidA.GetPolyhedron());
981 
982  for(G4int i=1; i<GetNumberOfSolids(); ++i)
983  {
984  G4VSolid* solidB = GetSolid(i);
985  const G4Transform3D transform=GetTransformation(i);
986  G4DisplacedSolid dispSolidB("placedB",solidB,transform);
987  G4Polyhedron* operand = dispSolidB.GetPolyhedron();
988  processor.push_back (operation, *operand);
989  }
990 
991  if (processor.execute(*top)) { return top; }
992  else { return 0; }
993 }
994 
995 //______________________________________________________________________________
997 {
998  if (!fpPolyhedron ||
1002  {
1003  G4AutoLock l(&polyhedronMutex);
1004  delete fpPolyhedron;
1006  fRebuildPolyhedron = false;
1007  l.unlock();
1008  }
1009  return fpPolyhedron;
1010 }
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
void set(double x, double y, double z)
CLHEP::HepRotation getRotation() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const XML_Char * name
Definition: expat.h:151
G4Polyhedron * GetPolyhedron() const
CLHEP::Hep3Vector getTranslation() const
CLHEP::Hep3Vector G4ThreeVector
G4VSolid * GetSolid(G4int index) const
void Voxelize(std::vector< G4VSolid * > &solids, std::vector< G4Transform3D > &transforms)
Definition: G4Voxelizer.cc:724
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double pos
G4Polyhedron * CreatePolyhedron() const
std::vector< G4Transform3D > fTransformObjs
virtual void AddSolid(const G4Box &)=0
static const G4double kInfinity
Definition: geomdefs.hh:42
G4ThreeVector GetGlobalPoint(const G4Transform3D &trans, const G4ThreeVector &lpoint) const
G4int GetCandidatesVoxelArray(const G4ThreeVector &point, std::vector< G4int > &list, G4SurfBits *crossed=0) const
Definition: G4Voxelizer.cc:984
void AddNode(G4VSolid &solid, G4Transform3D &trans)
Definition: G4MultiUnion.cc:77
G4bool fAccurate
G4ThreeVector SurfaceNormal(const G4ThreeVector &aPoint) const
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:625
std::ostream & StreamInfo(std::ostream &os) const
const char * p
Definition: xmltok.h:285
G4ThreeVector GetPointOnSurface() const
G4VSolid * Clone() const
Definition: G4MultiUnion.cc:84
G4double DistanceToInNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection) const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
G4Polyhedron * GetPolyhedron() const
G4double fSurfaceArea
double z() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double GetRadialTolerance() const
G4double DistanceToIn(const G4ThreeVector &aPoint) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
EInside InsideNoVoxels(const G4ThreeVector &aPoint) const
void push_back(Operation, const HepPolyhedron &)
G4double GetSurfaceArea()
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
void setX(double)
G4ThreeVector GetLocalPoint(const G4Transform3D &trans, const G4ThreeVector &gpoint) const
G4ThreeVector GetGlobalVector(const G4Transform3D &trans, const G4ThreeVector &lvec) const
virtual EInside Inside(const G4ThreeVector &p) const =0
G4double kRadTolerance
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
const XML_Char * s
Definition: expat.h:262
G4double DistanceToFirst(const G4ThreeVector &point, const G4ThreeVector &direction) const
G4double DistanceToOutNoVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:88
void BoundingLimits(G4ThreeVector &aMin, G4ThreeVector &aMax) const
void TransformLimits(G4ThreeVector &min, G4ThreeVector &max, const G4Transform3D &transformation) const
long G4long
Definition: G4Types.hh:80
void setZ(double)
const G4Transform3D & GetTransformation(G4int index) const
#define G4UniformRand()
Definition: Randomize.hh:53
void SetName(const G4String &name)
Float_t d
bool execute(HepPolyhedron &)
virtual std::ostream & StreamInfo(std::ostream &os) const =0
const std::vector< G4VoxelBox > & GetBoxes() const
G4double GetCubicVolume()
Hep3Vector unit() const
static G4int GetNumberOfRotationSteps()
G4bool fRebuildPolyhedron
EInside InsideWithExclusion(const G4ThreeVector &aPoint, G4SurfBits *bits=0) const
#define processor
Definition: xmlparse.cc:617
G4ThreeVector GetLocalVector(const G4Transform3D &trans, const G4ThreeVector &gvec) const
G4Polyhedron * fpPolyhedron
G4double DistanceToNext(const G4ThreeVector &point, const G4ThreeVector &direction, std::vector< G4int > &curVoxel) const
G4double fCubicVolume
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
G4double DistanceToOutVoxels(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, G4ThreeVector *aNormalVector) const
G4String GetName() const
std::vector< G4VSolid * > fSolids
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:263
double x() const
G4int GetBitsPerSlice() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4MultiUnion & operator=(const G4MultiUnion &rhs)
static G4GeometryTolerance * GetInstance()
Char_t n[5]
void SetBitNumber(unsigned int bitnumber, G4bool value=true)
Definition: G4SurfBits.hh:123
G4double DistanceToBoundingBox(const G4ThreeVector &point) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void Extent(EAxis aAxis, G4double &aMin, G4double &aMax) const
double y() const
EInside Inside(const G4ThreeVector &aPoint) const
void GetVoxel(std::vector< G4int > &curVoxel, const G4ThreeVector &point) const
G4Voxelizer fVoxels
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double DistanceToOut(const G4ThreeVector &aPoint) const
G4int GetNumberOfSolids() const
void setY(double)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static int operand(pchar begin, pchar end, double &result, pchar &endp, const dic_type &dictionary)
Definition: Evaluator.cc:164
G4double DistanceToInCandidates(const G4ThreeVector &aPoint, const G4ThreeVector &aDirection, std::vector< G4int > &candidates, G4SurfBits &bits) const
G4int SafetyFromOutsideNumberNode(const G4ThreeVector &aPoint, G4double &safety) const
std::mutex G4Mutex
Definition: G4Threading.hh:84