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