source: trunk/source/geometry/solids/specific/src/G4ExtrudedSolid.cc@ 1202

Last change on this file since 1202 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 23.0 KB
RevLine 
[831]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//
[921]27// $Id: G4ExtrudedSolid.cc,v 1.18 2008/10/30 11:47:45 ivana Exp $
[1058]28// GEANT4 tag $Name: geant4-09-02-ref-02 $
[831]29//
30//
31// --------------------------------------------------------------------
32// GEANT 4 class source file
33//
34// G4ExtrudedSolid.cc
35//
36// Author: Ivana Hrivnacova, IPN Orsay
37// --------------------------------------------------------------------
38
39#include <set>
40#include <algorithm>
41#include <cmath>
42
43#include "G4ExtrudedSolid.hh"
44#include "G4TriangularFacet.hh"
45#include "G4QuadrangularFacet.hh"
46
47//_____________________________________________________________________________
48
49G4ExtrudedSolid::G4ExtrudedSolid( const G4String& pName,
50 std::vector<G4TwoVector> polygon,
51 std::vector<ZSection> zsections)
52 : G4TessellatedSolid(pName),
53 fNv(polygon.size()),
54 fNz(zsections.size()),
55 fPolygon(),
56 fZSections(),
57 fTriangles(),
58 fIsConvex(false),
59 fGeometryType("G4ExtrudedSolid")
60
61{
62 // General constructor
63
[850]64 G4String errorDescription = "InvalidSetup in \"";
65 errorDescription += pName;
66 errorDescription += "\"";
67
[831]68 // First check input parameters
69
[850]70 if ( fNv < 3 )
71 {
72 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
73 FatalException, "Number of polygon vertices < 3");
[831]74 }
75
[850]76 if ( fNz < 2 )
77 {
78 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
79 FatalException, "Number of z-sides < 2");
[831]80 }
81
82 for ( G4int i=0; i<fNz-1; ++i )
83 {
84 if ( zsections[i].fZ > zsections[i+1].fZ )
85 {
[850]86 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
[831]87 FatalException,
88 "Z-sections have to be ordered by z value (z0 < z1 < z2 ...)");
89 }
[850]90 if ( std::fabs( zsections[i+1].fZ - zsections[i].fZ ) < kCarTolerance * 0.5 )
[831]91 {
[850]92 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
[831]93 FatalException,
94 "Z-sections with the same z position are not supported.");
95 }
96 }
97
98 // Copy polygon
99 //
100 for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
101
102 // Copy z-sections
103 //
104 for ( G4int i=0; i<fNz; ++i ) { fZSections.push_back(zsections[i]); }
105
106
107 G4bool result = MakeFacets();
108 if (!result)
109 {
[850]110 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
[831]111 FatalException, "Making facets failed.");
112 }
113 fIsConvex = IsConvex();
114
115
116 ComputeProjectionParameters();
117}
118
119//_____________________________________________________________________________
120
121G4ExtrudedSolid::G4ExtrudedSolid( const G4String& pName,
122 std::vector<G4TwoVector> polygon,
123 G4double dz,
124 G4TwoVector off1, G4double scale1,
125 G4TwoVector off2, G4double scale2 )
126 : G4TessellatedSolid(pName),
127 fNv(polygon.size()),
128 fNz(2),
129 fPolygon(),
130 fZSections(),
131 fTriangles(),
132 fIsConvex(false),
133 fGeometryType("G4ExtrudedSolid")
134
135{
136 // Special constructor for solid with 2 z-sections
137
[850]138 G4String errorDescription = "InvalidSetup in \"";
139 errorDescription += pName;
140 errorDescription += "\"";
141
[831]142 // First check input parameters
143 //
144 if ( fNv < 3 )
145 {
[850]146 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
[831]147 FatalException, "Number of polygon vertices < 3");
148 }
149
150 // Copy polygon
151 //
152 for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
153
154 // Copy z-sections
155 //
156 fZSections.push_back(ZSection(-dz, off1, scale1));
157 fZSections.push_back(ZSection( dz, off2, scale2));
158
159 G4bool result = MakeFacets();
160 if (!result)
161 {
[850]162 G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", errorDescription,
[831]163 FatalException, "Making facets failed.");
164 }
165 fIsConvex = IsConvex();
166
167 ComputeProjectionParameters();
168}
169
170//_____________________________________________________________________________
171
172G4ExtrudedSolid::G4ExtrudedSolid( __void__& a )
173 : G4TessellatedSolid(a)
174{
175 // Fake default constructor - sets only member data and allocates memory
176 // for usage restricted to object persistency.
177}
178
179
180//_____________________________________________________________________________
181
182G4ExtrudedSolid::~G4ExtrudedSolid()
183{
184 // Destructor
185}
186
187//_____________________________________________________________________________
188
189void G4ExtrudedSolid::ComputeProjectionParameters()
190{
191 // Compute parameters for point projections p(z)
192 // to the polygon scale & offset:
193 // scale(z) = k*z + scale0
194 // offset(z) = l*z + offset0
195 // p(z) = scale(z)*p0 + offset(z)
196 // p0 = (p(z) - offset(z))/scale(z);
197 //
198
199 for ( G4int iz=0; iz<fNz-1; ++iz)
200 {
201 G4double z1 = fZSections[iz].fZ;
202 G4double z2 = fZSections[iz+1].fZ;
203 G4double scale1 = fZSections[iz].fScale;
204 G4double scale2 = fZSections[iz+1].fScale;
205 G4TwoVector off1 = fZSections[iz].fOffset;
206 G4TwoVector off2 = fZSections[iz+1].fOffset;
207
208 G4double kscale = (scale2 - scale1)/(z2 - z1);
209 G4double scale0 = scale2 - kscale*(z2 - z1)/2.0;
210 G4TwoVector koff = (off2 - off1)/(z2 - z1);
211 G4TwoVector off0 = off2 - koff*(z2 - z1)/2.0;
212
213 fKScales.push_back(kscale);
214 fScale0s.push_back(scale0);
215 fKOffsets.push_back(koff);
216 fOffset0s.push_back(off0);
217 }
218}
219
220
221//_____________________________________________________________________________
222
223G4ThreeVector G4ExtrudedSolid::GetVertex(G4int iz, G4int ind) const
224{
225 // Shift and scale vertices
226
227 return G4ThreeVector( fPolygon[ind].x() * fZSections[iz].fScale
228 + fZSections[iz].fOffset.x(),
229 fPolygon[ind].y() * fZSections[iz].fScale
230 + fZSections[iz].fOffset.y(), fZSections[iz].fZ);
231}
232
233//_____________________________________________________________________________
234
235
236G4TwoVector G4ExtrudedSolid::ProjectPoint(const G4ThreeVector& point) const
237{
238 // Project point in the polygon scale
239 // scale(z) = k*z + scale0
240 // offset(z) = l*z + offset0
241 // p(z) = scale(z)*p0 + offset(z)
242 // p0 = (p(z) - offset(z))/scale(z);
243
244 // Select projection (z-segment of the solid) according to p.z()
245 //
246 G4int iz = 0;
247 while ( point.z() > fZSections[iz+1].fZ && iz < fNz-2 ) { ++iz; }
248
249 G4double z0 = ( fZSections[iz+1].fZ + fZSections[iz].fZ )/2.0;
250 G4TwoVector p2(point.x(), point.y());
251 G4double pscale = fKScales[iz]*(point.z()-z0) + fScale0s[iz];
252 G4TwoVector poffset = fKOffsets[iz]*(point.z()-z0) + fOffset0s[iz];
253
254 // G4cout << point << " projected to "
255 // << iz << "-th z-segment polygon as "
256 // << (p2 - poffset)/pscale << G4endl;
257
258 // pscale is always >0 as it is an interpolation between two
259 // positive scale values
260 //
261 return (p2 - poffset)/pscale;
262}
263
264//_____________________________________________________________________________
265
266G4bool G4ExtrudedSolid::IsSameLine(G4TwoVector p,
267 G4TwoVector l1, G4TwoVector l2) const
268{
269 // Return true if p is on the line through l1, l2
270
271 if ( l1.x() == l2.x() )
272 {
[850]273 return std::fabs(p.x() - l1.x()) < kCarTolerance * 0.5;
[831]274 }
275
276 return std::fabs (p.y() - l1.y() - ((l2.y() - l1.y())/(l2.x() - l1.x()))
[850]277 *(p.x() - l1.x())) < kCarTolerance * 0.5;
[831]278 }
279
280//_____________________________________________________________________________
281
282G4bool G4ExtrudedSolid::IsSameLineSegment(G4TwoVector p,
283 G4TwoVector l1, G4TwoVector l2) const
284{
285 // Return true if p is on the line through l1, l2 and lies between
286 // l1 and l2
287
[850]288 if ( p.x() < std::min(l1.x(), l2.x()) - kCarTolerance * 0.5 ||
289 p.x() > std::max(l1.x(), l2.x()) + kCarTolerance * 0.5 ||
290 p.y() < std::min(l1.y(), l2.y()) - kCarTolerance * 0.5 ||
291 p.y() > std::max(l1.y(), l2.y()) + kCarTolerance * 0.5 )
[831]292 {
293 return false;
294 }
295
296 return IsSameLine(p, l1, l2);
297}
298
299//_____________________________________________________________________________
300
301G4bool G4ExtrudedSolid::IsSameSide(G4TwoVector p1, G4TwoVector p2,
302 G4TwoVector l1, G4TwoVector l2) const
303{
304 // Return true if p1 and p2 are on the same side of the line through l1, l2
305
306 return ( (p1.x() - l1.x()) * (l2.y() - l1.y())
307 - (l2.x() - l1.x()) * (p1.y() - l1.y()) )
308 * ( (p2.x() - l1.x()) * (l2.y() - l1.y())
309 - (l2.x() - l1.x()) * (p2.y() - l1.y()) ) > 0;
310}
311
312//_____________________________________________________________________________
313
314G4bool G4ExtrudedSolid::IsPointInside(G4TwoVector a, G4TwoVector b,
315 G4TwoVector c, G4TwoVector p) const
316{
[850]317 // Return true if p is inside of triangle abc or on its edges,
318 // else returns false
[831]319
320 // Check extent first
321 //
322 if ( ( p.x() < a.x() && p.x() < b.x() && p.x() < c.x() ) ||
323 ( p.x() > a.x() && p.x() > b.x() && p.x() > c.x() ) ||
324 ( p.y() < a.y() && p.y() < b.y() && p.y() < c.y() ) ||
325 ( p.y() > a.y() && p.y() > b.y() && p.y() > c.y() ) ) return false;
326
[850]327 G4bool inside
328 = IsSameSide(p, a, b, c)
329 && IsSameSide(p, b, a, c)
330 && IsSameSide(p, c, a, b);
331
332 G4bool onEdge
333 = IsSameLineSegment(p, a, b)
334 || IsSameLineSegment(p, b, c)
335 || IsSameLineSegment(p, c, a);
336
337 return inside || onEdge;
[831]338}
339
340//_____________________________________________________________________________
341
342G4double
343G4ExtrudedSolid::GetAngle(G4TwoVector po, G4TwoVector pa, G4TwoVector pb) const
344{
345 // Return the angle of the vertex in po
346
347 G4TwoVector t1 = pa - po;
348 G4TwoVector t2 = pb - po;
349
[850]350 G4double result = (std::atan2(t1.y(), t1.x()) - std::atan2(t2.y(), t2.x()));
[831]351
352 if ( result < 0 ) result += 2*pi;
353
354 return result;
355}
356
357//_____________________________________________________________________________
358
359G4VFacet*
360G4ExtrudedSolid::MakeDownFacet(G4int ind1, G4int ind2, G4int ind3) const
361{
362 // Create a triangular facet from the polygon points given by indices
363 // forming the down side ( the normal goes in -z)
364
365 std::vector<G4ThreeVector> vertices;
366 vertices.push_back(GetVertex(0, ind1));
367 vertices.push_back(GetVertex(0, ind2));
368 vertices.push_back(GetVertex(0, ind3));
369
370 // first vertex most left
371 //
372 G4ThreeVector cross
373 = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
374
375 if ( cross.z() > 0.0 )
376 {
377 // vertices ardered clock wise has to be reordered
378
379 // G4cout << "G4ExtrudedSolid::MakeDownFacet: reordering vertices "
380 // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
381
382 G4ThreeVector tmp = vertices[1];
383 vertices[1] = vertices[2];
384 vertices[2] = tmp;
385 }
386
387 return new G4TriangularFacet(vertices[0], vertices[1],
388 vertices[2], ABSOLUTE);
389}
390
391//_____________________________________________________________________________
392
393G4VFacet*
394G4ExtrudedSolid::MakeUpFacet(G4int ind1, G4int ind2, G4int ind3) const
395{
396 // Creates a triangular facet from the polygon points given by indices
397 // forming the upper side ( z>0 )
398
399 std::vector<G4ThreeVector> vertices;
400 vertices.push_back(GetVertex(fNz-1, ind1));
401 vertices.push_back(GetVertex(fNz-1, ind2));
402 vertices.push_back(GetVertex(fNz-1, ind3));
403
404 // first vertex most left
405 //
406 G4ThreeVector cross
407 = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
408
409 if ( cross.z() < 0.0 )
410 {
411 // vertices ordered clock wise has to be reordered
412
413 // G4cout << "G4ExtrudedSolid::MakeUpFacet: reordering vertices "
414 // << ind1 << ", " << ind2 << ", " << ind3 << G4endl;
415
416 G4ThreeVector tmp = vertices[1];
417 vertices[1] = vertices[2];
418 vertices[2] = tmp;
419 }
420
421 return new G4TriangularFacet(vertices[0], vertices[1],
422 vertices[2], ABSOLUTE);
423}
424
425//_____________________________________________________________________________
426
427G4bool G4ExtrudedSolid::AddGeneralPolygonFacets()
428{
429 // Decompose polygonal sides in triangular facets
430
431 typedef std::pair < G4TwoVector, G4int > Vertex;
432
433 // Fill one more vector
434 //
435 std::vector< Vertex > verticesToBeDone;
436 for ( G4int i=0; i<fNv; ++i )
437 {
438 verticesToBeDone.push_back(Vertex(fPolygon[i], i));
439 }
440 std::vector< Vertex > ears;
441
442 std::vector< Vertex >::iterator c1 = verticesToBeDone.begin();
443 std::vector< Vertex >::iterator c2 = c1+1;
444 std::vector< Vertex >::iterator c3 = c1+2;
445 while ( verticesToBeDone.size()>2 )
446 {
447
448 // G4cout << "Looking at triangle : "
449 // << c1->second << " " << c2->second
450 // << " " << c3->second << G4endl;
451
452 // skip concave vertices
453 //
454 G4double angle = GetAngle(c2->first, c3->first, c1->first);
[850]455
456 if ( angle > pi )
457 {
[831]458 // G4cout << "Skipping concave vertex " << c2->second << G4endl;
459
460 // try next three consecutive vertices
461 //
462 c1 = c2;
463 c2 = c3;
464 ++c3;
465 if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
[850]466
467 // G4cout << "Looking at triangle : "
468 // << c1->second << " " << c2->second
469 // << " " << c3->second << G4endl;
470
[831]471 }
472
473 G4bool good = true;
474 std::vector< Vertex >::iterator it;
475 for ( it=verticesToBeDone.begin(); it != verticesToBeDone.end(); ++it )
476 {
477 // skip vertices of tested triangle
478 //
479 if ( it == c1 || it == c2 || it == c3 ) { continue; }
480
481 if ( IsPointInside(c1->first, c2->first, c3->first, it->first) )
482 {
483 // G4cout << "Point " << it->second << " is inside" << G4endl;
484 good = false;
485
486 // try next three consecutive vertices
487 //
488 c1 = c2;
489 c2 = c3;
490 ++c3;
491 if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
492 break;
493 }
494 // else
495 // { G4cout << "Point " << it->second << " is outside" << G4endl; }
496 }
497 if ( good )
498 {
499 // all points are outside triangle, we can make a facet
500
501 // G4cout << "Found triangle : "
502 // << c1->second << " " << c2->second
503 // << " " << c3->second << G4endl;
504
505 G4bool result;
506 result = AddFacet( MakeDownFacet(c1->second, c2->second, c3->second) );
507 if ( ! result ) { return false; }
508
509 result = AddFacet( MakeUpFacet(c1->second, c2->second, c3->second) );
510 if ( ! result ) { return false; }
511
512 std::vector<G4int> triangle(3);
513 triangle[0] = c1->second;
514 triangle[1] = c2->second;
515 triangle[2] = c3->second;
516 fTriangles.push_back(triangle);
517
518 // remove the ear point from verticesToBeDone
519 //
520 verticesToBeDone.erase(c2);
521 c1 = verticesToBeDone.begin();
522 c2 = c1+1;
523 c3 = c1+2;
524 }
525 }
526 return true;
527}
528
529//_____________________________________________________________________________
530
531G4bool G4ExtrudedSolid::MakeFacets()
532{
533 // Define facets
534
535 G4bool good;
536
537 // Decomposition of polygonal sides in the facets
538 //
539 if ( fNv == 3 )
540 {
541 good = AddFacet( new G4TriangularFacet( GetVertex(0, 0), GetVertex(0, 1),
542 GetVertex(0, 2), ABSOLUTE) );
543 if ( ! good ) { return false; }
544
545 good = AddFacet( new G4TriangularFacet( GetVertex(fNz-1, 2), GetVertex(fNz-1, 1),
546 GetVertex(fNz-1, 0), ABSOLUTE) );
547 if ( ! good ) { return false; }
548
549 std::vector<G4int> triangle(3);
550 triangle[0] = 0;
551 triangle[1] = 1;
552 triangle[2] = 2;
553 fTriangles.push_back(triangle);
554 }
555
556 else if ( fNv == 4 )
557 {
558 good = AddFacet( new G4QuadrangularFacet( GetVertex(0, 0),GetVertex(0, 1),
559 GetVertex(0, 2),GetVertex(0, 3),
560 ABSOLUTE) );
561 if ( ! good ) { return false; }
562
563 good = AddFacet( new G4QuadrangularFacet( GetVertex(fNz-1, 3), GetVertex(fNz-1, 2),
[921]564 GetVertex(fNz-1, 1), GetVertex(fNz-1, 0),
[831]565 ABSOLUTE) );
566 if ( ! good ) { return false; }
567
568 std::vector<G4int> triangle1(3);
569 triangle1[0] = 0;
570 triangle1[1] = 1;
571 triangle1[2] = 2;
572 fTriangles.push_back(triangle1);
573
574 std::vector<G4int> triangle2(3);
575 triangle2[0] = 0;
576 triangle2[1] = 2;
577 triangle2[2] = 3;
578 fTriangles.push_back(triangle2);
579 }
580 else
581 {
582 good = AddGeneralPolygonFacets();
583 if ( ! good ) { return false; }
584 }
585
586 // The quadrangular sides
587 //
588 for ( G4int iz = 0; iz < fNz-1; ++iz )
589 {
590 for ( G4int i = 0; i < fNv; ++i )
591 {
592 G4int j = (i+1) % fNv;
593 good = AddFacet( new G4QuadrangularFacet
594 ( GetVertex(iz, j), GetVertex(iz, i),
595 GetVertex(iz+1, i), GetVertex(iz+1, j), ABSOLUTE) );
596 if ( ! good ) { return false; }
597 }
598 }
599
600 SetSolidClosed(true);
601
602 return good;
603}
604
605//_____________________________________________________________________________
606
607G4bool G4ExtrudedSolid::IsConvex() const
608{
609 // Get polygon convexity (polygon is convex if all vertex angles are < pi )
610
611 for ( G4int i=0; i< fNv; ++i )
612 {
613 G4int j = ( i + 1 ) % fNv;
614 G4int k = ( i + 2 ) % fNv;
615 G4TwoVector v1 = fPolygon[i]-fPolygon[j];
616 G4TwoVector v2 = fPolygon[k]-fPolygon[j];
617 G4double dphi = v2.phi() - v1.phi();
618 if ( dphi < 0. ) { dphi += 2.*pi; }
619
620 if ( dphi >= pi ) { return false; }
621 }
622
623 return true;
624}
625
626//_____________________________________________________________________________
627
628G4GeometryType G4ExtrudedSolid::GetEntityType () const
629{
630 // Return entity type
631
632 return fGeometryType;
633}
634
635//_____________________________________________________________________________
636
637EInside G4ExtrudedSolid::Inside (const G4ThreeVector &p) const
638{
639 // Override the base class function as it fails in case of concave polygon.
640 // Project the point in the original polygon scale and check if it is inside
641 // for each triangle.
642
643 // Check first if outside extent
644 //
[850]645 if ( p.x() < GetMinXExtent() - kCarTolerance * 0.5 ||
646 p.x() > GetMaxXExtent() + kCarTolerance * 0.5 ||
647 p.y() < GetMinYExtent() - kCarTolerance * 0.5 ||
648 p.y() > GetMaxYExtent() + kCarTolerance * 0.5 ||
649 p.z() < GetMinZExtent() - kCarTolerance * 0.5 ||
650 p.z() > GetMaxZExtent() + kCarTolerance * 0.5 )
[831]651 {
652 // G4cout << "G4ExtrudedSolid::Outside extent: " << p << G4endl;
653 return kOutside;
654 }
655
656 // Project point p(z) to the polygon scale p0
657 //
658 G4TwoVector pscaled = ProjectPoint(p);
659
660 // Check if on surface of polygon
661 //
662 for ( G4int i=0; i<fNv; ++i )
663 {
664 G4int j = (i+1) % fNv;
665 if ( IsSameLine(pscaled, fPolygon[i], fPolygon[j]) )
666 {
667 // G4cout << "G4ExtrudedSolid::Inside return Surface (on polygon) "
668 // << G4endl;
669
670 return kSurface;
671 }
672 }
673
674 // Now check if inside triangles
675 //
676 std::vector< std::vector<G4int> >::const_iterator it = fTriangles.begin();
677 G4bool inside = false;
678 do
679 {
680 if ( IsPointInside(fPolygon[(*it)[0]], fPolygon[(*it)[1]],
681 fPolygon[(*it)[2]], pscaled) ) { inside = true; }
682 ++it;
683 } while ( (inside == false) && (it != fTriangles.end()) );
684
685 if ( inside )
686 {
687 // Check if on surface of z sides
688 //
[850]689 if ( std::fabs( p.z() - fZSections[0].fZ ) < kCarTolerance * 0.5 ||
690 std::fabs( p.z() - fZSections[fNz-1].fZ ) < kCarTolerance * 0.5 )
[831]691 {
692 // G4cout << "G4ExtrudedSolid::Inside return Surface (on z side)"
693 // << G4endl;
694
695 return kSurface;
696 }
697
698 // G4cout << "G4ExtrudedSolid::Inside return Inside" << G4endl;
699
700 return kInside;
701 }
702
703 // G4cout << "G4ExtrudedSolid::Inside return Outside " << G4endl;
704
705 return kOutside;
706}
707
708//_____________________________________________________________________________
709
710G4double G4ExtrudedSolid::DistanceToOut (const G4ThreeVector &p,
711 const G4ThreeVector &v,
712 const G4bool calcNorm,
713 G4bool *validNorm,
714 G4ThreeVector *n) const
715{
716 // Override the base class function to redefine validNorm
717 // (the solid can be concave)
718
719 G4double distOut =
720 G4TessellatedSolid::DistanceToOut(p, v, calcNorm, validNorm, n);
721 if (validNorm) { *validNorm = fIsConvex; }
722
723 return distOut;
724}
725
726
727//_____________________________________________________________________________
728
729G4double G4ExtrudedSolid::DistanceToOut (const G4ThreeVector &p) const
730{
731 // Override the overloaded base class function
732
733 return G4TessellatedSolid::DistanceToOut(p);
734}
735
736
737//_____________________________________________________________________________
738
739std::ostream& G4ExtrudedSolid::StreamInfo(std::ostream &os) const
740{
741 os << "-----------------------------------------------------------\n"
742 << " *** Dump for solid - " << GetName() << " ***\n"
743 << " ===================================================\n"
744 << " Solid geometry type: " << fGeometryType << G4endl;
745
746 if ( fIsConvex)
747 { os << " Convex polygon; list of vertices:" << G4endl; }
748 else
749 { os << " Concave polygon; list of vertices:" << G4endl; }
750
751 for ( G4int i=0; i<fNv; ++i )
752 {
753 os << " vx = " << fPolygon[i].x()/mm << " mm"
754 << " vy = " << fPolygon[i].y()/mm << " mm" << G4endl;
755 }
756
757 os << " Sections:" << G4endl;
758 for ( G4int iz=0; iz<fNz; ++iz )
759 {
760 os << " z = " << fZSections[iz].fZ/mm << " mm "
761 << " x0= " << fZSections[iz].fOffset.x()/mm << " mm "
762 << " y0= " << fZSections[iz].fOffset.y()/mm << " mm "
763 << " scale= " << fZSections[iz].fScale << G4endl;
764 }
765
766 return os;
767}
Note: See TracBrowser for help on using the repository browser.