Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GeomTestVolume.hh
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 //
27 // $Id$
28 //
29 // --------------------------------------------------------------------
30 // GEANT 4 class header file
31 //
32 // G4GeomTestVolume
33 //
34 // Class description:
35 //
36 // Checks for inconsistencies in the geometric boundaries of a physical
37 // volume and the boundaries of all its immediate daughters.
38 
39 // Author: D.C.Williams, UCSC (davidw@scipp.ucsc.edu)
40 // --------------------------------------------------------------------
41 #ifndef G4GeomTestVolume_hh
42 #define G4GeomTestVolume_hh
43 
44 #include "G4ThreeVector.hh"
45 #include "G4VisExtent.hh"
46 #include "G4GeomTestOverlapList.hh"
48 
49 #include <map>
50 
51 class G4VPhysicalVolume;
52 class G4GeomTestLogger;
53 
55 {
56  public: // with description
57 
59  G4GeomTestLogger *theLogger,
60  G4double theTolerance=1E-4 ); // mm
62  // Constructor and destructor
63 
64  G4double GetTolerance() const;
66  // Get/Set error tolerance (default set to 1E-4*mm)
67 
68  void TestCartGridXYZ( G4int nx=100, G4int ny=100, G4int nz=100 );
69  void TestCartGridX( G4int ny=100, G4int nz=100 );
70  void TestCartGridY( G4int nz=100, G4int nx=100 );
71  void TestCartGridZ( G4int nx=100, G4int ny=100 );
72  // Test using a grid of lines parallel to a cartesian axis
73 
74  void TestCartGrid( const G4ThreeVector &g1,
75  const G4ThreeVector &g2,
76  const G4ThreeVector &v,
77  G4int n1,
78  G4int n2 );
79  // Test using a grid of parallel lines
80  // g1 = First grid axis
81  // g2 = Second grid axis
82  // v = Direction of lines
83  // n1 = Number of grid points along g1
84  // n2 = Number of grid points along g2
85  // The spread of the grid points are automatically calculated
86  // based on the extent of the solid
87 
88  void TestRecursiveCartGrid( G4int nx=100, G4int ny=100, G4int nz=100,
89  G4int sLevel=0, G4int depth=-1 );
90  // Test using a grid, propagating recursively to the daughters, with
91  // possibility of specifying the initial level in the volume tree and
92  // the depth (default is the whole tree).
93  // Be careful: depending on the complexity of the geometry, this
94  // could require long computational time
95 
96  void TestCylinder( G4int nPhi=90, G4int nZ=50, G4int nRho=50,
97  G4double fracZ=0.8, G4double fracRho=0.8,
98  G4bool usePhi=false );
99  // Test using a set of lines in a cylindrical
100  // pattern of gradually increasing mesh size
101  // nPhi = Number lines per phi
102  // nZ = Number of z points
103  // nRho = Number of rho points
104  // fracZ = Fraction scale for points along z
105  // fracRho = Fraction scale for points along rho
106  // usePhi = Include phi set of lines
107  // Define a set of rho values such that:
108  // rho0 = Size of volume
109  // rho1 = frac*rho0
110  // rho2 = frac*rho1
111  // ... etc
112  // And define a set of z values
113  // z0 = z size of volume
114  // z1 = fracZ*z0
115  // z2 = fracZ*z1
116  // .... etc
117  // And define a set of nPhi phi values, evenly
118  // distributed in phi
119  //
120  // Three sets of lines are tested:
121  // * Imagine the set of lines parallel to the z axis
122  // through each rho point, at a phi angle taken the
123  // set of phi angles
124  // * Imagine the set of lines running perpendicular
125  // to the z axis and through a point on the z axis
126  // at +/- each z point and at an angle taken from the
127  // set of phi values
128  // * If usePhi==true, now take each pair of lines from the
129  // above two sets and imagine the line through the
130  // intersection and perpendicular to both
131 
132  void TestRecursiveCylinder( G4int nPhi=90, G4int nZ=50, G4int nRho=50,
133  G4double fracZ=0.8, G4double fracRho=0.8,
134  G4bool usePhi=false,
135  G4int sLevel=0, G4int depth=-1 );
136  // Test using a set of lines in a cylindrical pattern of gradually
137  // increasing mesh size, propagating recursively to the daughters, with
138  // possibility of specifying the initial level in the volume tree and
139  // the depth (default is the whole tree).
140  // Be careful: depending on the complexity of the geometry, this
141  // could require long computational time
142 
143  void TestOneLine( const G4ThreeVector &p, const G4ThreeVector &v );
144  // Test using a single line, specified by a point and direction,
145  // in the coordinate system of the target volume
146 
147  void TestRecursiveLine( const G4ThreeVector &p, const G4ThreeVector &v,
148  G4int sLevel=0, G4int depth=-1 );
149  // Test using a single line, specified by a point and direction,
150  // propagating recursively to the daughters, with possibility of
151  // specifying the initial level in the volume tree and
152  // the depth (default is the whole tree).
153 
154  void ReportErrors();
155  // Tabulate and report all errors so far encountered
156 
157  void ClearErrors();
158  // Clear list of errors
159 
160  protected:
161 
162  const G4VPhysicalVolume *target; // Target volume
163  G4GeomTestLogger *logger; // Error logger
164  G4double tolerance; // Error tolerance
165  G4VisExtent extent; // Geometric extent of volume
166 
167  std::map<G4long,G4GeomTestOverlapList> overlaps;
168  // A list of overlap errors, keyed by the
169  // daughter1*numDaughter+daughter2, where daughter1
170  // is the smaller of the two daughter indices
171 
172  std::map<G4long,G4GeomTestOvershootList> overshoots;
173  // A list of overshoot errors, keyed by the
174  // daughter number
175 };
176 
177 #endif