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: testG4Trd.cc,v 1.10 2006/06/29 18:46:29 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ |
---|
29 | // |
---|
30 | |
---|
31 | // testG4Trd |
---|
32 | // Ensure asserts are compiled in |
---|
33 | |
---|
34 | #include <assert.h> |
---|
35 | #include <cmath> |
---|
36 | |
---|
37 | #include "globals.hh" |
---|
38 | #include "geomdefs.hh" |
---|
39 | |
---|
40 | #include "ApproxEqual.hh" |
---|
41 | |
---|
42 | #include "G4ThreeVector.hh" |
---|
43 | #include "G4Trd.hh" |
---|
44 | #include "G4Box.hh" |
---|
45 | #include "G4RotationMatrix.hh" |
---|
46 | #include "G4AffineTransform.hh" |
---|
47 | #include "G4VoxelLimits.hh" |
---|
48 | |
---|
49 | /////////////////////////////////////////////////////////////////// |
---|
50 | // |
---|
51 | // Dave's auxiliary function |
---|
52 | |
---|
53 | const G4String OutputInside(const EInside a) |
---|
54 | { |
---|
55 | switch(a) |
---|
56 | { |
---|
57 | case kInside : return "kInside" ; |
---|
58 | case kOutside: return "kOutside" ; |
---|
59 | case kSurface: return "kSurface" ; |
---|
60 | } |
---|
61 | return "????"; |
---|
62 | } |
---|
63 | |
---|
64 | |
---|
65 | G4bool testG4Trd() |
---|
66 | { |
---|
67 | EInside inside ; |
---|
68 | G4ThreeVector pzero(0,0,0); |
---|
69 | G4ThreeVector ponxside(20,0,0),ponyside(0,30,0),ponzside(0,0,40); |
---|
70 | G4ThreeVector ponmxside(-20,0,0),ponmyside(0,-30,0),ponmzside(0,0,-40); |
---|
71 | G4ThreeVector ponzsidey(0,25,40),ponmzsidey(0,25,-40); |
---|
72 | |
---|
73 | G4ThreeVector pbigx(100,0,0),pbigy(0,100,0),pbigz(0,0,100); |
---|
74 | G4ThreeVector pbigmx(-100,0,0),pbigmy(0,-100,0),pbigmz(0,0,-100); |
---|
75 | |
---|
76 | G4ThreeVector vx(1,0,0),vy(0,1,0),vz(0,0,1); |
---|
77 | G4ThreeVector vmx(-1,0,0),vmy(0,-1,0),vmz(0,0,-1); |
---|
78 | G4ThreeVector vxy(1/std::sqrt(2.0),1/std::sqrt(2.0),0); |
---|
79 | G4ThreeVector vmxy(-1/std::sqrt(2.0),1/std::sqrt(2.0),0); |
---|
80 | G4ThreeVector vmxmy(-1/std::sqrt(2.0),-1/std::sqrt(2.0),0); |
---|
81 | G4ThreeVector vxmy(1/std::sqrt(2.0),-1/std::sqrt(2.0),0); |
---|
82 | |
---|
83 | G4double Dist, vol, volCheck; |
---|
84 | G4ThreeVector *pNorm,norm; |
---|
85 | G4bool *pgoodNorm,goodNorm,calcNorm=true; |
---|
86 | |
---|
87 | pNorm=&norm; |
---|
88 | pgoodNorm=&goodNorm; |
---|
89 | |
---|
90 | G4Trd trd1("Test Box #1",20,20,30,30,40); |
---|
91 | G4Trd trd2("Test Trd",10,30,20,40,40); |
---|
92 | G4Trd trd3("BABAR Trd",0.14999999999999999,0.14999999999999999, |
---|
93 | 24.707000000000001, 24.707000000000001, |
---|
94 | 22.699999999999999) ; |
---|
95 | |
---|
96 | |
---|
97 | // Check name |
---|
98 | assert(trd1.GetName()=="Test Box #1"); |
---|
99 | assert(trd2.GetName()=="Test Trd"); |
---|
100 | |
---|
101 | // check cubic volume |
---|
102 | |
---|
103 | vol = trd1.GetCubicVolume(); |
---|
104 | volCheck = 8*20*30*40; |
---|
105 | assert(ApproxEqual(vol,volCheck)); |
---|
106 | |
---|
107 | // Check Inside |
---|
108 | |
---|
109 | assert(trd1.Inside(pzero)==kInside); |
---|
110 | assert(trd1.Inside(pbigz)==kOutside); |
---|
111 | assert(trd1.Inside(ponxside)==kSurface); |
---|
112 | assert(trd1.Inside(ponyside)==kSurface); |
---|
113 | assert(trd1.Inside(ponzside)==kSurface); |
---|
114 | |
---|
115 | inside = trd1.Inside(G4ThreeVector(20,30,40)) ; |
---|
116 | // G4cout << "trd1.Inside((20,30,40)) = " << OutputInside(inside) << G4endl ; |
---|
117 | assert(inside == kSurface); |
---|
118 | |
---|
119 | inside = trd1.Inside(G4ThreeVector(-20,30,40)) ; |
---|
120 | // G4cout << "trd1.Inside((-20,30,40)) = " << OutputInside(inside) << G4endl ; |
---|
121 | assert(inside == kSurface); |
---|
122 | |
---|
123 | inside = trd1.Inside(G4ThreeVector(20,-30,40)) ; |
---|
124 | // G4cout << "trd1.Inside((20,-30,40)) = " << OutputInside(inside) << G4endl ; |
---|
125 | assert(inside == kSurface); |
---|
126 | |
---|
127 | inside = trd1.Inside(G4ThreeVector(20,30,-40)) ; |
---|
128 | // G4cout << "trd1.Inside((20,30,-40)) = " << OutputInside(inside) << G4endl ; |
---|
129 | assert(inside == kSurface); |
---|
130 | |
---|
131 | inside = trd1.Inside(G4ThreeVector(20,30,0)) ; |
---|
132 | // G4cout << "trd1.Inside((20,30,0)) = " << OutputInside(inside) << G4endl ; |
---|
133 | assert(inside == kSurface); |
---|
134 | |
---|
135 | inside = trd1.Inside(G4ThreeVector(0,30,40)) ; |
---|
136 | // G4cout << "trd1.Inside((0,30,40)) = " << OutputInside(inside) << G4endl ; |
---|
137 | assert(inside == kSurface); |
---|
138 | |
---|
139 | inside = trd1.Inside(G4ThreeVector(20,0,40)) ; |
---|
140 | // G4cout << "trd1.Inside((20,0,40)) = " << OutputInside(inside) << G4endl ; |
---|
141 | assert(inside == kSurface); |
---|
142 | |
---|
143 | inside = trd1.Inside(G4ThreeVector(-20,-30,-40)) ; |
---|
144 | // G4cout << "trd1.Inside((-20,-30,-40)) = " << OutputInside(inside) << G4endl ; |
---|
145 | assert(inside == kSurface); |
---|
146 | |
---|
147 | |
---|
148 | |
---|
149 | |
---|
150 | |
---|
151 | assert(trd2.Inside(pzero)==kInside); |
---|
152 | assert(trd2.Inside(pbigz)==kOutside); |
---|
153 | assert(trd2.Inside(ponxside)==kSurface); |
---|
154 | assert(trd2.Inside(ponyside)==kSurface); |
---|
155 | assert(trd2.Inside(ponzside)==kSurface); |
---|
156 | |
---|
157 | // Check Surface Normal |
---|
158 | G4ThreeVector normal; |
---|
159 | |
---|
160 | normal=trd1.SurfaceNormal(ponxside); |
---|
161 | assert(ApproxEqual(normal,G4ThreeVector(1,0,0))); |
---|
162 | normal=trd1.SurfaceNormal(ponmxside); |
---|
163 | assert(ApproxEqual(normal,G4ThreeVector(-1,0,0))); |
---|
164 | normal=trd1.SurfaceNormal(ponyside); |
---|
165 | assert(ApproxEqual(normal,G4ThreeVector(0,1,0))); |
---|
166 | normal=trd1.SurfaceNormal(ponmyside); |
---|
167 | assert(ApproxEqual(normal,G4ThreeVector(0,-1,0))); |
---|
168 | normal=trd1.SurfaceNormal(ponzside); |
---|
169 | assert(ApproxEqual(normal,G4ThreeVector(0,0,1))); |
---|
170 | normal=trd1.SurfaceNormal(ponmzside); |
---|
171 | assert(ApproxEqual(normal,G4ThreeVector(0,0,-1))); |
---|
172 | normal=trd1.SurfaceNormal(ponzsidey); |
---|
173 | assert(ApproxEqual(normal,G4ThreeVector(0,0,1))); |
---|
174 | normal=trd1.SurfaceNormal(ponmzsidey); |
---|
175 | assert(ApproxEqual(normal,G4ThreeVector(0,0,-1))); |
---|
176 | |
---|
177 | |
---|
178 | // Normals on Edges |
---|
179 | |
---|
180 | G4ThreeVector edgeXY( 20.0, 30., 0.0); |
---|
181 | G4ThreeVector edgemXmY( -20.0, -30., 0.0); |
---|
182 | G4ThreeVector edgeXmY( 20.0, -30., 0.0); |
---|
183 | G4ThreeVector edgemXY( -20.0, 30., 0.0); |
---|
184 | G4ThreeVector edgeXZ( 20.0, 0.0, 40.0); |
---|
185 | G4ThreeVector edgemXmZ( -20.0, 0.0, -40.0); |
---|
186 | G4ThreeVector edgeXmZ( 20.0, 0.0, -40.0); |
---|
187 | G4ThreeVector edgemXZ( -20.0, 0.0, 40.0); |
---|
188 | G4ThreeVector edgeYZ( 0.0, 30.0, 40.0); |
---|
189 | G4ThreeVector edgemYmZ( 0.0, -30.0, -40.0); |
---|
190 | G4ThreeVector edgeYmZ( 0.0, 30.0, -40.0); |
---|
191 | G4ThreeVector edgemYZ( 0.0, -30.0, 40.0); |
---|
192 | |
---|
193 | G4double invSqrt2 = 1.0 / std::sqrt( 2.0); |
---|
194 | G4double invSqrt3 = 1.0 / std::sqrt( 3.0); |
---|
195 | |
---|
196 | normal= trd1.SurfaceNormal( edgeXY ); |
---|
197 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt2, invSqrt2, 0.0) )); |
---|
198 | |
---|
199 | // G4cout << " Normal at " << edgeXY << " is " << normal |
---|
200 | // << " Expected is " << G4ThreeVector( invSqrt2, invSqrt2, 0.0) << G4endl; |
---|
201 | |
---|
202 | normal= trd1.SurfaceNormal( edgemXmY ); |
---|
203 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt2, -invSqrt2, 0.0) )); |
---|
204 | normal= trd1.SurfaceNormal( edgeXmY ); |
---|
205 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt2, -invSqrt2, 0.0) )); |
---|
206 | normal= trd1.SurfaceNormal( edgemXY ); |
---|
207 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt2, invSqrt2, 0.0) )); |
---|
208 | |
---|
209 | normal= trd1.SurfaceNormal( edgeXZ ); |
---|
210 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt2, 0.0, invSqrt2) )); |
---|
211 | normal= trd1.SurfaceNormal( edgemXmZ ); |
---|
212 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt2, 0.0, -invSqrt2) )); |
---|
213 | normal= trd1.SurfaceNormal( edgeXmZ ); |
---|
214 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt2, 0.0, -invSqrt2) )); |
---|
215 | normal= trd1.SurfaceNormal( edgemXZ ); |
---|
216 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt2, 0.0, invSqrt2) )); |
---|
217 | |
---|
218 | normal= trd1.SurfaceNormal( edgeYZ ); |
---|
219 | assert(ApproxEqual( normal, G4ThreeVector( 0.0, invSqrt2, invSqrt2) )); |
---|
220 | normal= trd1.SurfaceNormal( edgemYmZ ); |
---|
221 | assert(ApproxEqual( normal, G4ThreeVector( 0.0, -invSqrt2, -invSqrt2) )); |
---|
222 | normal= trd1.SurfaceNormal( edgeYmZ ); |
---|
223 | assert(ApproxEqual( normal, G4ThreeVector( 0.0, invSqrt2, -invSqrt2) )); |
---|
224 | normal= trd1.SurfaceNormal( edgemYZ ); |
---|
225 | assert(ApproxEqual( normal, G4ThreeVector( 0.0, -invSqrt2, invSqrt2) )); |
---|
226 | |
---|
227 | // Normals on corners |
---|
228 | |
---|
229 | G4ThreeVector cornerXYZ( 20.0, 30., 40.0); |
---|
230 | G4ThreeVector cornermXYZ( -20.0, 30., 40.0); |
---|
231 | G4ThreeVector cornerXmYZ( 20.0, -30., 40.0); |
---|
232 | G4ThreeVector cornermXmYZ( -20.0, -30., 40.0); |
---|
233 | G4ThreeVector cornerXYmZ( 20.0, 30., -40.0); |
---|
234 | G4ThreeVector cornermXYmZ( -20.0, 30., -40.0); |
---|
235 | G4ThreeVector cornerXmYmZ( 20.0, -30., -40.0); |
---|
236 | G4ThreeVector cornermXmYmZ( -20.0, -30., -40.0); |
---|
237 | |
---|
238 | normal= trd1.SurfaceNormal( cornerXYZ ); |
---|
239 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt3, invSqrt3, invSqrt3) )); |
---|
240 | normal= trd1.SurfaceNormal( cornermXYZ ); |
---|
241 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt3, invSqrt3, invSqrt3) )); |
---|
242 | normal= trd1.SurfaceNormal( cornerXmYZ ); |
---|
243 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt3, -invSqrt3, invSqrt3) )); |
---|
244 | normal= trd1.SurfaceNormal( cornermXmYZ ); |
---|
245 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt3, -invSqrt3, invSqrt3) )); |
---|
246 | normal= trd1.SurfaceNormal( cornerXYmZ ); |
---|
247 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt3, invSqrt3, -invSqrt3) )); |
---|
248 | normal= trd1.SurfaceNormal( cornermXYmZ ); |
---|
249 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt3, invSqrt3, -invSqrt3) )); |
---|
250 | normal= trd1.SurfaceNormal( cornerXmYmZ ); |
---|
251 | assert(ApproxEqual( normal, G4ThreeVector( invSqrt3, -invSqrt3, -invSqrt3) )); |
---|
252 | normal= trd1.SurfaceNormal( cornermXmYmZ ); |
---|
253 | assert(ApproxEqual( normal, G4ThreeVector( -invSqrt3, -invSqrt3, -invSqrt3) )); |
---|
254 | |
---|
255 | |
---|
256 | |
---|
257 | |
---|
258 | |
---|
259 | |
---|
260 | double cosa = 4/std::sqrt(17.), sina = 1/std::sqrt(17.), tanga = 1.0/4.0 ; |
---|
261 | |
---|
262 | normal=trd2.SurfaceNormal(ponxside); |
---|
263 | assert(ApproxEqual(normal,G4ThreeVector(cosa,0,-sina))); |
---|
264 | normal=trd2.SurfaceNormal(ponmxside); |
---|
265 | assert(ApproxEqual(normal,G4ThreeVector(-cosa,0,-sina))); |
---|
266 | normal=trd2.SurfaceNormal(ponyside); |
---|
267 | assert(ApproxEqual(normal,G4ThreeVector(0,cosa,-sina))); |
---|
268 | normal=trd2.SurfaceNormal(ponmyside); |
---|
269 | assert(ApproxEqual(normal,G4ThreeVector(0,-cosa,-sina))); |
---|
270 | normal=trd2.SurfaceNormal(ponzside); |
---|
271 | assert(ApproxEqual(normal,G4ThreeVector(0,0,1))); |
---|
272 | normal=trd2.SurfaceNormal(ponmzside); |
---|
273 | assert(ApproxEqual(normal,G4ThreeVector(0,0,-1))); |
---|
274 | normal=trd2.SurfaceNormal(ponzsidey); |
---|
275 | assert(ApproxEqual(normal,G4ThreeVector(0,0,1))); |
---|
276 | normal=trd2.SurfaceNormal(ponmzsidey); |
---|
277 | assert(ApproxEqual(normal,G4ThreeVector(0,0,-1))); // (0,cosa,-sina) ? |
---|
278 | |
---|
279 | // DistanceToOut(P) |
---|
280 | |
---|
281 | Dist=trd1.DistanceToOut(pzero); |
---|
282 | assert(ApproxEqual(Dist,20)); |
---|
283 | Dist=trd1.DistanceToOut(vx); |
---|
284 | assert(ApproxEqual(Dist,19)); |
---|
285 | Dist=trd1.DistanceToOut(vy); |
---|
286 | assert(ApproxEqual(Dist,20)); |
---|
287 | Dist=trd1.DistanceToOut(vz); |
---|
288 | assert(ApproxEqual(Dist,20)); |
---|
289 | |
---|
290 | Dist=trd2.DistanceToOut(pzero); |
---|
291 | assert(ApproxEqual(Dist,20*cosa)); |
---|
292 | Dist=trd2.DistanceToOut(vx); |
---|
293 | assert(ApproxEqual(Dist,19*cosa)); |
---|
294 | Dist=trd2.DistanceToOut(vy); |
---|
295 | assert(ApproxEqual(Dist,20*cosa)); |
---|
296 | Dist=trd2.DistanceToOut(vz); |
---|
297 | assert(ApproxEqual(Dist,20*cosa+sina)); |
---|
298 | |
---|
299 | |
---|
300 | // DistanceToOut(P,V) |
---|
301 | |
---|
302 | Dist=trd1.DistanceToOut(pzero,vx,calcNorm,pgoodNorm,pNorm); |
---|
303 | assert(ApproxEqual(Dist,20)&&ApproxEqual(*pNorm,vx)&&*pgoodNorm); |
---|
304 | Dist=trd1.DistanceToOut(pzero,vmx,calcNorm,pgoodNorm,pNorm); |
---|
305 | assert(ApproxEqual(Dist,20)&&ApproxEqual(norm,vmx)&&*pgoodNorm); |
---|
306 | Dist=trd1.DistanceToOut(pzero,vy,calcNorm,pgoodNorm,pNorm); |
---|
307 | assert(ApproxEqual(Dist,30)&&ApproxEqual(norm,vy)&&*pgoodNorm); |
---|
308 | Dist=trd1.DistanceToOut(pzero,vmy,calcNorm,pgoodNorm,pNorm); |
---|
309 | assert(ApproxEqual(Dist,30)&&ApproxEqual(norm,vmy)&&*pgoodNorm); |
---|
310 | Dist=trd1.DistanceToOut(pzero,vz,calcNorm,pgoodNorm,pNorm); |
---|
311 | assert(ApproxEqual(Dist,40)&&ApproxEqual(norm,vz)&&*pgoodNorm); |
---|
312 | Dist=trd1.DistanceToOut(pzero,vmz,calcNorm,pgoodNorm,pNorm); |
---|
313 | assert(ApproxEqual(Dist,40)&&ApproxEqual(norm,vmz)&&*pgoodNorm); |
---|
314 | Dist=trd1.DistanceToOut(pzero,vxy,calcNorm,pgoodNorm,pNorm); |
---|
315 | assert(ApproxEqual(Dist,std::sqrt(800.))&&*pgoodNorm); |
---|
316 | |
---|
317 | Dist=trd1.DistanceToOut(ponxside,vx,calcNorm,pgoodNorm,pNorm); |
---|
318 | assert(ApproxEqual(Dist,0)&&ApproxEqual(*pNorm,vx)&&*pgoodNorm); |
---|
319 | Dist=trd1.DistanceToOut(ponmxside,vmx,calcNorm,pgoodNorm,pNorm); |
---|
320 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vmx)&&*pgoodNorm); |
---|
321 | Dist=trd1.DistanceToOut(ponyside,vy,calcNorm,pgoodNorm,pNorm); |
---|
322 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vy)&&*pgoodNorm); |
---|
323 | Dist=trd1.DistanceToOut(ponmyside,vmy,calcNorm,pgoodNorm,pNorm); |
---|
324 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vmy)&&*pgoodNorm); |
---|
325 | Dist=trd1.DistanceToOut(ponzside,vz,calcNorm,pgoodNorm,pNorm); |
---|
326 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vz)&&*pgoodNorm); |
---|
327 | Dist=trd1.DistanceToOut(ponmzside,vmz,calcNorm,pgoodNorm,pNorm); |
---|
328 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vmz)&&*pgoodNorm); |
---|
329 | |
---|
330 | Dist=trd2.DistanceToOut(pzero,vx,calcNorm,pgoodNorm,pNorm); |
---|
331 | assert(ApproxEqual(Dist,20)&&ApproxEqual(*pNorm,G4ThreeVector(cosa,0,-sina))&&*pgoodNorm); |
---|
332 | Dist=trd2.DistanceToOut(pzero,vmx,calcNorm,pgoodNorm,pNorm); |
---|
333 | assert(ApproxEqual(Dist,20)&&ApproxEqual(norm,G4ThreeVector(-cosa,0,-sina))&&*pgoodNorm); |
---|
334 | Dist=trd2.DistanceToOut(pzero,vy,calcNorm,pgoodNorm,pNorm); |
---|
335 | assert(ApproxEqual(Dist,30)&&ApproxEqual(norm,G4ThreeVector(0,cosa,-sina))&&*pgoodNorm); |
---|
336 | Dist=trd2.DistanceToOut(pzero,vmy,calcNorm,pgoodNorm,pNorm); |
---|
337 | assert(ApproxEqual(Dist,30)&&ApproxEqual(norm,G4ThreeVector(0,-cosa,-sina))&&*pgoodNorm); |
---|
338 | Dist=trd2.DistanceToOut(pzero,vz,calcNorm,pgoodNorm,pNorm); |
---|
339 | assert(ApproxEqual(Dist,40)&&ApproxEqual(norm,vz)&&*pgoodNorm); |
---|
340 | Dist=trd2.DistanceToOut(pzero,vmz,calcNorm,pgoodNorm,pNorm); |
---|
341 | assert(ApproxEqual(Dist,40)&&ApproxEqual(norm,vmz)&&*pgoodNorm); |
---|
342 | Dist=trd2.DistanceToOut(pzero,vxy,calcNorm,pgoodNorm,pNorm); |
---|
343 | assert(ApproxEqual(Dist,std::sqrt(800.))&&*pgoodNorm); |
---|
344 | |
---|
345 | Dist=trd2.DistanceToOut(ponxside,vx,calcNorm,pgoodNorm,pNorm); |
---|
346 | assert(ApproxEqual(Dist,0)&&ApproxEqual(*pNorm,G4ThreeVector(cosa,0,-sina))&&*pgoodNorm); |
---|
347 | Dist=trd2.DistanceToOut(ponmxside,vmx,calcNorm,pgoodNorm,pNorm); |
---|
348 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,G4ThreeVector(-cosa,0,-sina))&&*pgoodNorm); |
---|
349 | Dist=trd2.DistanceToOut(ponyside,vy,calcNorm,pgoodNorm,pNorm); |
---|
350 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,G4ThreeVector(0,cosa,-sina))&&*pgoodNorm); |
---|
351 | Dist=trd2.DistanceToOut(ponmyside,vmy,calcNorm,pgoodNorm,pNorm); |
---|
352 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,G4ThreeVector(0,-cosa,-sina))&&*pgoodNorm); |
---|
353 | Dist=trd2.DistanceToOut(ponzside,vz,calcNorm,pgoodNorm,pNorm); |
---|
354 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vz)&&*pgoodNorm); |
---|
355 | Dist=trd2.DistanceToOut(ponmzside,vmz,calcNorm,pgoodNorm,pNorm); |
---|
356 | assert(ApproxEqual(Dist,0)&&ApproxEqual(norm,vmz)&&*pgoodNorm); |
---|
357 | |
---|
358 | |
---|
359 | //DistanceToIn(P) |
---|
360 | |
---|
361 | Dist=trd1.DistanceToIn(pbigx); |
---|
362 | assert(ApproxEqual(Dist,80)); |
---|
363 | Dist=trd1.DistanceToIn(pbigmx); |
---|
364 | assert(ApproxEqual(Dist,80)); |
---|
365 | Dist=trd1.DistanceToIn(pbigy); |
---|
366 | assert(ApproxEqual(Dist,70)); |
---|
367 | Dist=trd1.DistanceToIn(pbigmy); |
---|
368 | assert(ApproxEqual(Dist,70)); |
---|
369 | Dist=trd1.DistanceToIn(pbigz); |
---|
370 | assert(ApproxEqual(Dist,60)); |
---|
371 | Dist=trd1.DistanceToIn(pbigmz); |
---|
372 | assert(ApproxEqual(Dist,60)); |
---|
373 | |
---|
374 | Dist=trd2.DistanceToIn(pbigx); |
---|
375 | assert(ApproxEqual(Dist,80*cosa)); |
---|
376 | Dist=trd2.DistanceToIn(pbigmx); |
---|
377 | assert(ApproxEqual(Dist,80*cosa)); |
---|
378 | Dist=trd2.DistanceToIn(pbigy); |
---|
379 | assert(ApproxEqual(Dist,70*cosa)); |
---|
380 | Dist=trd2.DistanceToIn(pbigmy); |
---|
381 | assert(ApproxEqual(Dist,70*cosa)); |
---|
382 | Dist=trd2.DistanceToIn(pbigz); |
---|
383 | assert(ApproxEqual(Dist,60)); |
---|
384 | Dist=trd2.DistanceToIn(pbigmz); |
---|
385 | assert(ApproxEqual(Dist,60)); |
---|
386 | |
---|
387 | |
---|
388 | // DistanceToIn(P,V) |
---|
389 | |
---|
390 | Dist=trd1.DistanceToIn(pbigx,vmx); |
---|
391 | assert(ApproxEqual(Dist,80)); |
---|
392 | Dist=trd1.DistanceToIn(pbigmx,vx); |
---|
393 | assert(ApproxEqual(Dist,80)); |
---|
394 | Dist=trd1.DistanceToIn(pbigy,vmy); |
---|
395 | assert(ApproxEqual(Dist,70)); |
---|
396 | Dist=trd1.DistanceToIn(pbigmy,vy); |
---|
397 | assert(ApproxEqual(Dist,70)); |
---|
398 | Dist=trd1.DistanceToIn(pbigz,vmz); |
---|
399 | assert(ApproxEqual(Dist,60)); |
---|
400 | Dist=trd1.DistanceToIn(pbigmz,vz); |
---|
401 | assert(ApproxEqual(Dist,60)); |
---|
402 | Dist=trd1.DistanceToIn(pbigx,vxy); |
---|
403 | assert(ApproxEqual(Dist,kInfinity)); |
---|
404 | Dist=trd1.DistanceToIn(pbigmx,vxy); |
---|
405 | assert(ApproxEqual(Dist,kInfinity)); |
---|
406 | |
---|
407 | Dist=trd2.DistanceToIn(pbigx,vmx); |
---|
408 | assert(ApproxEqual(Dist,80)); |
---|
409 | Dist=trd2.DistanceToIn(pbigmx,vx); |
---|
410 | assert(ApproxEqual(Dist,80)); |
---|
411 | Dist=trd2.DistanceToIn(pbigy,vmy); |
---|
412 | assert(ApproxEqual(Dist,70)); |
---|
413 | Dist=trd2.DistanceToIn(pbigmy,vy); |
---|
414 | assert(ApproxEqual(Dist,70)); |
---|
415 | Dist=trd2.DistanceToIn(pbigz,vmz); |
---|
416 | assert(ApproxEqual(Dist,60)); |
---|
417 | Dist=trd2.DistanceToIn(pbigmz,vz); |
---|
418 | assert(ApproxEqual(Dist,60)); |
---|
419 | Dist=trd2.DistanceToIn(pbigx,vxy); |
---|
420 | assert(ApproxEqual(Dist,kInfinity)); |
---|
421 | Dist=trd2.DistanceToIn(pbigmx,vxy); |
---|
422 | assert(ApproxEqual(Dist,kInfinity)); |
---|
423 | |
---|
424 | Dist=trd3.DistanceToIn(G4ThreeVector( 0.15000000000000185, |
---|
425 | -22.048743592955137, |
---|
426 | 2.4268539333219472), |
---|
427 | G4ThreeVector(-0.76165597579890043, |
---|
428 | 0.64364445891356026, |
---|
429 | -0.074515708658524193)) ; |
---|
430 | |
---|
431 | // G4cout<<"BABAR trd distance = "<<Dist<<G4endl ; |
---|
432 | assert(ApproxEqual(Dist,0.0)); |
---|
433 | |
---|
434 | // return-value = 2.4415531753644804e-15 |
---|
435 | |
---|
436 | |
---|
437 | |
---|
438 | |
---|
439 | // CalculateExtent |
---|
440 | |
---|
441 | G4VoxelLimits limit; // Unlimited |
---|
442 | G4AffineTransform origin(pzero); |
---|
443 | G4double min,max; |
---|
444 | |
---|
445 | assert(trd1.CalculateExtent(kXAxis,limit,origin,min,max)); |
---|
446 | assert(ApproxEqual(min,-20)&&ApproxEqual(max,20)); |
---|
447 | assert(trd1.CalculateExtent(kYAxis,limit,origin,min,max)); |
---|
448 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
449 | assert(trd1.CalculateExtent(kZAxis,limit,origin,min,max)); |
---|
450 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
451 | |
---|
452 | assert(trd2.CalculateExtent(kXAxis,limit,origin,min,max)); |
---|
453 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
454 | assert(trd2.CalculateExtent(kYAxis,limit,origin,min,max)); |
---|
455 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
456 | assert(trd2.CalculateExtent(kZAxis,limit,origin,min,max)); |
---|
457 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
458 | |
---|
459 | |
---|
460 | G4ThreeVector pmxmymz(-100,-110,-120); |
---|
461 | G4AffineTransform tPosOnly(pmxmymz); |
---|
462 | |
---|
463 | assert(trd1.CalculateExtent(kXAxis,limit,tPosOnly,min,max)); |
---|
464 | assert(ApproxEqual(min,-120)&&ApproxEqual(max,-80)); |
---|
465 | assert(trd1.CalculateExtent(kYAxis,limit,tPosOnly,min,max)); |
---|
466 | assert(ApproxEqual(min,-140)&&ApproxEqual(max,-80)); |
---|
467 | assert(trd1.CalculateExtent(kZAxis,limit,tPosOnly,min,max)); |
---|
468 | assert(ApproxEqual(min,-160)&&ApproxEqual(max,-80)); |
---|
469 | |
---|
470 | assert(trd2.CalculateExtent(kXAxis,limit,tPosOnly,min,max)); |
---|
471 | assert(ApproxEqual(min,-130)&&ApproxEqual(max,-70)); |
---|
472 | assert(trd2.CalculateExtent(kYAxis,limit,tPosOnly,min,max)); |
---|
473 | assert(ApproxEqual(min,-150)&&ApproxEqual(max,-70)); |
---|
474 | assert(trd2.CalculateExtent(kZAxis,limit,tPosOnly,min,max)); |
---|
475 | assert(ApproxEqual(min,-160)&&ApproxEqual(max,-80)); |
---|
476 | |
---|
477 | |
---|
478 | G4RotationMatrix r90Z; |
---|
479 | r90Z.rotateZ(halfpi); |
---|
480 | G4AffineTransform tRotZ(r90Z,pzero); |
---|
481 | |
---|
482 | assert(trd1.CalculateExtent(kXAxis,limit,tRotZ,min,max)); |
---|
483 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
484 | assert(trd1.CalculateExtent(kYAxis,limit,tRotZ,min,max)); |
---|
485 | assert(ApproxEqual(min,-20)&&ApproxEqual(max,20)); |
---|
486 | assert(trd1.CalculateExtent(kZAxis,limit,tRotZ,min,max)); |
---|
487 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
488 | |
---|
489 | assert(trd2.CalculateExtent(kXAxis,limit,tRotZ,min,max)); |
---|
490 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
491 | assert(trd2.CalculateExtent(kYAxis,limit,tRotZ,min,max)); |
---|
492 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
493 | assert(trd2.CalculateExtent(kZAxis,limit,tRotZ,min,max)); |
---|
494 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
495 | |
---|
496 | |
---|
497 | // Check that clipped away |
---|
498 | |
---|
499 | G4VoxelLimits xClip; |
---|
500 | xClip.AddLimit(kXAxis,-100,-50); |
---|
501 | |
---|
502 | assert(!trd1.CalculateExtent(kXAxis,xClip,origin,min,max)); |
---|
503 | |
---|
504 | assert(!trd2.CalculateExtent(kXAxis,xClip,origin,min,max)); |
---|
505 | |
---|
506 | |
---|
507 | // Assert clipped to volume |
---|
508 | |
---|
509 | G4VoxelLimits allClip; |
---|
510 | allClip.AddLimit(kXAxis,-5,+5); |
---|
511 | allClip.AddLimit(kYAxis,-5,+5); |
---|
512 | allClip.AddLimit(kZAxis,-5,+5); |
---|
513 | |
---|
514 | G4RotationMatrix genRot; |
---|
515 | genRot.rotateX(pi/6); |
---|
516 | genRot.rotateY(pi/6); |
---|
517 | genRot.rotateZ(pi/6); |
---|
518 | |
---|
519 | G4AffineTransform tGen(genRot,vx); |
---|
520 | |
---|
521 | assert(trd1.CalculateExtent(kXAxis,allClip,tGen,min,max)); |
---|
522 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
523 | assert(trd1.CalculateExtent(kYAxis,allClip,tGen,min,max)); |
---|
524 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
525 | assert(trd1.CalculateExtent(kZAxis,allClip,tGen,min,max)); |
---|
526 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
527 | |
---|
528 | assert(trd2.CalculateExtent(kXAxis,allClip,tGen,min,max)); |
---|
529 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
530 | assert(trd2.CalculateExtent(kYAxis,allClip,tGen,min,max)); |
---|
531 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
532 | assert(trd2.CalculateExtent(kZAxis,allClip,tGen,min,max)); |
---|
533 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
534 | |
---|
535 | |
---|
536 | G4VoxelLimits buggyClip2; |
---|
537 | buggyClip2.AddLimit(kXAxis,5,15); |
---|
538 | |
---|
539 | assert(trd1.CalculateExtent(kXAxis,buggyClip2,origin,min,max)); |
---|
540 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
541 | assert(trd1.CalculateExtent(kYAxis,buggyClip2,origin,min,max)); |
---|
542 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
543 | assert(trd1.CalculateExtent(kZAxis,buggyClip2,origin,min,max)); |
---|
544 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
545 | |
---|
546 | assert(trd2.CalculateExtent(kXAxis,buggyClip2,origin,min,max)); |
---|
547 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
548 | assert(trd2.CalculateExtent(kYAxis,buggyClip2,origin,min,max)); |
---|
549 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
550 | assert(trd2.CalculateExtent(kZAxis,buggyClip2,origin,min,max)); |
---|
551 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
552 | |
---|
553 | |
---|
554 | buggyClip2.AddLimit(kYAxis,5,15); |
---|
555 | |
---|
556 | assert(trd1.CalculateExtent(kXAxis,buggyClip2,origin,min,max)); |
---|
557 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
558 | assert(trd1.CalculateExtent(kYAxis,buggyClip2,origin,min,max)); |
---|
559 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
560 | assert(trd1.CalculateExtent(kZAxis,buggyClip2,origin,min,max)); |
---|
561 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
562 | |
---|
563 | assert(trd2.CalculateExtent(kXAxis,buggyClip2,origin,min,max)); |
---|
564 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
565 | assert(trd2.CalculateExtent(kYAxis,buggyClip2,origin,min,max)); |
---|
566 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
567 | assert(trd2.CalculateExtent(kZAxis,buggyClip2,origin,min,max)); |
---|
568 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
569 | |
---|
570 | |
---|
571 | G4VoxelLimits buggyClip1; |
---|
572 | buggyClip1.AddLimit(kXAxis,-5,+5); |
---|
573 | |
---|
574 | assert(trd1.CalculateExtent(kXAxis,buggyClip1,origin,min,max)); |
---|
575 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
576 | assert(trd1.CalculateExtent(kYAxis,buggyClip1,origin,min,max)); |
---|
577 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
578 | assert(trd1.CalculateExtent(kZAxis,buggyClip1,origin,min,max)); |
---|
579 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
580 | |
---|
581 | assert(trd2.CalculateExtent(kXAxis,buggyClip1,origin,min,max)); |
---|
582 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
583 | assert(trd2.CalculateExtent(kYAxis,buggyClip1,origin,min,max)); |
---|
584 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
585 | assert(trd2.CalculateExtent(kZAxis,buggyClip1,origin,min,max)); |
---|
586 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
587 | |
---|
588 | |
---|
589 | buggyClip1.AddLimit(kYAxis,-5,+5); |
---|
590 | |
---|
591 | assert(trd1.CalculateExtent(kXAxis,buggyClip1,origin,min,max)); |
---|
592 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
593 | assert(trd1.CalculateExtent(kYAxis,buggyClip1,origin,min,max)); |
---|
594 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
595 | assert(trd1.CalculateExtent(kZAxis,buggyClip1,origin,min,max)); |
---|
596 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
597 | |
---|
598 | assert(trd2.CalculateExtent(kXAxis,buggyClip1,origin,min,max)); |
---|
599 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
600 | assert(trd2.CalculateExtent(kYAxis,buggyClip1,origin,min,max)); |
---|
601 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
602 | assert(trd2.CalculateExtent(kZAxis,buggyClip1,origin,min,max)); |
---|
603 | assert(ApproxEqual(min,-40)&&ApproxEqual(max,40)); |
---|
604 | |
---|
605 | |
---|
606 | G4VoxelLimits newvoxlim; |
---|
607 | newvoxlim.AddLimit(kZAxis,-5,+5); |
---|
608 | |
---|
609 | assert(trd1.CalculateExtent(kXAxis,newvoxlim,origin,min,max)); |
---|
610 | assert(ApproxEqual(min,-20)&&ApproxEqual(max,20)); |
---|
611 | assert(trd1.CalculateExtent(kYAxis,newvoxlim,origin,min,max)); |
---|
612 | assert(ApproxEqual(min,-30)&&ApproxEqual(max,30)); |
---|
613 | assert(trd1.CalculateExtent(kZAxis,newvoxlim,origin,min,max)); |
---|
614 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
615 | |
---|
616 | assert(trd2.CalculateExtent(kXAxis,newvoxlim,origin,min,max)); |
---|
617 | assert(ApproxEqual(min,-(20+5*tanga))&&ApproxEqual(max,20+5*tanga)); |
---|
618 | assert(trd2.CalculateExtent(kYAxis,newvoxlim,origin,min,max)); |
---|
619 | assert(ApproxEqual(min,-(30+5*tanga))&&ApproxEqual(max,30+5*tanga)); |
---|
620 | assert(trd2.CalculateExtent(kZAxis,newvoxlim,origin,min,max)); |
---|
621 | assert(ApproxEqual(min,-5)&&ApproxEqual(max,5)); |
---|
622 | |
---|
623 | G4VoxelLimits nonsymvox; |
---|
624 | nonsymvox.AddLimit(kZAxis,5,15); |
---|
625 | |
---|
626 | assert(trd2.CalculateExtent(kXAxis,nonsymvox,origin,min,max)); |
---|
627 | assert(ApproxEqual(min,-(20+15*tanga))&&ApproxEqual(max,20+15*tanga)); |
---|
628 | assert(trd2.CalculateExtent(kYAxis,nonsymvox,origin,min,max)); |
---|
629 | assert(ApproxEqual(min,-(30+15*tanga))&&ApproxEqual(max,30+15*tanga)); |
---|
630 | assert(trd2.CalculateExtent(kZAxis,nonsymvox,origin,min,max)); |
---|
631 | assert(ApproxEqual(min,5)&&ApproxEqual(max,15)); |
---|
632 | |
---|
633 | return true; |
---|
634 | } |
---|
635 | |
---|
636 | int main() |
---|
637 | { |
---|
638 | #ifdef NDEBUG |
---|
639 | G4Exception("FAIL: *** Assertions must be compiled in! ***"); |
---|
640 | #endif |
---|
641 | assert(testG4Trd()); |
---|
642 | return 0; |
---|
643 | } |
---|
644 | |
---|
645 | |
---|
646 | |
---|
647 | |
---|
648 | |
---|