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: G4Polyhedra.cc,v 1.42 2008/05/15 13:45:15 gcosmo Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
29 | // |
---|
30 | // |
---|
31 | // -------------------------------------------------------------------- |
---|
32 | // GEANT 4 class source file |
---|
33 | // |
---|
34 | // |
---|
35 | // G4Polyhedra.cc |
---|
36 | // |
---|
37 | // Implementation of a CSG polyhedra, as an inherited class of G4VCSGfaceted. |
---|
38 | // |
---|
39 | // To be done: |
---|
40 | // * Cracks: there are probably small cracks in the seams between the |
---|
41 | // phi face (G4PolyPhiFace) and sides (G4PolyhedraSide) that are not |
---|
42 | // entirely leakproof. Also, I am not sure all vertices are leak proof. |
---|
43 | // * Many optimizations are possible, but not implemented. |
---|
44 | // * Visualization needs to be updated outside of this routine. |
---|
45 | // |
---|
46 | // Utility classes: |
---|
47 | // * G4EnclosingCylinder: I decided a quick check of geometry would be a |
---|
48 | // good idea (for CPU speed). If the quick check fails, the regular |
---|
49 | // full-blown G4VCSGfaceted version is invoked. |
---|
50 | // * G4ReduciblePolygon: Really meant as a check of input parameters, |
---|
51 | // this utility class also "converts" the GEANT3-like PGON/PCON |
---|
52 | // arguments into the newer ones. |
---|
53 | // Both these classes are implemented outside this file because they are |
---|
54 | // shared with G4Polycone. |
---|
55 | // |
---|
56 | // -------------------------------------------------------------------- |
---|
57 | |
---|
58 | #include "G4Polyhedra.hh" |
---|
59 | |
---|
60 | #include "G4PolyhedraSide.hh" |
---|
61 | #include "G4PolyPhiFace.hh" |
---|
62 | |
---|
63 | #include "Randomize.hh" |
---|
64 | |
---|
65 | #include "G4Polyhedron.hh" |
---|
66 | #include "G4EnclosingCylinder.hh" |
---|
67 | #include "G4ReduciblePolygon.hh" |
---|
68 | #include "G4VPVParameterisation.hh" |
---|
69 | |
---|
70 | #include <sstream> |
---|
71 | |
---|
72 | using namespace CLHEP; |
---|
73 | |
---|
74 | // |
---|
75 | // Constructor (GEANT3 style parameters) |
---|
76 | // |
---|
77 | // GEANT3 PGON radii are specified in the distance to the norm of each face. |
---|
78 | // |
---|
79 | G4Polyhedra::G4Polyhedra( const G4String& name, |
---|
80 | G4double phiStart, |
---|
81 | G4double thePhiTotal, |
---|
82 | G4int theNumSide, |
---|
83 | G4int numZPlanes, |
---|
84 | const G4double zPlane[], |
---|
85 | const G4double rInner[], |
---|
86 | const G4double rOuter[] ) |
---|
87 | : G4VCSGfaceted( name ), genericPgon(false) |
---|
88 | { |
---|
89 | if (theNumSide <= 0) |
---|
90 | { |
---|
91 | G4cerr << "ERROR - G4Polyhedra::G4Polyhedra(): " << GetName() << G4endl |
---|
92 | << " No sides specified !" |
---|
93 | << G4endl; |
---|
94 | G4Exception("G4Polyhedra::G4Polyhedra()", "InvalidSetup", |
---|
95 | FatalException, "Solid must have at least one side."); |
---|
96 | } |
---|
97 | |
---|
98 | // |
---|
99 | // Calculate conversion factor from G3 radius to G4 radius |
---|
100 | // |
---|
101 | G4double phiTotal = thePhiTotal; |
---|
102 | if ( (phiTotal <=0) || (phiTotal >= twopi*(1-DBL_EPSILON)) ) |
---|
103 | { phiTotal = twopi; } |
---|
104 | G4double convertRad = std::cos(0.5*phiTotal/theNumSide); |
---|
105 | |
---|
106 | // |
---|
107 | // Some historical stuff |
---|
108 | // |
---|
109 | original_parameters = new G4PolyhedraHistorical; |
---|
110 | |
---|
111 | original_parameters->numSide = theNumSide; |
---|
112 | original_parameters->Start_angle = phiStart; |
---|
113 | original_parameters->Opening_angle = phiTotal; |
---|
114 | original_parameters->Num_z_planes = numZPlanes; |
---|
115 | original_parameters->Z_values = new G4double[numZPlanes]; |
---|
116 | original_parameters->Rmin = new G4double[numZPlanes]; |
---|
117 | original_parameters->Rmax = new G4double[numZPlanes]; |
---|
118 | |
---|
119 | G4int i; |
---|
120 | for (i=0; i<numZPlanes; i++) |
---|
121 | { |
---|
122 | if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] )) |
---|
123 | { |
---|
124 | if( (rInner[i] > rOuter[i+1]) |
---|
125 | ||(rInner[i+1] > rOuter[i]) ) |
---|
126 | { |
---|
127 | DumpInfo(); |
---|
128 | G4cerr << "ERROR - G4Polyhedra::G4Polyhedra()" |
---|
129 | << G4endl |
---|
130 | << " Segments are not contiguous !" << G4endl |
---|
131 | << " rMin[" << i << "] = " << rInner[i] |
---|
132 | << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl |
---|
133 | << " rMin[" << i+1 << "] = " << rInner[i+1] |
---|
134 | << " -- rMax[" << i << "] = " << rOuter[i] << G4endl; |
---|
135 | G4Exception("G4Polyhedra::G4Polyhedra()","InvalidSetup",FatalException, |
---|
136 | "Cannot create a Polyhedra with no contiguous segments."); |
---|
137 | } |
---|
138 | } |
---|
139 | original_parameters->Z_values[i] = zPlane[i]; |
---|
140 | original_parameters->Rmin[i] = rInner[i]/convertRad; |
---|
141 | original_parameters->Rmax[i] = rOuter[i]/convertRad; |
---|
142 | } |
---|
143 | |
---|
144 | |
---|
145 | // |
---|
146 | // Build RZ polygon using special PCON/PGON GEANT3 constructor |
---|
147 | // |
---|
148 | G4ReduciblePolygon *rz = |
---|
149 | new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes ); |
---|
150 | rz->ScaleA( 1/convertRad ); |
---|
151 | |
---|
152 | // |
---|
153 | // Do the real work |
---|
154 | // |
---|
155 | Create( phiStart, phiTotal, theNumSide, rz ); |
---|
156 | |
---|
157 | delete rz; |
---|
158 | } |
---|
159 | |
---|
160 | |
---|
161 | // |
---|
162 | // Constructor (generic parameters) |
---|
163 | // |
---|
164 | G4Polyhedra::G4Polyhedra( const G4String& name, |
---|
165 | G4double phiStart, |
---|
166 | G4double phiTotal, |
---|
167 | G4int theNumSide, |
---|
168 | G4int numRZ, |
---|
169 | const G4double r[], |
---|
170 | const G4double z[] ) |
---|
171 | : G4VCSGfaceted( name ), genericPgon(true) |
---|
172 | { |
---|
173 | G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ ); |
---|
174 | |
---|
175 | Create( phiStart, phiTotal, theNumSide, rz ); |
---|
176 | |
---|
177 | // Set original_parameters struct for consistency |
---|
178 | // |
---|
179 | SetOriginalParameters(); |
---|
180 | |
---|
181 | delete rz; |
---|
182 | } |
---|
183 | |
---|
184 | |
---|
185 | // |
---|
186 | // Create |
---|
187 | // |
---|
188 | // Generic create routine, called by each constructor |
---|
189 | // after conversion of arguments |
---|
190 | // |
---|
191 | void G4Polyhedra::Create( G4double phiStart, |
---|
192 | G4double phiTotal, |
---|
193 | G4int theNumSide, |
---|
194 | G4ReduciblePolygon *rz ) |
---|
195 | { |
---|
196 | // |
---|
197 | // Perform checks of rz values |
---|
198 | // |
---|
199 | if (rz->Amin() < 0.0) |
---|
200 | { |
---|
201 | G4cerr << "ERROR - G4Polyhedra::Create() " << GetName() << G4endl |
---|
202 | << " All R values must be >= 0 !" |
---|
203 | << G4endl; |
---|
204 | G4Exception("G4Polyhedra::Create()", "InvalidSetup", |
---|
205 | FatalException, "Illegal input parameters."); |
---|
206 | } |
---|
207 | |
---|
208 | G4double rzArea = rz->Area(); |
---|
209 | if (rzArea < -kCarTolerance) |
---|
210 | rz->ReverseOrder(); |
---|
211 | |
---|
212 | else if (rzArea < -kCarTolerance) |
---|
213 | { |
---|
214 | G4cerr << "ERROR - G4Polyhedra::Create() " << GetName() << G4endl |
---|
215 | << " R/Z cross section is zero or near zero: " |
---|
216 | << rzArea << G4endl; |
---|
217 | G4Exception("G4Polyhedra::Create()", "InvalidSetup", |
---|
218 | FatalException, "Illegal input parameters."); |
---|
219 | } |
---|
220 | |
---|
221 | if ( (!rz->RemoveDuplicateVertices( kCarTolerance )) |
---|
222 | || (!rz->RemoveRedundantVertices( kCarTolerance )) ) |
---|
223 | { |
---|
224 | G4cerr << "ERROR - G4Polyhedra::Create() " << GetName() << G4endl |
---|
225 | << " Too few unique R/Z values !" |
---|
226 | << G4endl; |
---|
227 | G4Exception("G4Polyhedra::Create()", "InvalidSetup", |
---|
228 | FatalException, "Illegal input parameters."); |
---|
229 | } |
---|
230 | |
---|
231 | if (rz->CrossesItself( 1/kInfinity )) |
---|
232 | { |
---|
233 | G4cerr << "ERROR - G4Polyhedra::Create() " << GetName() << G4endl |
---|
234 | << " R/Z segments cross !" |
---|
235 | << G4endl; |
---|
236 | G4Exception("G4Polyhedra::Create()", "InvalidSetup", |
---|
237 | FatalException, "Illegal input parameters."); |
---|
238 | } |
---|
239 | |
---|
240 | numCorner = rz->NumVertices(); |
---|
241 | |
---|
242 | |
---|
243 | startPhi = phiStart; |
---|
244 | while( startPhi < 0 ) startPhi += twopi; |
---|
245 | // |
---|
246 | // Phi opening? Account for some possible roundoff, and interpret |
---|
247 | // nonsense value as representing no phi opening |
---|
248 | // |
---|
249 | if ( (phiTotal <= 0) || (phiTotal > twopi*(1-DBL_EPSILON)) ) |
---|
250 | { |
---|
251 | phiIsOpen = false; |
---|
252 | endPhi = phiStart+twopi; |
---|
253 | } |
---|
254 | else |
---|
255 | { |
---|
256 | phiIsOpen = true; |
---|
257 | |
---|
258 | // |
---|
259 | // Convert phi into our convention |
---|
260 | // |
---|
261 | endPhi = phiStart+phiTotal; |
---|
262 | while( endPhi < startPhi ) endPhi += twopi; |
---|
263 | } |
---|
264 | |
---|
265 | // |
---|
266 | // Save number sides |
---|
267 | // |
---|
268 | numSide = theNumSide; |
---|
269 | |
---|
270 | // |
---|
271 | // Allocate corner array. |
---|
272 | // |
---|
273 | corners = new G4PolyhedraSideRZ[numCorner]; |
---|
274 | |
---|
275 | // |
---|
276 | // Copy corners |
---|
277 | // |
---|
278 | G4ReduciblePolygonIterator iterRZ(rz); |
---|
279 | |
---|
280 | G4PolyhedraSideRZ *next = corners; |
---|
281 | iterRZ.Begin(); |
---|
282 | do |
---|
283 | { |
---|
284 | next->r = iterRZ.GetA(); |
---|
285 | next->z = iterRZ.GetB(); |
---|
286 | } while( ++next, iterRZ.Next() ); |
---|
287 | |
---|
288 | // |
---|
289 | // Allocate face pointer array |
---|
290 | // |
---|
291 | numFace = phiIsOpen ? numCorner+2 : numCorner; |
---|
292 | faces = new G4VCSGface*[numFace]; |
---|
293 | |
---|
294 | // |
---|
295 | // Construct side faces |
---|
296 | // |
---|
297 | // To do so properly, we need to keep track of four successive RZ |
---|
298 | // corners. |
---|
299 | // |
---|
300 | // But! Don't construct a face if both points are at zero radius! |
---|
301 | // |
---|
302 | G4PolyhedraSideRZ *corner = corners, |
---|
303 | *prev = corners + numCorner-1, |
---|
304 | *nextNext; |
---|
305 | G4VCSGface **face = faces; |
---|
306 | do |
---|
307 | { |
---|
308 | next = corner+1; |
---|
309 | if (next >= corners+numCorner) next = corners; |
---|
310 | nextNext = next+1; |
---|
311 | if (nextNext >= corners+numCorner) nextNext = corners; |
---|
312 | |
---|
313 | if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue; |
---|
314 | |
---|
315 | // |
---|
316 | // We must decide here if we can dare declare one of our faces |
---|
317 | // as having a "valid" normal (i.e. allBehind = true). This |
---|
318 | // is never possible if the face faces "inward" in r *unless* |
---|
319 | // we have only one side |
---|
320 | // |
---|
321 | G4bool allBehind; |
---|
322 | if ((corner->z > next->z) && (numSide > 1)) |
---|
323 | { |
---|
324 | allBehind = false; |
---|
325 | } |
---|
326 | else |
---|
327 | { |
---|
328 | // |
---|
329 | // Otherwise, it is only true if the line passing |
---|
330 | // through the two points of the segment do not |
---|
331 | // split the r/z cross section |
---|
332 | // |
---|
333 | allBehind = !rz->BisectedBy( corner->r, corner->z, |
---|
334 | next->r, next->z, kCarTolerance ); |
---|
335 | } |
---|
336 | |
---|
337 | *face++ = new G4PolyhedraSide( prev, corner, next, nextNext, |
---|
338 | numSide, startPhi, endPhi-startPhi, phiIsOpen ); |
---|
339 | } while( prev=corner, corner=next, corner > corners ); |
---|
340 | |
---|
341 | if (phiIsOpen) |
---|
342 | { |
---|
343 | // |
---|
344 | // Construct phi open edges |
---|
345 | // |
---|
346 | *face++ = new G4PolyPhiFace( rz, startPhi, phiTotal/numSide, endPhi ); |
---|
347 | *face++ = new G4PolyPhiFace( rz, endPhi, phiTotal/numSide, startPhi ); |
---|
348 | } |
---|
349 | |
---|
350 | // |
---|
351 | // We might have dropped a face or two: recalculate numFace |
---|
352 | // |
---|
353 | numFace = face-faces; |
---|
354 | |
---|
355 | // |
---|
356 | // Make enclosingCylinder |
---|
357 | // |
---|
358 | enclosingCylinder = |
---|
359 | new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal ); |
---|
360 | } |
---|
361 | |
---|
362 | |
---|
363 | // |
---|
364 | // Fake default constructor - sets only member data and allocates memory |
---|
365 | // for usage restricted to object persistency. |
---|
366 | // |
---|
367 | G4Polyhedra::G4Polyhedra( __void__& a ) |
---|
368 | : G4VCSGfaceted(a), genericPgon(false), corners(0), |
---|
369 | original_parameters(0), enclosingCylinder(0) |
---|
370 | { |
---|
371 | } |
---|
372 | |
---|
373 | |
---|
374 | // |
---|
375 | // Destructor |
---|
376 | // |
---|
377 | G4Polyhedra::~G4Polyhedra() |
---|
378 | { |
---|
379 | delete [] corners; |
---|
380 | if (original_parameters) delete original_parameters; |
---|
381 | |
---|
382 | delete enclosingCylinder; |
---|
383 | } |
---|
384 | |
---|
385 | |
---|
386 | // |
---|
387 | // Copy constructor |
---|
388 | // |
---|
389 | G4Polyhedra::G4Polyhedra( const G4Polyhedra &source ) |
---|
390 | : G4VCSGfaceted( source ) |
---|
391 | { |
---|
392 | CopyStuff( source ); |
---|
393 | } |
---|
394 | |
---|
395 | |
---|
396 | // |
---|
397 | // Assignment operator |
---|
398 | // |
---|
399 | const G4Polyhedra &G4Polyhedra::operator=( const G4Polyhedra &source ) |
---|
400 | { |
---|
401 | if (this == &source) return *this; |
---|
402 | |
---|
403 | G4VCSGfaceted::operator=( source ); |
---|
404 | |
---|
405 | delete [] corners; |
---|
406 | if (original_parameters) delete original_parameters; |
---|
407 | |
---|
408 | delete enclosingCylinder; |
---|
409 | |
---|
410 | CopyStuff( source ); |
---|
411 | |
---|
412 | return *this; |
---|
413 | } |
---|
414 | |
---|
415 | |
---|
416 | // |
---|
417 | // CopyStuff |
---|
418 | // |
---|
419 | void G4Polyhedra::CopyStuff( const G4Polyhedra &source ) |
---|
420 | { |
---|
421 | // |
---|
422 | // Simple stuff |
---|
423 | // |
---|
424 | numSide = source.numSide; |
---|
425 | startPhi = source.startPhi; |
---|
426 | endPhi = source.endPhi; |
---|
427 | phiIsOpen = source.phiIsOpen; |
---|
428 | numCorner = source.numCorner; |
---|
429 | genericPgon= source.genericPgon; |
---|
430 | |
---|
431 | // |
---|
432 | // The corner array |
---|
433 | // |
---|
434 | corners = new G4PolyhedraSideRZ[numCorner]; |
---|
435 | |
---|
436 | G4PolyhedraSideRZ *corn = corners, |
---|
437 | *sourceCorn = source.corners; |
---|
438 | do |
---|
439 | { |
---|
440 | *corn = *sourceCorn; |
---|
441 | } while( ++sourceCorn, ++corn < corners+numCorner ); |
---|
442 | |
---|
443 | // |
---|
444 | // Original parameters |
---|
445 | // |
---|
446 | if (source.original_parameters) |
---|
447 | { |
---|
448 | original_parameters = |
---|
449 | new G4PolyhedraHistorical( *source.original_parameters ); |
---|
450 | } |
---|
451 | |
---|
452 | // |
---|
453 | // Enclosing cylinder |
---|
454 | // |
---|
455 | enclosingCylinder = new G4EnclosingCylinder( *source.enclosingCylinder ); |
---|
456 | } |
---|
457 | |
---|
458 | |
---|
459 | // |
---|
460 | // Reset |
---|
461 | // |
---|
462 | // Recalculates and reshapes the solid, given pre-assigned scaled |
---|
463 | // original_parameters. |
---|
464 | // |
---|
465 | G4bool G4Polyhedra::Reset() |
---|
466 | { |
---|
467 | if (genericPgon) |
---|
468 | { |
---|
469 | G4cerr << "Solid " << GetName() << " built using generic construct." |
---|
470 | << G4endl << "Not applicable to the generic construct !" << G4endl; |
---|
471 | G4Exception("G4Polyhedra::Reset()", "NotApplicableConstruct", |
---|
472 | JustWarning, "Parameters NOT resetted."); |
---|
473 | return 1; |
---|
474 | } |
---|
475 | |
---|
476 | // |
---|
477 | // Clear old setup |
---|
478 | // |
---|
479 | G4VCSGfaceted::DeleteStuff(); |
---|
480 | delete [] corners; |
---|
481 | delete enclosingCylinder; |
---|
482 | |
---|
483 | // |
---|
484 | // Rebuild polyhedra |
---|
485 | // |
---|
486 | G4ReduciblePolygon *rz = |
---|
487 | new G4ReduciblePolygon( original_parameters->Rmin, |
---|
488 | original_parameters->Rmax, |
---|
489 | original_parameters->Z_values, |
---|
490 | original_parameters->Num_z_planes ); |
---|
491 | Create( original_parameters->Start_angle, |
---|
492 | original_parameters->Opening_angle, |
---|
493 | original_parameters->numSide, rz ); |
---|
494 | delete rz; |
---|
495 | |
---|
496 | return 0; |
---|
497 | } |
---|
498 | |
---|
499 | |
---|
500 | // |
---|
501 | // Inside |
---|
502 | // |
---|
503 | // This is an override of G4VCSGfaceted::Inside, created in order |
---|
504 | // to speed things up by first checking with G4EnclosingCylinder. |
---|
505 | // |
---|
506 | EInside G4Polyhedra::Inside( const G4ThreeVector &p ) const |
---|
507 | { |
---|
508 | // |
---|
509 | // Quick test |
---|
510 | // |
---|
511 | if (enclosingCylinder->MustBeOutside(p)) return kOutside; |
---|
512 | |
---|
513 | // |
---|
514 | // Long answer |
---|
515 | // |
---|
516 | return G4VCSGfaceted::Inside(p); |
---|
517 | } |
---|
518 | |
---|
519 | |
---|
520 | // |
---|
521 | // DistanceToIn |
---|
522 | // |
---|
523 | // This is an override of G4VCSGfaceted::Inside, created in order |
---|
524 | // to speed things up by first checking with G4EnclosingCylinder. |
---|
525 | // |
---|
526 | G4double G4Polyhedra::DistanceToIn( const G4ThreeVector &p, |
---|
527 | const G4ThreeVector &v ) const |
---|
528 | { |
---|
529 | // |
---|
530 | // Quick test |
---|
531 | // |
---|
532 | if (enclosingCylinder->ShouldMiss(p,v)) |
---|
533 | return kInfinity; |
---|
534 | |
---|
535 | // |
---|
536 | // Long answer |
---|
537 | // |
---|
538 | return G4VCSGfaceted::DistanceToIn( p, v ); |
---|
539 | } |
---|
540 | |
---|
541 | |
---|
542 | // |
---|
543 | // DistanceToIn |
---|
544 | // |
---|
545 | G4double G4Polyhedra::DistanceToIn( const G4ThreeVector &p ) const |
---|
546 | { |
---|
547 | return G4VCSGfaceted::DistanceToIn(p); |
---|
548 | } |
---|
549 | |
---|
550 | |
---|
551 | // |
---|
552 | // ComputeDimensions |
---|
553 | // |
---|
554 | void G4Polyhedra::ComputeDimensions( G4VPVParameterisation* p, |
---|
555 | const G4int n, |
---|
556 | const G4VPhysicalVolume* pRep ) |
---|
557 | { |
---|
558 | p->ComputeDimensions(*this,n,pRep); |
---|
559 | } |
---|
560 | |
---|
561 | |
---|
562 | // |
---|
563 | // GetEntityType |
---|
564 | // |
---|
565 | G4GeometryType G4Polyhedra::GetEntityType() const |
---|
566 | { |
---|
567 | return G4String("G4Polyhedra"); |
---|
568 | } |
---|
569 | |
---|
570 | |
---|
571 | // |
---|
572 | // Stream object contents to an output stream |
---|
573 | // |
---|
574 | std::ostream& G4Polyhedra::StreamInfo( std::ostream& os ) const |
---|
575 | { |
---|
576 | os << "-----------------------------------------------------------\n" |
---|
577 | << " *** Dump for solid - " << GetName() << " ***\n" |
---|
578 | << " ===================================================\n" |
---|
579 | << " Solid type: G4Polyhedra\n" |
---|
580 | << " Parameters: \n" |
---|
581 | << " starting phi angle : " << startPhi/degree << " degrees \n" |
---|
582 | << " ending phi angle : " << endPhi/degree << " degrees \n"; |
---|
583 | G4int i=0; |
---|
584 | if (!genericPgon) |
---|
585 | { |
---|
586 | G4int numPlanes = original_parameters->Num_z_planes; |
---|
587 | os << " number of Z planes: " << numPlanes << "\n" |
---|
588 | << " Z values: \n"; |
---|
589 | for (i=0; i<numPlanes; i++) |
---|
590 | { |
---|
591 | os << " Z plane " << i << ": " |
---|
592 | << original_parameters->Z_values[i] << "\n"; |
---|
593 | } |
---|
594 | os << " Tangent distances to inner surface (Rmin): \n"; |
---|
595 | for (i=0; i<numPlanes; i++) |
---|
596 | { |
---|
597 | os << " Z plane " << i << ": " |
---|
598 | << original_parameters->Rmin[i] << "\n"; |
---|
599 | } |
---|
600 | os << " Tangent distances to outer surface (Rmax): \n"; |
---|
601 | for (i=0; i<numPlanes; i++) |
---|
602 | { |
---|
603 | os << " Z plane " << i << ": " |
---|
604 | << original_parameters->Rmax[i] << "\n"; |
---|
605 | } |
---|
606 | } |
---|
607 | os << " number of RZ points: " << numCorner << "\n" |
---|
608 | << " RZ values (corners): \n"; |
---|
609 | for (i=0; i<numCorner; i++) |
---|
610 | { |
---|
611 | os << " " |
---|
612 | << corners[i].r << ", " << corners[i].z << "\n"; |
---|
613 | } |
---|
614 | os << "-----------------------------------------------------------\n"; |
---|
615 | |
---|
616 | return os; |
---|
617 | } |
---|
618 | |
---|
619 | |
---|
620 | // |
---|
621 | // GetPointOnPlane |
---|
622 | // |
---|
623 | // Auxiliary method for get point on surface |
---|
624 | // |
---|
625 | G4ThreeVector G4Polyhedra::GetPointOnPlane(G4ThreeVector p0, G4ThreeVector p1, |
---|
626 | G4ThreeVector p2, G4ThreeVector p3) const |
---|
627 | { |
---|
628 | G4double lambda1, lambda2, chose,aOne,aTwo; |
---|
629 | G4ThreeVector t, u, v, w, Area, normal; |
---|
630 | aOne = 1.; |
---|
631 | aTwo = 1.; |
---|
632 | |
---|
633 | t = p1 - p0; |
---|
634 | u = p2 - p1; |
---|
635 | v = p3 - p2; |
---|
636 | w = p0 - p3; |
---|
637 | |
---|
638 | chose = RandFlat::shoot(0.,aOne+aTwo); |
---|
639 | if( (chose>=0.) && (chose < aOne) ) |
---|
640 | { |
---|
641 | lambda1 = RandFlat::shoot(0.,1.); |
---|
642 | lambda2 = RandFlat::shoot(0.,lambda1); |
---|
643 | return (p2+lambda1*v+lambda2*w); |
---|
644 | } |
---|
645 | |
---|
646 | lambda1 = RandFlat::shoot(0.,1.); |
---|
647 | lambda2 = RandFlat::shoot(0.,lambda1); |
---|
648 | return (p0+lambda1*t+lambda2*u); |
---|
649 | } |
---|
650 | |
---|
651 | |
---|
652 | // |
---|
653 | // GetPointOnTriangle |
---|
654 | // |
---|
655 | // Auxiliary method for get point on surface |
---|
656 | // |
---|
657 | G4ThreeVector G4Polyhedra::GetPointOnTriangle(G4ThreeVector p1, |
---|
658 | G4ThreeVector p2, |
---|
659 | G4ThreeVector p3) const |
---|
660 | { |
---|
661 | G4double lambda1,lambda2; |
---|
662 | G4ThreeVector v=p3-p1, w=p1-p2; |
---|
663 | |
---|
664 | lambda1 = RandFlat::shoot(0.,1.); |
---|
665 | lambda2 = RandFlat::shoot(0.,lambda1); |
---|
666 | |
---|
667 | return (p2 + lambda1*w + lambda2*v); |
---|
668 | } |
---|
669 | |
---|
670 | |
---|
671 | // |
---|
672 | // GetPointOnSurface |
---|
673 | // |
---|
674 | G4ThreeVector G4Polyhedra::GetPointOnSurface() const |
---|
675 | { |
---|
676 | if( !genericPgon ) // Polyhedra by faces |
---|
677 | { |
---|
678 | G4int j, numPlanes = original_parameters->Num_z_planes, Flag=0; |
---|
679 | G4double chose, totArea=0., Achose1, Achose2, |
---|
680 | rad1, rad2, sinphi1, sinphi2, cosphi1, cosphi2; |
---|
681 | G4double a, b, l2, rang, totalPhi, ksi, |
---|
682 | area, aTop=0., aBottom=0., zVal=0.; |
---|
683 | |
---|
684 | G4ThreeVector p0, p1, p2, p3; |
---|
685 | std::vector<G4double> aVector1; |
---|
686 | std::vector<G4double> aVector2; |
---|
687 | std::vector<G4double> aVector3; |
---|
688 | |
---|
689 | totalPhi= (phiIsOpen) ? (endPhi-startPhi) : twopi; |
---|
690 | ksi = totalPhi/numSide; |
---|
691 | G4double cosksi = std::cos(ksi/2.); |
---|
692 | |
---|
693 | // Below we generate the areas relevant to our solid |
---|
694 | // |
---|
695 | for(j=0; j<numPlanes-1; j++) |
---|
696 | { |
---|
697 | a = original_parameters->Rmax[j+1]; |
---|
698 | b = original_parameters->Rmax[j]; |
---|
699 | l2 = sqr(original_parameters->Z_values[j] |
---|
700 | -original_parameters->Z_values[j+1]) + sqr(b-a); |
---|
701 | area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi; |
---|
702 | aVector1.push_back(area); |
---|
703 | } |
---|
704 | |
---|
705 | for(j=0; j<numPlanes-1; j++) |
---|
706 | { |
---|
707 | a = original_parameters->Rmin[j+1];//*cosksi; |
---|
708 | b = original_parameters->Rmin[j];//*cosksi; |
---|
709 | l2 = sqr(original_parameters->Z_values[j] |
---|
710 | -original_parameters->Z_values[j+1]) + sqr(b-a); |
---|
711 | area = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi; |
---|
712 | aVector2.push_back(area); |
---|
713 | } |
---|
714 | |
---|
715 | for(j=0; j<numPlanes-1; j++) |
---|
716 | { |
---|
717 | if(phiIsOpen == true) |
---|
718 | { |
---|
719 | aVector3.push_back(0.5*(original_parameters->Rmax[j] |
---|
720 | -original_parameters->Rmin[j] |
---|
721 | +original_parameters->Rmax[j+1] |
---|
722 | -original_parameters->Rmin[j+1]) |
---|
723 | *std::fabs(original_parameters->Z_values[j+1] |
---|
724 | -original_parameters->Z_values[j])); |
---|
725 | } |
---|
726 | else { aVector3.push_back(0.); } |
---|
727 | } |
---|
728 | |
---|
729 | for(j=0; j<numPlanes-1; j++) |
---|
730 | { |
---|
731 | totArea += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j]; |
---|
732 | } |
---|
733 | |
---|
734 | // Must include top and bottom areas |
---|
735 | // |
---|
736 | if(original_parameters->Rmax[numPlanes-1] != 0.) |
---|
737 | { |
---|
738 | a = original_parameters->Rmax[numPlanes-1]; |
---|
739 | b = original_parameters->Rmin[numPlanes-1]; |
---|
740 | l2 = sqr(a-b); |
---|
741 | aTop = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi; |
---|
742 | } |
---|
743 | |
---|
744 | if(original_parameters->Rmax[0] != 0.) |
---|
745 | { |
---|
746 | a = original_parameters->Rmax[0]; |
---|
747 | b = original_parameters->Rmin[0]; |
---|
748 | l2 = sqr(a-b); |
---|
749 | aBottom = std::sqrt(l2-sqr((a-b)*cosksi))*(a+b)*cosksi; |
---|
750 | } |
---|
751 | |
---|
752 | Achose1 = 0.; |
---|
753 | Achose2 = numSide*(aVector1[0]+aVector2[0])+2.*aVector3[0]; |
---|
754 | |
---|
755 | chose = RandFlat::shoot(0.,totArea+aTop+aBottom); |
---|
756 | if( (chose >= 0.) && (chose < aTop + aBottom) ) |
---|
757 | { |
---|
758 | chose = RandFlat::shoot(startPhi,startPhi+totalPhi); |
---|
759 | rang = std::floor((chose-startPhi)/ksi-0.01); |
---|
760 | if(rang<0) { rang=0; } |
---|
761 | rang = std::fabs(rang); |
---|
762 | sinphi1 = std::sin(startPhi+rang*ksi); |
---|
763 | sinphi2 = std::sin(startPhi+(rang+1)*ksi); |
---|
764 | cosphi1 = std::cos(startPhi+rang*ksi); |
---|
765 | cosphi2 = std::cos(startPhi+(rang+1)*ksi); |
---|
766 | chose = RandFlat::shoot(0., aTop + aBottom); |
---|
767 | if(chose>=0. && chose<aTop) |
---|
768 | { |
---|
769 | rad1 = original_parameters->Rmin[numPlanes-1]; |
---|
770 | rad2 = original_parameters->Rmax[numPlanes-1]; |
---|
771 | zVal = original_parameters->Z_values[numPlanes-1]; |
---|
772 | } |
---|
773 | else |
---|
774 | { |
---|
775 | rad1 = original_parameters->Rmin[0]; |
---|
776 | rad2 = original_parameters->Rmax[0]; |
---|
777 | zVal = original_parameters->Z_values[0]; |
---|
778 | } |
---|
779 | p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal); |
---|
780 | p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal); |
---|
781 | p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal); |
---|
782 | p3 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal); |
---|
783 | return GetPointOnPlane(p0,p1,p2,p3); |
---|
784 | } |
---|
785 | else |
---|
786 | { |
---|
787 | for (j=0; j<numPlanes-1; j++) |
---|
788 | { |
---|
789 | if( ((chose >= Achose1) && (chose < Achose2)) || (j == numPlanes-2) ) |
---|
790 | { |
---|
791 | Flag = j; break; |
---|
792 | } |
---|
793 | Achose1 += numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j]; |
---|
794 | Achose2 = Achose1 + numSide*(aVector1[j+1]+aVector2[j+1]) |
---|
795 | + 2.*aVector3[j+1]; |
---|
796 | } |
---|
797 | } |
---|
798 | |
---|
799 | // At this point we have chosen a subsection |
---|
800 | // between to adjacent plane cuts... |
---|
801 | |
---|
802 | j = Flag; |
---|
803 | |
---|
804 | totArea = numSide*(aVector1[j]+aVector2[j])+2.*aVector3[j]; |
---|
805 | chose = RandFlat::shoot(0.,totArea); |
---|
806 | |
---|
807 | if( (chose>=0.) && (chose<numSide*aVector1[j]) ) |
---|
808 | { |
---|
809 | chose = RandFlat::shoot(startPhi,startPhi+totalPhi); |
---|
810 | rang = std::floor((chose-startPhi)/ksi-0.01); |
---|
811 | if(rang<0) { rang=0; } |
---|
812 | rang = std::fabs(rang); |
---|
813 | rad1 = original_parameters->Rmax[j]; |
---|
814 | rad2 = original_parameters->Rmax[j+1]; |
---|
815 | sinphi1 = std::sin(startPhi+rang*ksi); |
---|
816 | sinphi2 = std::sin(startPhi+(rang+1)*ksi); |
---|
817 | cosphi1 = std::cos(startPhi+rang*ksi); |
---|
818 | cosphi2 = std::cos(startPhi+(rang+1)*ksi); |
---|
819 | zVal = original_parameters->Z_values[j]; |
---|
820 | |
---|
821 | p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal); |
---|
822 | p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal); |
---|
823 | |
---|
824 | zVal = original_parameters->Z_values[j+1]; |
---|
825 | |
---|
826 | p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal); |
---|
827 | p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal); |
---|
828 | return GetPointOnPlane(p0,p1,p2,p3); |
---|
829 | } |
---|
830 | else if ( (chose >= numSide*aVector1[j]) |
---|
831 | && (chose <= numSide*(aVector1[j]+aVector2[j])) ) |
---|
832 | { |
---|
833 | chose = RandFlat::shoot(startPhi,startPhi+totalPhi); |
---|
834 | rang = std::floor((chose-startPhi)/ksi-0.01); |
---|
835 | if(rang<0) { rang=0; } |
---|
836 | rang = std::fabs(rang); |
---|
837 | rad1 = original_parameters->Rmin[j]; |
---|
838 | rad2 = original_parameters->Rmin[j+1]; |
---|
839 | sinphi1 = std::sin(startPhi+rang*ksi); |
---|
840 | sinphi2 = std::sin(startPhi+(rang+1)*ksi); |
---|
841 | cosphi1 = std::cos(startPhi+rang*ksi); |
---|
842 | cosphi2 = std::cos(startPhi+(rang+1)*ksi); |
---|
843 | zVal = original_parameters->Z_values[j]; |
---|
844 | |
---|
845 | p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1,zVal); |
---|
846 | p1 = G4ThreeVector(rad1*cosphi2,rad1*sinphi2,zVal); |
---|
847 | |
---|
848 | zVal = original_parameters->Z_values[j+1]; |
---|
849 | |
---|
850 | p2 = G4ThreeVector(rad2*cosphi2,rad2*sinphi2,zVal); |
---|
851 | p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1,zVal); |
---|
852 | return GetPointOnPlane(p0,p1,p2,p3); |
---|
853 | } |
---|
854 | |
---|
855 | chose = RandFlat::shoot(0.,2.2); |
---|
856 | if( (chose>=0.) && (chose < 1.) ) |
---|
857 | { |
---|
858 | rang = startPhi; |
---|
859 | } |
---|
860 | else |
---|
861 | { |
---|
862 | rang = endPhi; |
---|
863 | } |
---|
864 | |
---|
865 | cosphi1 = std::cos(rang); rad1 = original_parameters->Rmin[j]; |
---|
866 | sinphi1 = std::sin(rang); rad2 = original_parameters->Rmax[j]; |
---|
867 | |
---|
868 | p0 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1, |
---|
869 | original_parameters->Z_values[j]); |
---|
870 | p1 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1, |
---|
871 | original_parameters->Z_values[j]); |
---|
872 | |
---|
873 | rad1 = original_parameters->Rmax[j+1]; |
---|
874 | rad2 = original_parameters->Rmin[j+1]; |
---|
875 | |
---|
876 | p2 = G4ThreeVector(rad1*cosphi1,rad1*sinphi1, |
---|
877 | original_parameters->Z_values[j+1]); |
---|
878 | p3 = G4ThreeVector(rad2*cosphi1,rad2*sinphi1, |
---|
879 | original_parameters->Z_values[j+1]); |
---|
880 | return GetPointOnPlane(p0,p1,p2,p3); |
---|
881 | } |
---|
882 | else // Generic polyhedra |
---|
883 | { |
---|
884 | return GetPointOnSurfaceGeneric(); |
---|
885 | } |
---|
886 | } |
---|
887 | |
---|
888 | // |
---|
889 | // CreatePolyhedron |
---|
890 | // |
---|
891 | G4Polyhedron* G4Polyhedra::CreatePolyhedron() const |
---|
892 | { |
---|
893 | if (!genericPgon) |
---|
894 | { |
---|
895 | return new G4PolyhedronPgon( original_parameters->Start_angle, |
---|
896 | original_parameters->Opening_angle, |
---|
897 | original_parameters->numSide, |
---|
898 | original_parameters->Num_z_planes, |
---|
899 | original_parameters->Z_values, |
---|
900 | original_parameters->Rmin, |
---|
901 | original_parameters->Rmax); |
---|
902 | } |
---|
903 | else |
---|
904 | { |
---|
905 | // The following code prepares for: |
---|
906 | // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces, |
---|
907 | // const double xyz[][3], |
---|
908 | // const int faces_vec[][4]) |
---|
909 | // Here is an extract from the header file HepPolyhedron.h: |
---|
910 | /** |
---|
911 | * Creates user defined polyhedron. |
---|
912 | * This function allows to the user to define arbitrary polyhedron. |
---|
913 | * The faces of the polyhedron should be either triangles or planar |
---|
914 | * quadrilateral. Nodes of a face are defined by indexes pointing to |
---|
915 | * the elements in the xyz array. Numeration of the elements in the |
---|
916 | * array starts from 1 (like in fortran). The indexes can be positive |
---|
917 | * or negative. Negative sign means that the corresponding edge is |
---|
918 | * invisible. The normal of the face should be directed to exterior |
---|
919 | * of the polyhedron. |
---|
920 | * |
---|
921 | * @param Nnodes number of nodes |
---|
922 | * @param Nfaces number of faces |
---|
923 | * @param xyz nodes |
---|
924 | * @param faces_vec faces (quadrilaterals or triangles) |
---|
925 | * @return status of the operation - is non-zero in case of problem |
---|
926 | */ |
---|
927 | G4int nNodes; |
---|
928 | G4int nFaces; |
---|
929 | typedef G4double double3[3]; |
---|
930 | double3* xyz; |
---|
931 | typedef G4int int4[4]; |
---|
932 | int4* faces_vec; |
---|
933 | if (phiIsOpen) |
---|
934 | { |
---|
935 | // Triangulate open ends. Simple ear-chopping algorithm... |
---|
936 | // I'm not sure how robust this algorithm is (J.Allison). |
---|
937 | // |
---|
938 | std::vector<G4bool> chopped(numCorner, false); |
---|
939 | std::vector<G4int*> triQuads; |
---|
940 | G4int remaining = numCorner; |
---|
941 | G4int iStarter = 0; |
---|
942 | while (remaining >= 3) |
---|
943 | { |
---|
944 | // Find unchopped corners... |
---|
945 | // |
---|
946 | G4int A = -1, B = -1, C = -1; |
---|
947 | G4int iStepper = iStarter; |
---|
948 | do |
---|
949 | { |
---|
950 | if (A < 0) { A = iStepper; } |
---|
951 | else if (B < 0) { B = iStepper; } |
---|
952 | else if (C < 0) { C = iStepper; } |
---|
953 | do |
---|
954 | { |
---|
955 | if (++iStepper >= numCorner) iStepper = 0; |
---|
956 | } |
---|
957 | while (chopped[iStepper]); |
---|
958 | } |
---|
959 | while (C < 0 && iStepper != iStarter); |
---|
960 | |
---|
961 | // Check triangle at B is pointing outward (an "ear"). |
---|
962 | // Sign of z cross product determines... |
---|
963 | |
---|
964 | G4double BAr = corners[A].r - corners[B].r; |
---|
965 | G4double BAz = corners[A].z - corners[B].z; |
---|
966 | G4double BCr = corners[C].r - corners[B].r; |
---|
967 | G4double BCz = corners[C].z - corners[B].z; |
---|
968 | if (BAr * BCz - BAz * BCr < kCarTolerance) |
---|
969 | { |
---|
970 | G4int* tq = new G4int[3]; |
---|
971 | tq[0] = A + 1; |
---|
972 | tq[1] = B + 1; |
---|
973 | tq[2] = C + 1; |
---|
974 | triQuads.push_back(tq); |
---|
975 | chopped[B] = true; |
---|
976 | --remaining; |
---|
977 | } |
---|
978 | else |
---|
979 | { |
---|
980 | do |
---|
981 | { |
---|
982 | if (++iStarter >= numCorner) { iStarter = 0; } |
---|
983 | } |
---|
984 | while (chopped[iStarter]); |
---|
985 | } |
---|
986 | } |
---|
987 | |
---|
988 | // Transfer to faces... |
---|
989 | |
---|
990 | nNodes = (numSide + 1) * numCorner; |
---|
991 | nFaces = numSide * numCorner + 2 * triQuads.size(); |
---|
992 | faces_vec = new int4[nFaces]; |
---|
993 | G4int iface = 0; |
---|
994 | G4int addition = numCorner * numSide; |
---|
995 | G4int d = numCorner - 1; |
---|
996 | for (G4int iEnd = 0; iEnd < 2; ++iEnd) |
---|
997 | { |
---|
998 | for (size_t i = 0; i < triQuads.size(); ++i) |
---|
999 | { |
---|
1000 | // Negative for soft/auxiliary/normally invisible edges... |
---|
1001 | // |
---|
1002 | G4int a, b, c; |
---|
1003 | if (iEnd == 0) |
---|
1004 | { |
---|
1005 | a = triQuads[i][0]; |
---|
1006 | b = triQuads[i][1]; |
---|
1007 | c = triQuads[i][2]; |
---|
1008 | } |
---|
1009 | else |
---|
1010 | { |
---|
1011 | a = triQuads[i][0] + addition; |
---|
1012 | b = triQuads[i][2] + addition; |
---|
1013 | c = triQuads[i][1] + addition; |
---|
1014 | } |
---|
1015 | G4int ab = std::abs(b - a); |
---|
1016 | G4int bc = std::abs(c - b); |
---|
1017 | G4int ca = std::abs(a - c); |
---|
1018 | faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a; |
---|
1019 | faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b; |
---|
1020 | faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c; |
---|
1021 | faces_vec[iface][3] = 0; |
---|
1022 | ++iface; |
---|
1023 | } |
---|
1024 | } |
---|
1025 | |
---|
1026 | // Continue with sides... |
---|
1027 | |
---|
1028 | xyz = new double3[nNodes]; |
---|
1029 | const G4double dPhi = (endPhi - startPhi) / numSide; |
---|
1030 | G4double phi = startPhi; |
---|
1031 | G4int ixyz = 0; |
---|
1032 | for (G4int iSide = 0; iSide < numSide; ++iSide) |
---|
1033 | { |
---|
1034 | for (G4int iCorner = 0; iCorner < numCorner; ++iCorner) |
---|
1035 | { |
---|
1036 | xyz[ixyz][0] = corners[iCorner].r * std::cos(phi); |
---|
1037 | xyz[ixyz][1] = corners[iCorner].r * std::sin(phi); |
---|
1038 | xyz[ixyz][2] = corners[iCorner].z; |
---|
1039 | if (iCorner < numCorner - 1) |
---|
1040 | { |
---|
1041 | faces_vec[iface][0] = ixyz + 1; |
---|
1042 | faces_vec[iface][1] = ixyz + numCorner + 1; |
---|
1043 | faces_vec[iface][2] = ixyz + numCorner + 2; |
---|
1044 | faces_vec[iface][3] = ixyz + 2; |
---|
1045 | } |
---|
1046 | else |
---|
1047 | { |
---|
1048 | faces_vec[iface][0] = ixyz + 1; |
---|
1049 | faces_vec[iface][1] = ixyz + numCorner + 1; |
---|
1050 | faces_vec[iface][2] = ixyz + 2; |
---|
1051 | faces_vec[iface][3] = ixyz - numCorner + 2; |
---|
1052 | } |
---|
1053 | ++iface; |
---|
1054 | ++ixyz; |
---|
1055 | } |
---|
1056 | phi += dPhi; |
---|
1057 | } |
---|
1058 | |
---|
1059 | // Last corners... |
---|
1060 | |
---|
1061 | for (G4int iCorner = 0; iCorner < numCorner; ++iCorner) |
---|
1062 | { |
---|
1063 | xyz[ixyz][0] = corners[iCorner].r * std::cos(phi); |
---|
1064 | xyz[ixyz][1] = corners[iCorner].r * std::sin(phi); |
---|
1065 | xyz[ixyz][2] = corners[iCorner].z; |
---|
1066 | ++ixyz; |
---|
1067 | } |
---|
1068 | } |
---|
1069 | else // !phiIsOpen - i.e., a complete 360 degrees. |
---|
1070 | { |
---|
1071 | nNodes = numSide * numCorner; |
---|
1072 | nFaces = numSide * numCorner;; |
---|
1073 | xyz = new double3[nNodes]; |
---|
1074 | faces_vec = new int4[nFaces]; |
---|
1075 | // const G4double dPhi = (endPhi - startPhi) / numSide; |
---|
1076 | const G4double dPhi = twopi / numSide; // !phiIsOpen endPhi-startPhi = 360 degrees. |
---|
1077 | G4double phi = startPhi; |
---|
1078 | G4int ixyz = 0, iface = 0; |
---|
1079 | for (G4int iSide = 0; iSide < numSide; ++iSide) |
---|
1080 | { |
---|
1081 | for (G4int iCorner = 0; iCorner < numCorner; ++iCorner) |
---|
1082 | { |
---|
1083 | xyz[ixyz][0] = corners[iCorner].r * std::cos(phi); |
---|
1084 | xyz[ixyz][1] = corners[iCorner].r * std::sin(phi); |
---|
1085 | xyz[ixyz][2] = corners[iCorner].z; |
---|
1086 | if (iSide < numSide - 1) |
---|
1087 | { |
---|
1088 | if (iCorner < numCorner - 1) |
---|
1089 | { |
---|
1090 | faces_vec[iface][0] = ixyz + 1; |
---|
1091 | faces_vec[iface][1] = ixyz + numCorner + 1; |
---|
1092 | faces_vec[iface][2] = ixyz + numCorner + 2; |
---|
1093 | faces_vec[iface][3] = ixyz + 2; |
---|
1094 | } |
---|
1095 | else |
---|
1096 | { |
---|
1097 | faces_vec[iface][0] = ixyz + 1; |
---|
1098 | faces_vec[iface][1] = ixyz + numCorner + 1; |
---|
1099 | faces_vec[iface][2] = ixyz + 2; |
---|
1100 | faces_vec[iface][3] = ixyz - numCorner + 2; |
---|
1101 | } |
---|
1102 | } |
---|
1103 | else // Last side joins ends... |
---|
1104 | { |
---|
1105 | if (iCorner < numCorner - 1) |
---|
1106 | { |
---|
1107 | faces_vec[iface][0] = ixyz + 1; |
---|
1108 | faces_vec[iface][1] = ixyz + numCorner - nFaces + 1; |
---|
1109 | faces_vec[iface][2] = ixyz + numCorner - nFaces + 2; |
---|
1110 | faces_vec[iface][3] = ixyz + 2; |
---|
1111 | } |
---|
1112 | else |
---|
1113 | { |
---|
1114 | faces_vec[iface][0] = ixyz + 1; |
---|
1115 | faces_vec[iface][1] = ixyz - nFaces + numCorner + 1; |
---|
1116 | faces_vec[iface][2] = ixyz - nFaces + 2; |
---|
1117 | faces_vec[iface][3] = ixyz - numCorner + 2; |
---|
1118 | } |
---|
1119 | } |
---|
1120 | ++ixyz; |
---|
1121 | ++iface; |
---|
1122 | } |
---|
1123 | phi += dPhi; |
---|
1124 | } |
---|
1125 | } |
---|
1126 | G4Polyhedron* polyhedron = new G4Polyhedron; |
---|
1127 | G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec); |
---|
1128 | delete faces_vec; |
---|
1129 | delete xyz; |
---|
1130 | if (problem) |
---|
1131 | { |
---|
1132 | std::ostringstream oss; |
---|
1133 | oss << "Problem creating G4Polyhedron for: " << GetName(); |
---|
1134 | G4Exception("G4Polyhedra::CreatePolyhedron()", "BadPolyhedron", |
---|
1135 | JustWarning, oss.str().c_str()); |
---|
1136 | delete polyhedron; |
---|
1137 | return 0; |
---|
1138 | } |
---|
1139 | else |
---|
1140 | { |
---|
1141 | return polyhedron; |
---|
1142 | } |
---|
1143 | } |
---|
1144 | } |
---|
1145 | |
---|
1146 | // |
---|
1147 | // CreateNURBS |
---|
1148 | // |
---|
1149 | G4NURBS *G4Polyhedra::CreateNURBS() const |
---|
1150 | { |
---|
1151 | return 0; |
---|
1152 | } |
---|
1153 | |
---|
1154 | |
---|
1155 | // |
---|
1156 | // G4PolyhedraHistorical stuff |
---|
1157 | // |
---|
1158 | G4PolyhedraHistorical::G4PolyhedraHistorical() |
---|
1159 | : Z_values(0), Rmin(0), Rmax(0) |
---|
1160 | { |
---|
1161 | } |
---|
1162 | |
---|
1163 | G4PolyhedraHistorical::~G4PolyhedraHistorical() |
---|
1164 | { |
---|
1165 | delete [] Z_values; |
---|
1166 | delete [] Rmin; |
---|
1167 | delete [] Rmax; |
---|
1168 | } |
---|
1169 | |
---|
1170 | G4PolyhedraHistorical:: |
---|
1171 | G4PolyhedraHistorical( const G4PolyhedraHistorical& source ) |
---|
1172 | { |
---|
1173 | Start_angle = source.Start_angle; |
---|
1174 | Opening_angle = source.Opening_angle; |
---|
1175 | numSide = source.numSide; |
---|
1176 | Num_z_planes = source.Num_z_planes; |
---|
1177 | |
---|
1178 | Z_values = new G4double[Num_z_planes]; |
---|
1179 | Rmin = new G4double[Num_z_planes]; |
---|
1180 | Rmax = new G4double[Num_z_planes]; |
---|
1181 | |
---|
1182 | for( G4int i = 0; i < Num_z_planes; i++) |
---|
1183 | { |
---|
1184 | Z_values[i] = source.Z_values[i]; |
---|
1185 | Rmin[i] = source.Rmin[i]; |
---|
1186 | Rmax[i] = source.Rmax[i]; |
---|
1187 | } |
---|
1188 | } |
---|
1189 | |
---|
1190 | G4PolyhedraHistorical& |
---|
1191 | G4PolyhedraHistorical::operator=( const G4PolyhedraHistorical& right ) |
---|
1192 | { |
---|
1193 | if ( &right == this ) return *this; |
---|
1194 | |
---|
1195 | if (&right) |
---|
1196 | { |
---|
1197 | Start_angle = right.Start_angle; |
---|
1198 | Opening_angle = right.Opening_angle; |
---|
1199 | numSide = right.numSide; |
---|
1200 | Num_z_planes = right.Num_z_planes; |
---|
1201 | |
---|
1202 | delete [] Z_values; |
---|
1203 | delete [] Rmin; |
---|
1204 | delete [] Rmax; |
---|
1205 | Z_values = new G4double[Num_z_planes]; |
---|
1206 | Rmin = new G4double[Num_z_planes]; |
---|
1207 | Rmax = new G4double[Num_z_planes]; |
---|
1208 | |
---|
1209 | for( G4int i = 0; i < Num_z_planes; i++) |
---|
1210 | { |
---|
1211 | Z_values[i] = right.Z_values[i]; |
---|
1212 | Rmin[i] = right.Rmin[i]; |
---|
1213 | Rmax[i] = right.Rmax[i]; |
---|
1214 | } |
---|
1215 | } |
---|
1216 | return *this; |
---|
1217 | } |
---|