Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4GeomTestVolume.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: G4GeomTestVolume.cc 101727 2016-11-23 08:49:05Z gcosmo $
28 //
29 // --------------------------------------------------------------------
30 // GEANT 4 class source file
31 //
32 // G4GeomTestVolume
33 //
34 // Author: G.Cosmo, CERN
35 // --------------------------------------------------------------------
36 
37 #include <set>
38 
39 #include "G4GeomTestVolume.hh"
40 
41 #include "G4PhysicalConstants.hh"
42 #include "G4VPhysicalVolume.hh"
43 #include "G4LogicalVolume.hh"
44 #include "G4VSolid.hh"
45 
46 //
47 // Constructor
48 //
50  G4double theTolerance,
51  G4int numberOfPoints,
52  G4bool theVerbosity )
53  : target(theTarget), tolerance(theTolerance),
54  resolution(numberOfPoints), maxErr(1), verbosity(theVerbosity)
55 {;}
56 
57 //
58 // Destructor
59 //
61 
62 //
63 // Get error tolerance
64 //
66 {
67  return tolerance;
68 }
69 
70 //
71 // Set error tolerance
72 //
74 {
75  tolerance = tol;
76 }
77 
78 //
79 // Get number of points to check (resolution)
80 //
82 {
83  return resolution;
84 }
85 
86 //
87 // Set number of points to check (resolution)
88 //
90 {
91  resolution = np;
92 }
93 
94 //
95 // Get verbosity
96 //
98 {
99  return verbosity;
100 }
101 
102 //
103 // Set verbosity
104 //
106 {
107  verbosity = verb;
108 }
109 
110 //
111 // Get errors reporting threshold
112 //
114 {
115  return maxErr;
116 }
117 
118 //
119 // Set maximum number of errors to report
120 //
122 {
123  maxErr = max;
124 }
125 
126 //
127 // TestRecursiveOverlap
128 //
130 {
131  // If reached requested level of depth (i.e. set to 0), exit.
132  // If not depth specified (i.e. set to -1), visit the whole tree.
133  // If requested initial level of depth is not zero, visit from beginning
134  //
135  if (depth == 0) return;
136  if (depth != -1) depth--;
137  if (slevel != 0) slevel--;
138 
139  //
140  // As long as we reached the requested
141  // initial level of depth, test ourselves
142  //
143  if ( slevel==0 )
144  {
145  target->CheckOverlaps(resolution, tolerance, verbosity, maxErr);
146  }
147 
148  //
149  // Loop over unique daughters
150  //
151  std::set<const G4LogicalVolume *> tested;
152 
153  const G4LogicalVolume *logical = target->GetLogicalVolume();
154  G4int nDaughter = logical->GetNoDaughters();
155  G4int iDaughter;
156  for( iDaughter=0; iDaughter<nDaughter; ++iDaughter )
157  {
158  G4VPhysicalVolume *daughter = logical->GetDaughter(iDaughter);
159 
160  // Tested already?
161  //
162  // const G4LogicalVolume *daughterLogical =
163  // daughter->GetLogicalVolume();
164  // std::pair<std::set<const G4LogicalVolume *>::iterator, G4bool>
165  // there = tested.insert(daughterLogical);
166  // if (!there.second) continue;
167 
168  //
169  // Recurse
170  //
171  G4GeomTestVolume vTest( daughter, tolerance, resolution, verbosity );
172  vTest.SetErrorsThreshold(maxErr);
173  vTest.TestRecursiveOverlap( slevel,depth );
174  }
175 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetResolution(G4int points)
G4bool GetVerbosity() const
const XML_Char * target
Definition: expat.h:268
G4GeomTestVolume(G4VPhysicalVolume *theTarget, G4double theTolerance=0.0, G4int numberOfPoints=10000, G4bool theVerbosity=true)
void TestRecursiveOverlap(G4int sLevel=0, G4int depth=-1)
G4int GetNoDaughters() const
void SetTolerance(G4double tolerance)
void SetErrorsThreshold(G4int max)
G4double GetTolerance() const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ErrorTarget * theTarget
Definition: errprop.cc:59
G4int GetResolution() const
G4int GetErrorsThreshold() const
int G4int
Definition: G4Types.hh:78
void SetVerbosity(G4bool verbosity)
G4VPhysicalVolume * GetDaughter(const G4int i) const