source: trunk/source/geometry/solids/specific/src/G4VCSGfaceted.cc@ 837

Last change on this file since 837 was 831, checked in by garnier, 17 years ago

import all except CVS

File size: 12.6 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// the GEANT4 collaboration.
27//
28// By copying, distributing or modifying the Program (or any work
29// based on the Program) you indicate your acceptance of this statement,
30// and all its terms.
31//
32// $Id: G4VCSGfaceted.cc,v 1.20 2006/10/20 14:21:36 gcosmo Exp $
33// GEANT4 tag $Name: $
34//
35//
36// --------------------------------------------------------------------
37// GEANT 4 class source file
38//
39//
40// G4VCSGfaceted.cc
41//
42// Implementation of the virtual class of a CSG type shape that is built
43// entirely out of G4VCSGface faces.
44//
45// --------------------------------------------------------------------
46
47#include "G4VCSGfaceted.hh"
48#include "G4VCSGface.hh"
49#include "G4SolidExtentList.hh"
50
51#include "G4VoxelLimits.hh"
52#include "G4AffineTransform.hh"
53
54#include "G4Polyhedron.hh"
55#include "G4VGraphicsScene.hh"
56#include "G4NURBS.hh"
57#include "G4NURBSbox.hh"
58#include "G4VisExtent.hh"
59
60//
61// Constructor
62//
63G4VCSGfaceted::G4VCSGfaceted( const G4String& name )
64 : G4VSolid(name),
65 numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
66 fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
67{
68}
69
70
71//
72// Fake default constructor - sets only member data and allocates memory
73// for usage restricted to object persistency.
74//
75G4VCSGfaceted::G4VCSGfaceted( __void__& a )
76 : G4VSolid(a),
77 numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
78 fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
79{
80}
81
82//
83// Destructor
84//
85G4VCSGfaceted::~G4VCSGfaceted()
86{
87 DeleteStuff();
88 delete fpPolyhedron;
89}
90
91
92//
93// Copy constructor
94//
95G4VCSGfaceted::G4VCSGfaceted( const G4VCSGfaceted &source )
96 : G4VSolid( source )
97{
98 CopyStuff( source );
99}
100
101
102//
103// Assignment operator
104//
105const G4VCSGfaceted &G4VCSGfaceted::operator=( const G4VCSGfaceted &source )
106{
107 if (&source == this) return *this;
108
109 DeleteStuff();
110 CopyStuff( source );
111
112 return *this;
113}
114
115
116//
117// CopyStuff (protected)
118//
119// Copy the contents of source
120//
121void G4VCSGfaceted::CopyStuff( const G4VCSGfaceted &source )
122{
123 numFace = source.numFace;
124 if (numFace == 0) return; // odd, but permissable?
125
126 faces = new G4VCSGface*[numFace];
127
128 G4VCSGface **face = faces,
129 **sourceFace = source.faces;
130 do {
131 *face = (*sourceFace)->Clone();
132 } while( ++sourceFace, ++face < faces+numFace );
133 fCubicVolume = source.fCubicVolume;
134 fpPolyhedron = source.fpPolyhedron;
135}
136
137
138//
139// DeleteStuff (protected)
140//
141// Delete all allocated objects
142//
143void G4VCSGfaceted::DeleteStuff()
144{
145 if (numFace)
146 {
147 G4VCSGface **face = faces;
148 do {
149 delete *face;
150 } while( ++face < faces + numFace );
151
152 delete [] faces;
153 }
154}
155
156
157//
158// CalculateExtent
159//
160G4bool G4VCSGfaceted::CalculateExtent( const EAxis axis,
161 const G4VoxelLimits &voxelLimit,
162 const G4AffineTransform &transform,
163 G4double &min,
164 G4double &max ) const
165{
166 G4SolidExtentList extentList( axis, voxelLimit );
167
168 //
169 // Loop over all faces, checking min/max extent as we go.
170 //
171 G4VCSGface **face = faces;
172 do {
173 (*face)->CalculateExtent( axis, voxelLimit, transform, extentList );
174 } while( ++face < faces + numFace );
175
176 //
177 // Return min/max value
178 //
179 return extentList.GetExtent( min, max );
180}
181
182
183//
184// Inside
185//
186// It could be a good idea to override this virtual
187// member to add first a simple test (such as spherical
188// test or whatnot) and to call this version only if
189// the simplier test fails.
190//
191EInside G4VCSGfaceted::Inside( const G4ThreeVector &p ) const
192{
193 EInside answer=kOutside;
194 G4VCSGface **face = faces;
195 G4double best = kInfinity;
196 do {
197 G4double distance;
198 EInside result = (*face)->Inside( p, kCarTolerance/2, &distance );
199 if (result == kSurface) return kSurface;
200 if (distance < best)
201 {
202 best = distance;
203 answer = result;
204 }
205 } while( ++face < faces + numFace );
206
207 return answer;
208}
209
210
211//
212// SurfaceNormal
213//
214G4ThreeVector G4VCSGfaceted::SurfaceNormal( const G4ThreeVector& p ) const
215{
216 G4ThreeVector answer;
217 G4VCSGface **face = faces;
218 G4double best = kInfinity;
219 do {
220 G4double distance;
221 G4ThreeVector normal = (*face)->Normal( p, &distance );
222 if (distance < best)
223 {
224 best = distance;
225 answer = normal;
226 }
227 } while( ++face < faces + numFace );
228
229 return answer;
230}
231
232
233//
234// DistanceToIn(p,v)
235//
236G4double G4VCSGfaceted::DistanceToIn( const G4ThreeVector &p,
237 const G4ThreeVector &v ) const
238{
239 G4double distance = kInfinity;
240 G4double distFromSurface = kInfinity;
241 G4VCSGface *bestFace=0;
242 G4VCSGface **face = faces;
243 do {
244 G4double faceDistance,
245 faceDistFromSurface;
246 G4ThreeVector faceNormal;
247 G4bool faceAllBehind;
248 if ((*face)->Intersect( p, v, false, kCarTolerance/2,
249 faceDistance, faceDistFromSurface,
250 faceNormal, faceAllBehind ) )
251 {
252 //
253 // Intersecting face
254 //
255 if (faceDistance < distance)
256 {
257 distance = faceDistance;
258 distFromSurface = faceDistFromSurface;
259 bestFace = *face;
260 if (distFromSurface <= 0) return 0;
261 }
262 }
263 } while( ++face < faces + numFace );
264
265 if (distance < kInfinity && distFromSurface<kCarTolerance/2)
266 {
267 if (bestFace->Distance(p,false) < kCarTolerance/2) distance = 0;
268 }
269
270 return distance;
271}
272
273
274//
275// DistanceToIn(p)
276//
277G4double G4VCSGfaceted::DistanceToIn( const G4ThreeVector &p ) const
278{
279 return DistanceTo( p, false );
280}
281
282
283//
284// DistanceToOut(p,v)
285//
286G4double G4VCSGfaceted::DistanceToOut( const G4ThreeVector &p,
287 const G4ThreeVector &v,
288 const G4bool calcNorm,
289 G4bool *validNorm,
290 G4ThreeVector *n ) const
291{
292 G4bool allBehind = true;
293 G4double distance = kInfinity;
294 G4double distFromSurface = kInfinity;
295 G4ThreeVector normal;
296 G4VCSGface *bestFace=0;
297
298 G4VCSGface **face = faces;
299 do {
300 G4double faceDistance,
301 faceDistFromSurface;
302 G4ThreeVector faceNormal;
303 G4bool faceAllBehind;
304 if ((*face)->Intersect( p, v, true, kCarTolerance/2,
305 faceDistance, faceDistFromSurface,
306 faceNormal, faceAllBehind ) )
307 {
308 //
309 // Intersecting face
310 //
311 if ( (distance < kInfinity) || (!faceAllBehind) ) allBehind = false;
312 if (faceDistance < distance)
313 {
314 distance = faceDistance;
315 distFromSurface = faceDistFromSurface;
316 normal = faceNormal;
317 bestFace = *face;
318 if (distFromSurface <= 0) break;
319 }
320 }
321 } while( ++face < faces + numFace );
322
323 if (distance < kInfinity)
324 {
325 if (distFromSurface <= 0)
326 distance = 0;
327 else if (distFromSurface<kCarTolerance/2)
328 {
329 if (bestFace->Distance(p,true) < kCarTolerance/2) distance = 0;
330 }
331
332 if (calcNorm)
333 {
334 *validNorm = allBehind;
335 *n = normal;
336 }
337 }
338 else
339 {
340 if (Inside(p) == kSurface) distance = 0;
341 if (calcNorm) *validNorm = false;
342 }
343
344 return distance;
345}
346
347
348//
349// DistanceToOut(p)
350//
351G4double G4VCSGfaceted::DistanceToOut( const G4ThreeVector &p ) const
352{
353 return DistanceTo( p, true );
354}
355
356
357//
358// DistanceTo
359//
360// Protected routine called by DistanceToIn and DistanceToOut
361//
362G4double G4VCSGfaceted::DistanceTo( const G4ThreeVector &p,
363 const G4bool outgoing ) const
364{
365 G4VCSGface **face = faces;
366 G4double best = kInfinity;
367 do {
368 G4double distance = (*face)->Distance( p, outgoing );
369 if (distance < best) best = distance;
370 } while( ++face < faces + numFace );
371
372 return (best < 0.5*kCarTolerance) ? 0 : best;
373}
374
375
376//
377// DescribeYourselfTo
378//
379void G4VCSGfaceted::DescribeYourselfTo( G4VGraphicsScene& scene ) const
380{
381 scene.AddSolid( *this );
382}
383
384
385//
386// GetExtent
387//
388// Define the sides of the box into which our solid instance would fit.
389//
390G4VisExtent G4VCSGfaceted::GetExtent() const
391{
392 static const G4ThreeVector xMax(1,0,0), xMin(-1,0,0),
393 yMax(0,1,0), yMin(0,-1,0),
394 zMax(0,0,1), zMin(0,0,-1);
395 static const G4ThreeVector *axes[6] =
396 { &xMin, &xMax, &yMin, &yMax, &zMin, &zMax };
397
398 G4double answers[6] =
399 {-kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity};
400
401 G4VCSGface **face = faces;
402 do {
403 const G4ThreeVector **axis = axes+5 ;
404 G4double *answer = answers+5;
405 do {
406 G4double testFace = (*face)->Extent( **axis );
407 if (testFace > *answer) *answer = testFace;
408 }
409 while( --axis, --answer >= answers );
410
411 } while( ++face < faces + numFace );
412
413 return G4VisExtent( -answers[0], answers[1],
414 -answers[2], answers[3],
415 -answers[4], answers[5] );
416}
417
418
419//
420// GetEntityType
421//
422G4GeometryType G4VCSGfaceted::GetEntityType() const
423{
424 return G4String("G4CSGfaceted");
425}
426
427
428//
429// Stream object contents to an output stream
430//
431std::ostream& G4VCSGfaceted::StreamInfo( std::ostream& os ) const
432{
433 os << "-----------------------------------------------------------\n"
434 << " *** Dump for solid - " << GetName() << " ***\n"
435 << " ===================================================\n"
436 << " Solid type: G4VCSGfaceted\n"
437 << " Parameters: \n"
438 << " number of faces: " << numFace << "\n"
439 << "-----------------------------------------------------------\n";
440
441 return os;
442}
443
444
445//
446// GetCubVolStatistics
447//
448G4int G4VCSGfaceted::GetCubVolStatistics() const
449{
450 return fStatistics;
451}
452
453
454//
455// GetCubVolEpsilon
456//
457G4double G4VCSGfaceted::GetCubVolEpsilon() const
458{
459 return fCubVolEpsilon;
460}
461
462
463//
464// SetCubVolStatistics
465//
466void G4VCSGfaceted::SetCubVolStatistics(G4int st)
467{
468 fCubicVolume=0.;
469 fStatistics=st;
470}
471
472
473//
474// SetCubVolEpsilon
475//
476void G4VCSGfaceted::SetCubVolEpsilon(G4double ep)
477{
478 fCubicVolume=0.;
479 fCubVolEpsilon=ep;
480}
481
482
483//
484// GetAreaStatistics
485//
486G4int G4VCSGfaceted::GetAreaStatistics() const
487{
488 return fStatistics;
489}
490
491
492//
493// GetAreaAccuracy
494//
495G4double G4VCSGfaceted::GetAreaAccuracy() const
496{
497 return fAreaAccuracy;
498}
499
500
501//
502// SetAreaStatistics
503//
504void G4VCSGfaceted::SetAreaStatistics(G4int st)
505{
506 fSurfaceArea=0.;
507 fStatistics=st;
508}
509
510
511//
512// SetAreaAccuracy
513//
514void G4VCSGfaceted::SetAreaAccuracy(G4double ep)
515{
516 fSurfaceArea=0.;
517 fAreaAccuracy=ep;
518}
519
520
521//
522// GetCubicVolume
523//
524G4double G4VCSGfaceted::GetCubicVolume()
525{
526 if(fCubicVolume != 0.) ;
527 else fCubicVolume = EstimateCubicVolume(fStatistics,fCubVolEpsilon);
528 return fCubicVolume;
529}
530
531
532//
533// GetSurfaceArea
534//
535G4double G4VCSGfaceted::GetSurfaceArea()
536{
537 if(fSurfaceArea != 0.) ;
538 else fSurfaceArea = EstimateCubicVolume(fStatistics,fAreaAccuracy);
539 return fSurfaceArea;
540}
541
542
543//
544// GetPolyhedron
545//
546G4Polyhedron* G4VCSGfaceted::GetPolyhedron () const
547{
548 if (!fpPolyhedron ||
549 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
550 fpPolyhedron->GetNumberOfRotationSteps())
551 {
552 delete fpPolyhedron;
553 fpPolyhedron = CreatePolyhedron();
554 }
555 return fpPolyhedron;
556}
Note: See TracBrowser for help on using the repository browser.