Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ClippablePolygon.cc
Go to the documentation of this file.
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// G4ClippablePolygon implementation
27//
28// Includes code from G4VSolid (P.Kent, V.Grichine, J.Allison)
29// --------------------------------------------------------------------
30
31#include "G4ClippablePolygon.hh"
32
33#include "G4VoxelLimits.hh"
35
36// Constructor
37//
43
44// AddVertexInOrder
45//
47{
48 vertices.push_back( vertex );
49}
50
51// ClearAllVertices
52//
54{
55 vertices.clear();
56}
57
58// Clip
59//
61{
62 if (voxelLimit.IsLimited())
63 {
64 ClipAlongOneAxis( voxelLimit, kXAxis );
65 ClipAlongOneAxis( voxelLimit, kYAxis );
66 ClipAlongOneAxis( voxelLimit, kZAxis );
67 }
68
69 return (!vertices.empty());
70}
71
72// PartialClip
73//
74// Clip, while ignoring the indicated axis
75//
77 const EAxis IgnoreMe )
78{
79 if (voxelLimit.IsLimited())
80 {
81 if (IgnoreMe != kXAxis) { ClipAlongOneAxis( voxelLimit, kXAxis ); }
82 if (IgnoreMe != kYAxis) { ClipAlongOneAxis( voxelLimit, kYAxis ); }
83 if (IgnoreMe != kZAxis) { ClipAlongOneAxis( voxelLimit, kZAxis ); }
84 }
85
86 return (!vertices.empty());
87}
88
89// GetExtent
90//
92 G4double& min,
93 G4double& max ) const
94{
95 //
96 // Okay, how many entries do we have?
97 //
98 std::size_t noLeft = vertices.size();
99
100 //
101 // Return false if nothing is left
102 //
103 if (noLeft == 0) { return false; }
104
105 //
106 // Initialize min and max to our first vertex
107 //
108 min = max = vertices[0].operator()( axis );
109
110 //
111 // Compare to the rest
112 //
113 for( std::size_t i=1; i<noLeft; ++i )
114 {
115 G4double component = vertices[i].operator()( axis );
116 if (component < min )
117 {
118 min = component;
119 }
120 else if (component > max )
121 {
122 max = component;
123 }
124 }
125
126 return true;
127}
128
129// GetMinPoint
130//
131// Returns pointer to minimum point along the specified axis.
132// Take care! Do not use pointer after destroying parent polygon.
133//
135{
136 std::size_t noLeft = vertices.size();
137 if (noLeft==0)
138 {
139 G4Exception("G4ClippablePolygon::GetMinPoint()",
140 "GeomSolids0002", FatalException, "Empty polygon.");
141 }
142
143 const G4ThreeVector *answer = &(vertices[0]);
144 G4double min = answer->operator()(axis);
145
146 for( std::size_t i=1; i<noLeft; ++i )
147 {
148 G4double component = vertices[i].operator()( axis );
149 if (component < min)
150 {
151 answer = &(vertices[i]);
152 min = component;
153 }
154 }
155
156 return answer;
157}
158
159// GetMaxPoint
160//
161// Returns pointer to maximum point along the specified axis.
162// Take care! Do not use pointer after destroying parent polygon.
163//
165{
166 std::size_t noLeft = vertices.size();
167 if (noLeft==0)
168 {
169 G4Exception("G4ClippablePolygon::GetMaxPoint()",
170 "GeomSolids0002", FatalException, "Empty polygon.");
171 }
172
173 const G4ThreeVector *answer = &(vertices[0]);
174 G4double max = answer->operator()(axis);
175
176 for( std::size_t i=1; i<noLeft; ++i )
177 {
178 G4double component = vertices[i].operator()( axis );
179 if (component > max)
180 {
181 answer = &(vertices[i]);
182 max = component;
183 }
184 }
185
186 return answer;
187}
188
189// InFrontOf
190//
191// Decide if this polygon is in "front" of another when
192// viewed along the specified axis. For our purposes here,
193// it is sufficient to use the minimum extent of the
194// polygon along the axis to determine this.
195//
196// In case the minima of the two polygons are equal,
197// we use a more sophisticated test.
198//
199// Note that it is possible for the two following
200// statements to both return true or both return false:
201// polygon1.InFrontOf(polygon2)
202// polygon2.BehindOf(polygon1)
203//
205 EAxis axis ) const
206{
207 //
208 // If things are empty, do something semi-sensible
209 //
210 std::size_t noLeft = vertices.size();
211 if (noLeft==0) { return false; }
212
213 if (other.Empty()) { return true; }
214
215 //
216 // Get minimum of other polygon
217 //
218 const G4ThreeVector *minPointOther = other.GetMinPoint( axis );
219 const G4double minOther = minPointOther->operator()(axis);
220
221 //
222 // Get minimum of this polygon
223 //
224 const G4ThreeVector *minPoint = GetMinPoint( axis );
225 const G4double min = minPoint->operator()(axis);
226
227 //
228 // Easy decision
229 //
230 if (min < minOther-kCarTolerance) { return true; } // Clear winner
231
232 if (minOther < min-kCarTolerance) { return false; } // Clear loser
233
234 //
235 // We have a tie (this will not be all that rare since our
236 // polygons are connected)
237 //
238 // Check to see if there is a vertex in the other polygon
239 // that is behind this one (or vice versa)
240 //
241 G4bool answer;
242 G4ThreeVector normalOther = other.GetNormal();
243
244 if (std::fabs(normalOther(axis)) > std::fabs(normal(axis)))
245 {
246 G4double minP, maxP;
247 GetPlanerExtent( *minPointOther, normalOther, minP, maxP );
248
249 answer = (normalOther(axis) > 0) ? (minP < -kCarTolerance)
250 : (maxP > +kCarTolerance);
251 }
252 else
253 {
254 G4double minP, maxP;
255 other.GetPlanerExtent( *minPoint, normal, minP, maxP );
256
257 answer = (normal(axis) > 0) ? (maxP > +kCarTolerance)
258 : (minP < -kCarTolerance);
259 }
260 return answer;
261}
262
263// BehindOf
264//
265// Decide if this polygon is behind another.
266// See notes in method "InFrontOf"
267//
269 EAxis axis ) const
270{
271 //
272 // If things are empty, do something semi-sensible
273 //
274 std::size_t noLeft = vertices.size();
275 if (noLeft==0) { return false; }
276
277 if (other.Empty()) { return true; }
278
279 //
280 // Get minimum of other polygon
281 //
282 const G4ThreeVector *maxPointOther = other.GetMaxPoint( axis );
283 const G4double maxOther = maxPointOther->operator()(axis);
284
285 //
286 // Get minimum of this polygon
287 //
288 const G4ThreeVector *maxPoint = GetMaxPoint( axis );
289 const G4double max = maxPoint->operator()(axis);
290
291 //
292 // Easy decision
293 //
294 if (max > maxOther+kCarTolerance) { return true; } // Clear winner
295
296 if (maxOther > max+kCarTolerance) { return false; } // Clear loser
297
298 //
299 // We have a tie (this will not be all that rare since our
300 // polygons are connected)
301 //
302 // Check to see if there is a vertex in the other polygon
303 // that is in front of this one (or vice versa)
304 //
305 G4bool answer;
306 G4ThreeVector normalOther = other.GetNormal();
307
308 if (std::fabs(normalOther(axis)) > std::fabs(normal(axis)))
309 {
310 G4double minP, maxP;
311 GetPlanerExtent( *maxPointOther, normalOther, minP, maxP );
312
313 answer = (normalOther(axis) > 0) ? (maxP > +kCarTolerance)
314 : (minP < -kCarTolerance);
315 }
316 else
317 {
318 G4double minP, maxP;
319 other.GetPlanerExtent( *maxPoint, normal, minP, maxP );
320
321 answer = (normal(axis) > 0) ? (minP < -kCarTolerance)
322 : (maxP > +kCarTolerance);
323 }
324 return answer;
325}
326
327// GetPlanerExtent
328//
329// Get min/max distance in or out of a plane
330//
332 const G4ThreeVector& planeNormal,
333 G4double& min,
334 G4double& max ) const
335{
336 //
337 // Okay, how many entries do we have?
338 //
339 std::size_t noLeft = vertices.size();
340
341 //
342 // Return false if nothing is left
343 //
344 if (noLeft == 0) { return false; }
345
346 //
347 // Initialize min and max to our first vertex
348 //
349 min = max = planeNormal.dot(vertices[0]-pointOnPlane);
350
351 //
352 // Compare to the rest
353 //
354 for( std::size_t i=1; i<noLeft; ++i )
355 {
356 G4double component = planeNormal.dot(vertices[i] - pointOnPlane);
357 if (component < min )
358 {
359 min = component;
360 }
361 else if (component > max )
362 {
363 max = component;
364 }
365 }
366
367 return true;
368}
369
370// ClipAlongOneAxis
371//
372// Clip along just one axis, as specified in voxelLimit
373//
375 const EAxis axis )
376{
377 if (!voxelLimit.IsLimited(axis)) { return; }
378
379 G4ThreeVectorList tempPolygon;
380
381 //
382 // Build a "simple" voxelLimit that includes only the min extent
383 // and apply this to our vertices, producing result in tempPolygon
384 //
385 G4VoxelLimits simpleLimit1;
386 simpleLimit1.AddLimit( axis, voxelLimit.GetMinExtent(axis), kInfinity );
387 ClipToSimpleLimits( vertices, tempPolygon, simpleLimit1 );
388
389 //
390 // If nothing is left from the above clip, we might as well return now
391 // (but with an empty vertices)
392 //
393 if (tempPolygon.empty())
394 {
395 vertices.clear();
396 return;
397 }
398
399 //
400 // Now do the same, but using a "simple" limit that includes only the max
401 // extent. Apply this to out tempPolygon, producing result in vertices.
402 //
403 G4VoxelLimits simpleLimit2;
404 simpleLimit2.AddLimit( axis, -kInfinity, voxelLimit.GetMaxExtent(axis) );
405 ClipToSimpleLimits( tempPolygon, vertices, simpleLimit2 );
406
407 //
408 // If nothing is left, return now
409 //
410 if (vertices.empty()) { return; }
411}
412
413// ClipToSimpleLimits
414//
415// pVoxelLimits must be only limited along one axis, and either the maximum
416// along the axis must be +kInfinity, or the minimum -kInfinity
417//
418void G4ClippablePolygon::ClipToSimpleLimits( G4ThreeVectorList& pPolygon,
419 G4ThreeVectorList& outputPolygon,
420 const G4VoxelLimits& pVoxelLimit )
421{
422 std::size_t noVertices = pPolygon.size();
423 G4ThreeVector vEnd,vStart;
424
425 outputPolygon.clear();
426
427 for (std::size_t i=0; i<noVertices; ++i)
428 {
429 vStart=pPolygon[i];
430 if (i==noVertices-1)
431 {
432 vEnd=pPolygon[0];
433 }
434 else
435 {
436 vEnd=pPolygon[i+1];
437 }
438
439 if (pVoxelLimit.Inside(vStart))
440 {
441 if (pVoxelLimit.Inside(vEnd))
442 {
443 // vStart and vEnd inside -> output end point
444 //
445 outputPolygon.push_back(vEnd);
446 }
447 else
448 {
449 // vStart inside, vEnd outside -> output crossing point
450 //
451 pVoxelLimit.ClipToLimits(vStart,vEnd);
452 outputPolygon.push_back(vEnd);
453 }
454 }
455 else
456 {
457 if (pVoxelLimit.Inside(vEnd))
458 {
459 // vStart outside, vEnd inside -> output inside section
460 //
461 pVoxelLimit.ClipToLimits(vStart,vEnd);
462 outputPolygon.push_back(vStart);
463 outputPolygon.push_back(vEnd);
464 }
465 else // Both point outside -> no output
466 {
467 }
468 }
469 }
470}
std::vector< G4ThreeVector > G4ThreeVectorList
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
double dot(const Hep3Vector &) const
G4bool GetPlanerExtent(const G4ThreeVector &pointOnPlane, const G4ThreeVector &planeNormal, G4double &min, G4double &max) const
G4bool Clip(const G4VoxelLimits &voxelLimit)
G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
void AddVertexInOrder(const G4ThreeVector &vertex)
const G4ThreeVector * GetMinPoint(const EAxis axis) const
G4bool GetExtent(const EAxis axis, G4double &min, G4double &max) const
void ClipAlongOneAxis(const G4VoxelLimits &voxelLimit, const EAxis axis)
G4bool InFrontOf(const G4ClippablePolygon &other, EAxis axis) const
G4bool BehindOf(const G4ClippablePolygon &other, EAxis axis) const
G4bool Empty() const
const G4ThreeVector GetNormal() const
const G4ThreeVector * GetMaxPoint(const EAxis axis) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
G4double GetMinExtent(const EAxis pAxis) const
G4bool ClipToLimits(G4ThreeVector &pStart, G4ThreeVector &pEnd) const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4double GetMaxExtent(const EAxis pAxis) const
G4bool Inside(const G4ThreeVector &pVec) const
G4bool IsLimited() const
EAxis
Definition geomdefs.hh:54
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668