source: trunk/source/geometry/volumes/test/testG4ReplicaNavigation.cc@ 1350

Last change on this file since 1350 was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

File size: 13.9 KB
RevLine 
[1316]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: testG4ReplicaNavigation.cc,v 1.11 2006/06/29 18:58:48 gunter Exp $
[1347]28// GEANT4 tag $Name: geant4-09-04-ref-00 $
[1316]29//
30//
31// Test private location & distance computation functions of
32// G4ReplicaNavigation Paul Kent Aug 96
33
34
35#include <assert.h>
36#include "ApproxEqual.hh"
37#include "globals.hh"
38#include "G4Box.hh"
39#include "G4Sphere.hh"
40#include "G4LogicalVolume.hh"
41#include "G4ReplicaNavigation.hh"
42#include "G4PVReplica.hh"
43#include "G4PVPlacement.hh"
44
45class G4ReplicaNavigationTester
46{
47public:
48 EInside Inside(const G4VPhysicalVolume *pVol,
49 const G4int replicaNo,
50 const G4ThreeVector &localPoint) const
51 {
52 return nav.Inside(pVol,replicaNo,localPoint);
53 }
54
55 G4double DistanceToOut(const G4VPhysicalVolume *pVol,
56 const G4int replicaNo,
57 const G4ThreeVector &localPoint) const
58 {
59 return nav.DistanceToOut(pVol,replicaNo,localPoint);
60 }
61
62 G4double DistanceToOut(const G4VPhysicalVolume *pVol,
63 const G4int replicaNo,
64 const G4ThreeVector &localPoint,
65 const G4ThreeVector &localDirection) const
66 {
67 return nav.DistanceToOut(pVol,replicaNo,localPoint,localDirection);
68 }
69
70private:
71 G4ReplicaNavigation nav;
72};
73
74G4bool testG4ReplicaNavigation()
75{
76 EInside in;
77 G4double Dist;
78
79 // Define two worlds
80 //
81 G4Box* hall_box =
82 new G4Box("expHall_box",3000,3000,3000);
83 G4LogicalVolume* hall_log1 =
84 new G4LogicalVolume(hall_box,0,"expHall_log",0,0,0);
85 G4VPhysicalVolume* hall_phys1 =
86 new G4PVPlacement(0,G4ThreeVector(),"expHall1",hall_log1,0,false,0);
87 G4LogicalVolume* hall_log2 =
88 new G4LogicalVolume(hall_box,0,"expHall_log",0,0,0);
89 G4VPhysicalVolume* hall_phys2 =
90 new G4PVPlacement(0,G4ThreeVector(),"expHall2",hall_log2,0,false,0);
91
92 // Define volumes to be sliced
93 //
94 G4Box* fBox = new G4Box("Test Box",120.,120.,120.);
95 G4LogicalVolume *pMotherVol1X= new G4LogicalVolume(fBox, 0, "lmoth1X", 0, 0, 0);
96 new G4PVPlacement(0,G4ThreeVector(),"pmoth1",pMotherVol1X,hall_phys1,false,0);
97 G4LogicalVolume *pMotherVol1Y= new G4LogicalVolume(fBox, 0, "lmoth1Y", 0, 0, 0);
98 new G4PVPlacement(0,G4ThreeVector(),"pmoth1",pMotherVol1Y,hall_phys1,false,0);
99 G4LogicalVolume *pMotherVol1Z= new G4LogicalVolume(fBox, 0, "lmoth1Z", 0, 0, 0);
100 new G4PVPlacement(0,G4ThreeVector(),"pmoth1",pMotherVol1Z,hall_phys1,false,0);
101
102 G4Sphere* fSphere = new G4Sphere("Test Sphere",0.,80.,0*deg,360*deg,0*deg,360*deg);
103 G4LogicalVolume *pMotherVol2P= new G4LogicalVolume(fSphere, 0, "lmoth2P", 0, 0, 0);
104 new G4PVPlacement(0,G4ThreeVector(),"pmoth2",pMotherVol2P,hall_phys2,false,0);
105 G4LogicalVolume *pMotherVol2R= new G4LogicalVolume(fSphere, 0, "lmoth2R", 0, 0, 0);
106 new G4PVPlacement(0,G4ThreeVector(),"pmoth2",pMotherVol2R,hall_phys2,false,0);
107
108 G4ReplicaNavigationTester repNav;
109
110 // Define the actual slices (cartesian axis)
111 //
112 G4Box* xBoxSlice = new G4Box("Sliced Box X",40.,120.,120.);
113 G4LogicalVolume* xBoxLog= new G4LogicalVolume(xBoxSlice, 0, "xBoxSlice", 0, 0, 0);
114 G4PVReplica xRep("TestX",xBoxLog,pMotherVol1X,kXAxis,3,40);
115
116 G4Box* yBoxSlice = new G4Box("Sliced Box Y",120.,40.,120.);
117 G4LogicalVolume* yBoxLog= new G4LogicalVolume(yBoxSlice, 0, "yBoxSlice", 0, 0, 0);
118 G4PVReplica yRep("TestY",yBoxLog,pMotherVol1Y,kYAxis,3,40);
119
120 G4Box* zBoxSlice = new G4Box("Sliced Box Z",120.,120.,40.);
121 G4LogicalVolume* zBoxLog= new G4LogicalVolume(zBoxSlice, 0, "zBoxSlice", 0, 0, 0);
122 G4PVReplica zRep("TestZ",zBoxLog,pMotherVol1Z,kZAxis,3,40);
123
124 // Define the actual slices (Phi and Rho)
125 //
126 G4Sphere* fSphereP = new G4Sphere("Sliced Sphere Phi",0.,80.,0*deg,90*deg,0*deg,360*deg);
127 G4LogicalVolume* phiSphereLog= new G4LogicalVolume(fSphereP, 0, "PhiSlice", 0, 0, 0);
128 G4PVReplica phiRep("TestPhi",phiSphereLog,pMotherVol2P,kPhi,4,pi*0.5);
129 G4Sphere* fSphereR = new G4Sphere("Sliced Sphere Rho",0.,20.,0*deg,360*deg,0*deg,360*deg);
130 G4LogicalVolume* rhoSphereLog= new G4LogicalVolume(fSphereR, 0, "RhoSlice", 0, 0, 0);
131 G4PVReplica radRep("TestRho",rhoSphereLog,pMotherVol2R,kRho,4,20);
132
133 in=repNav.Inside(&xRep,0,G4ThreeVector(21,0,0));
134 assert(in==kOutside);
135 in=repNav.Inside(&xRep,0,G4ThreeVector(20,0,0));
136 assert(in==kSurface);
137 in=repNav.Inside(&xRep,0,G4ThreeVector(19,0,0));
138 assert(in==kInside);
139 in=repNav.Inside(&xRep,0,G4ThreeVector(-20,0,0));
140 assert(in==kSurface);
141 in=repNav.Inside(&xRep,0,G4ThreeVector(-21,0,0));
142 assert(in==kOutside);
143
144 in=repNav.Inside(&yRep,0,G4ThreeVector(0,21,0));
145 assert(in==kOutside);
146 in=repNav.Inside(&yRep,0,G4ThreeVector(0,20,0));
147 assert(in==kSurface);
148 in=repNav.Inside(&yRep,0,G4ThreeVector(0,19,0));
149 assert(in==kInside);
150 in=repNav.Inside(&yRep,0,G4ThreeVector(0,-20,0));
151 assert(in==kSurface);
152 in=repNav.Inside(&yRep,0,G4ThreeVector(0,-21,0));
153 assert(in==kOutside);
154
155 in=repNav.Inside(&zRep,0,G4ThreeVector(0,0,21));
156 assert(in==kOutside);
157 in=repNav.Inside(&zRep,0,G4ThreeVector(0,0,20));
158 assert(in==kSurface);
159 in=repNav.Inside(&zRep,0,G4ThreeVector(0,0,19));
160 assert(in==kInside);
161 in=repNav.Inside(&zRep,0,G4ThreeVector(0,0,-20));
162 assert(in==kSurface);
163 in=repNav.Inside(&zRep,0,G4ThreeVector(0,0,-21));
164 assert(in==kOutside);
165
166 in=repNav.Inside(&phiRep,0,G4ThreeVector(0,0,0));
167 assert(in==kSurface);
168 in=repNav.Inside(&phiRep,0,G4ThreeVector(10,0,0));
169 assert(in==kInside);
170 in=repNav.Inside(&phiRep,0,G4ThreeVector(-10,0,0));
171 assert(in==kOutside);
172 in=repNav.Inside(&phiRep,0,G4ThreeVector(10,10,0));
173 assert(in==kSurface);
174 in=repNav.Inside(&phiRep,0,G4ThreeVector(10,10.1,0));
175 assert(in==kOutside);
176 in=repNav.Inside(&phiRep,0,G4ThreeVector(10,-10,0));
177 assert(in==kSurface);
178 in=repNav.Inside(&phiRep,0,G4ThreeVector(10,-10.1,0));
179 assert(in==kOutside);
180
181 in=repNav.Inside(&radRep,0,G4ThreeVector(0,0,0));
182 assert(in==kInside);
183 in=repNav.Inside(&radRep,0,G4ThreeVector(0,20,0));
184 assert(in==kSurface);
185 in=repNav.Inside(&radRep,0,G4ThreeVector(0,21,0));
186 assert(in==kOutside);
187
188 in=repNav.Inside(&radRep,1,G4ThreeVector(0,0,0));
189 assert(in==kOutside);
190 in=repNav.Inside(&radRep,1,G4ThreeVector(0,20,0));
191 assert(in==kSurface);
192 in=repNav.Inside(&radRep,1,G4ThreeVector(0,30,0));
193 assert(in==kInside);
194
195 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(0,0,0));
196 assert(ApproxEqual(Dist,20));
197 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(20,20,20));
198 assert(ApproxEqual(Dist,0));
199 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(-21,-21,-21));
200 assert(Dist==0);
201 Dist=repNav.DistanceToOut(&yRep,0,G4ThreeVector(0,0,0));
202 assert(ApproxEqual(Dist,20));
203 Dist=repNav.DistanceToOut(&yRep,0,G4ThreeVector(20,20,20));
204 assert(ApproxEqual(Dist,0));
205 Dist=repNav.DistanceToOut(&yRep,0,G4ThreeVector(-21,-21,-21));
206 assert(Dist==0);
207 Dist=repNav.DistanceToOut(&zRep,0,G4ThreeVector(0,0,0));
208 assert(ApproxEqual(Dist,20));
209 Dist=repNav.DistanceToOut(&zRep,0,G4ThreeVector(20,20,20));
210 assert(ApproxEqual(Dist,0));
211 Dist=repNav.DistanceToOut(&zRep,0,G4ThreeVector(-21,-21,-21));
212 assert(Dist==0);
213
214
215 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,0,0));
216 assert(ApproxEqual(Dist,0));
217 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0));
218 assert(ApproxEqual(Dist,10*std::sin(pi*0.25)));
219 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(-10,0,0));
220 assert(Dist==0);
221 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,10,0));
222 assert(ApproxEqual(Dist,0));
223 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,-10,0));
224 assert(ApproxEqual(Dist,0));
225 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,5,0));
226 assert(ApproxEqual(Dist,std::sqrt(125.)*std::sin(pi*0.25-std::atan(0.5))));
227
228 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,0,0));
229 assert(ApproxEqual(Dist,20));
230 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,20,0));
231 assert(ApproxEqual(Dist,0));
232 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,21,0));
233 assert(Dist==0);
234
235 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(0,0,0));
236 assert(Dist==0);
237 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(0,20,0));
238 assert(ApproxEqual(Dist,0));
239 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(0,21,0));
240 assert(Dist==1);
241 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(21,21,0));
242 std::cout.precision(8);
243 // G4cout << " Dist is " << Dist << " and expected= " << std::sqrt(2.*441.)-20. << G4endl;
244 // G4cout << " a difference of " << Dist-(std::sqrt(2.*441.)-20.) << G4endl;
245 assert( Dist - (std::sqrt(2.*441.)-20.) < 1.e-14 );
246 // assert(ApproxEqual(Dist, std::sqrt(2.*441.)-20.));
247
248
249 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(0,0,0),
250 G4ThreeVector(1,0,0));
251 assert(ApproxEqual(Dist,20));
252 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(0,0,0),
253 G4ThreeVector(-1,0,0));
254 assert(ApproxEqual(Dist,20));
255 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(20,0,0),
256 G4ThreeVector(1,0,0));
257 assert(ApproxEqual(Dist,0));
258 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(20,0,0),
259 G4ThreeVector(-1,0,0));
260 assert(ApproxEqual(Dist,40));
261 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(21,0,0),
262 G4ThreeVector(1,0,0));
263 assert(Dist==0);
264 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(20,0,0),
265 G4ThreeVector(-1/std::sqrt(2.),-1/std::sqrt(2.),0));
266 assert(ApproxEqual(Dist,40*std::sqrt(2.)));
267 Dist=repNav.DistanceToOut(&xRep,0,G4ThreeVector(20,0,0),
268 G4ThreeVector(0,1,0));
269 assert(Dist==kInfinity);
270
271
272 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,0,0),
273 G4ThreeVector(1,0,0));
274 assert(Dist==kInfinity);
275 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(-1,0,0),
276 G4ThreeVector(1,0,0));
277 assert(Dist==kInfinity);
278 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,-1,0),
279 G4ThreeVector(1,0,0));
280 assert(Dist==kInfinity);
281 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,1,0),
282 G4ThreeVector(1,0,0));
283 assert(Dist==kInfinity);
284 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(-1,0,0),
285 G4ThreeVector(-1,0,0));
286 assert(Dist==0);
287 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,-1,0),
288 G4ThreeVector(-1,0,0));
289// assert(Dist==0);
290 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,1,0),
291 G4ThreeVector(-1,0,0));
292 assert(Dist==0);
293 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(0,0,0),
294 G4ThreeVector(-1,0,0));
295 assert(ApproxEqual(Dist,0));
296 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0),
297 G4ThreeVector(-1,0,0));
298 assert(ApproxEqual(Dist,10));
299 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0),
300 G4ThreeVector(0,1,0));
301 assert(ApproxEqual(Dist,10));
302 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0),
303 G4ThreeVector(0,-1,0));
304 assert(ApproxEqual(Dist,10));
305 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0),
306 G4ThreeVector(-1/std::sqrt(2.),1/std::sqrt(2.),0));
307 assert(ApproxEqual(Dist,10*std::sin(pi*0.25)));
308 Dist=repNav.DistanceToOut(&phiRep,0,G4ThreeVector(10,0,0),
309 G4ThreeVector(-1/std::sqrt(2.),-1/std::sqrt(2.),0));
310 assert(ApproxEqual(Dist,10*std::sin(pi*0.25)));
311
312
313 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,0,0),
314 G4ThreeVector(1,0,0));
315 assert(ApproxEqual(Dist,20));
316 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,0,0),
317 G4ThreeVector(-1,0,0));
318 assert(ApproxEqual(Dist,20));
319 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(0,0,0),
320 G4ThreeVector(-1/std::sqrt(2.),-1/std::sqrt(2.),0));
321 assert(ApproxEqual(Dist,20));
322 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(std::sqrt(200.),std::sqrt(200.),0),
323 G4ThreeVector(-1/std::sqrt(2.),-1/std::sqrt(2.),0));
324 assert(ApproxEqual(Dist,40));
325 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(std::sqrt(200.),std::sqrt(200.),0),
326 G4ThreeVector(1/std::sqrt(2.),1/std::sqrt(2.),0));
327 assert(ApproxEqual(Dist,0));
328 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(std::sqrt(200.),std::sqrt(200.),0),
329 G4ThreeVector(0,0,1));
330 assert(Dist==kInfinity);
331 Dist=repNav.DistanceToOut(&radRep,0,G4ThreeVector(21,0,0),
332 G4ThreeVector(1,0,0));
333 assert(Dist==0);
334 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(20,0,0),
335 G4ThreeVector(1,0,0));
336 assert(ApproxEqual(Dist,20));
337 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(20,0,0),
338 G4ThreeVector(-1,0,0));
339 assert(ApproxEqual(Dist,0));
340 Dist=repNav.DistanceToOut(&radRep,1,G4ThreeVector(20,0,0),
341 G4ThreeVector(0,-1,0));
342 assert(ApproxEqual(Dist,std::sqrt(40.*40.-20.*20.)));
343
344
345 return true;
346}
347
348int main()
349{
350#ifdef NDEBUG
351 G4Exception("FAIL: *** Assertions must be compiled in! ***");
352#endif
353 assert(testG4ReplicaNavigation());
354 return 0;
355}
356
Note: See TracBrowser for help on using the repository browser.