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

Last change on this file since 1316 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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