source: trunk/source/geometry/management/src/G4VSolid.cc@ 900

Last change on this file since 900 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 20.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//
[850]27// $Id: G4VSolid.cc,v 1.37 2008/02/20 15:24:26 gcosmo Exp $
28// GEANT4 tag $Name: HEAD $
[831]29//
30// class G4VSolid
31//
32// Implementation for solid base class
33//
34// History:
35//
36// 06.12.02 V.Grichine, restored original conditions in ClipPolygon()
37// 10.05.02 V.Grichine, ClipPolygon(): clip only other axis and limited voxels
38// 15.04.02 V.Grichine, bug fixed in ClipPolygon(): clip only one axis
39// 13.03.02 V.Grichine, cosmetics of voxel limit functions
40// 15.11.00 D.Williams, V.Grichine, fix in CalculateClippedPolygonExtent()
41// 10.07.95 P.Kent, Added == operator, solid Store entry
42// 30.06.95 P.Kent, Created.
43// --------------------------------------------------------------------
44
45#include "G4VSolid.hh"
46#include "G4SolidStore.hh"
47#include "globals.hh"
48#include "Randomize.hh"
49#include "G4GeometryTolerance.hh"
50
51#include "G4VoxelLimits.hh"
52#include "G4AffineTransform.hh"
53#include "G4VisExtent.hh"
54
55//////////////////////////////////////////////////////////////////////////
56//
57// Constructor
58// - Copies name
59// - Add ourselves to solid Store
60
61G4VSolid::G4VSolid(const G4String& name)
62 : fshapeName(name)
63{
64 kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
65
66 // Register to store
67 //
68 G4SolidStore::GetInstance()->Register(this);
69}
70
71//////////////////////////////////////////////////////////////////////////
72//
73// Protected copy constructor
74//
75
76G4VSolid::G4VSolid(const G4VSolid& rhs)
77 : kCarTolerance(rhs.kCarTolerance), fshapeName(rhs.fshapeName)
78{
79 // Register to store
80 //
81 G4SolidStore::GetInstance()->Register(this);
82}
83
84//////////////////////////////////////////////////////////////////////////
85//
86// Fake default constructor - sets only member data and allocates memory
87// for usage restricted to object persistency.
88//
89G4VSolid::G4VSolid( __void__& )
90 : fshapeName("")
91{
92 // Register to store
93 //
94 G4SolidStore::GetInstance()->Register(this);
95}
96
97//////////////////////////////////////////////////////////////////////////
98//
99// Destructor (virtual)
100// - Remove ourselves from solid Store
101
102G4VSolid::~G4VSolid()
103{
104 G4SolidStore::GetInstance()->DeRegister(this);
105}
106
107//////////////////////////////////////////////////////////////////////////
108//
109// Assignment operator
110
111G4VSolid& G4VSolid::operator = (const G4VSolid& rhs)
112{
113 // Check assignment to self
114 //
115 if (this == &rhs) { return *this; }
116
117 // Copy data
118 //
119 kCarTolerance = rhs.kCarTolerance;
120 fshapeName = rhs.fshapeName;
121
122 return *this;
123}
124
125//////////////////////////////////////////////////////////////////////////
126//
127// Streaming operator dumping solid contents
128
129std::ostream& operator<< ( std::ostream& os, const G4VSolid& e )
130{
131 return e.StreamInfo(os);
132}
133
134//////////////////////////////////////////////////////////////////////////
135//
136// Throw exception if ComputeDimensions called for illegal derived class
137
138void G4VSolid::ComputeDimensions(G4VPVParameterisation*,
139 const G4int,
140 const G4VPhysicalVolume*)
141{
142 G4cerr << "ERROR - Illegal call to G4VSolid::ComputeDimensions()" << G4endl
143 << " Method not overloaded by derived class !" << G4endl;
144 G4Exception("G4VSolid::ComputeDimensions()", "NotApplicable",
145 FatalException, "Illegal call to case class.");
146}
147
148//////////////////////////////////////////////////////////////////////////
149//
150// Throw exception (warning) for solids not implementing the method
151
152G4ThreeVector G4VSolid::GetPointOnSurface() const
153{
154 G4cerr << "WARNING - G4VSolid::GetPointOnSurface()" << G4endl
155 << " Not implemented for solid: "
156 << this->GetEntityType() << " !" << G4endl;
157 G4Exception("G4VSolid::GetPointOnSurface()", "NotImplemented",
158 JustWarning, "Not implemented for this solid ! Returning origin.");
159 return G4ThreeVector(0,0,0);
160}
161
162///////////////////////////////////////////////////////////////////////////
163//
164// Calculate the maximum and minimum extents of the polygon described
165// by the vertices: pSectionIndex->pSectionIndex+1->
166// pSectionIndex+2->pSectionIndex+3->pSectionIndex
167// in the List pVertices
168//
169// If the minimum is <pMin pMin is set to the new minimum
170// If the maximum is >pMax pMax is set to the new maximum
171//
172// No modifications are made to pVertices
173//
174
175void G4VSolid::ClipCrossSection( G4ThreeVectorList* pVertices,
176 const G4int pSectionIndex,
177 const G4VoxelLimits& pVoxelLimit,
178 const EAxis pAxis,
179 G4double& pMin, G4double& pMax) const
180{
181
182 G4ThreeVectorList polygon;
183 polygon.push_back((*pVertices)[pSectionIndex]);
184 polygon.push_back((*pVertices)[pSectionIndex+1]);
185 polygon.push_back((*pVertices)[pSectionIndex+2]);
186 polygon.push_back((*pVertices)[pSectionIndex+3]);
187 // G4cout<<"ClipCrossSection: 0-1-2-3"<<G4endl;
188 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
189 return;
190}
191
192//////////////////////////////////////////////////////////////////////////////////
193//
194// Calculate the maximum and minimum extents of the polygons
195// joining the CrossSections at pSectionIndex->pSectionIndex+3 and
196// pSectionIndex+4->pSectionIndex7
197//
198// in the List pVertices, within the boundaries of the voxel limits pVoxelLimit
199//
200// If the minimum is <pMin pMin is set to the new minimum
201// If the maximum is >pMax pMax is set to the new maximum
202//
203// No modifications are made to pVertices
204
205void G4VSolid::ClipBetweenSections( G4ThreeVectorList* pVertices,
206 const G4int pSectionIndex,
207 const G4VoxelLimits& pVoxelLimit,
208 const EAxis pAxis,
209 G4double& pMin, G4double& pMax) const
210{
211 G4ThreeVectorList polygon;
212 polygon.push_back((*pVertices)[pSectionIndex]);
213 polygon.push_back((*pVertices)[pSectionIndex+4]);
214 polygon.push_back((*pVertices)[pSectionIndex+5]);
215 polygon.push_back((*pVertices)[pSectionIndex+1]);
216 // G4cout<<"ClipBetweenSections: 0-4-5-1"<<G4endl;
217 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
218 polygon.clear();
219
220 polygon.push_back((*pVertices)[pSectionIndex+1]);
221 polygon.push_back((*pVertices)[pSectionIndex+5]);
222 polygon.push_back((*pVertices)[pSectionIndex+6]);
223 polygon.push_back((*pVertices)[pSectionIndex+2]);
224 // G4cout<<"ClipBetweenSections: 1-5-6-2"<<G4endl;
225 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
226 polygon.clear();
227
228 polygon.push_back((*pVertices)[pSectionIndex+2]);
229 polygon.push_back((*pVertices)[pSectionIndex+6]);
230 polygon.push_back((*pVertices)[pSectionIndex+7]);
231 polygon.push_back((*pVertices)[pSectionIndex+3]);
232 // G4cout<<"ClipBetweenSections: 2-6-7-3"<<G4endl;
233 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
234 polygon.clear();
235
236 polygon.push_back((*pVertices)[pSectionIndex+3]);
237 polygon.push_back((*pVertices)[pSectionIndex+7]);
238 polygon.push_back((*pVertices)[pSectionIndex+4]);
239 polygon.push_back((*pVertices)[pSectionIndex]);
240 // G4cout<<"ClipBetweenSections: 3-7-4-0"<<G4endl;
241 CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
242 return;
243}
244
245
246///////////////////////////////////////////////////////////////////////////////
247//
248// Calculate the maximum and minimum extents of the convex polygon pPolygon
249// along the axis pAxis, within the limits pVoxelLimit
250//
251
252void
253G4VSolid::CalculateClippedPolygonExtent(G4ThreeVectorList& pPolygon,
254 const G4VoxelLimits& pVoxelLimit,
255 const EAxis pAxis,
256 G4double& pMin,
257 G4double& pMax) const
258{
259 G4int noLeft,i;
260 G4double component;
261 /*
262 G4cout<<G4endl;
263 for(i = 0 ; i < pPolygon.size() ; i++ )
264 {
265 G4cout << i << "\t"
266 << "p.x = " << pPolygon[i].operator()(pAxis) << "\t"
267 // << "p.y = " << pPolygon[i].y() << "\t"
268 // << "p.z = " << pPolygon[i].z() << "\t"
269 << G4endl;
270 }
271 G4cout<<G4endl;
272 */
273 ClipPolygon(pPolygon,pVoxelLimit,pAxis);
274 noLeft = pPolygon.size();
275
276 if ( noLeft )
277 {
278 // G4cout<<G4endl;
279 for (i=0;i<noLeft;i++)
280 {
281 component = pPolygon[i].operator()(pAxis);
282 // G4cout <<i<<"\t"<<component<<G4endl;
283
284 if (component < pMin)
285 {
286 // G4cout <<i<<"\t"<<"Pmin = "<<component<<G4endl;
287 pMin = component;
288 }
289 if (component > pMax)
290 {
291 // G4cout <<i<<"\t"<<"PMax = "<<component<<G4endl;
292 pMax = component;
293 }
294 }
295 // G4cout<<G4endl;
296 }
297 // G4cout<<"pMin = "<<pMin<<"\t"<<"pMax = "<<pMax<<G4endl;
298}
299
300/////////////////////////////////////////////////////////////////////////////
301//
302// Clip the convex polygon described by the vertices at
303// pSectionIndex ->pSectionIndex+3 within pVertices to the limits pVoxelLimit
304//
305// Set pMin to the smallest
306//
307// Calculate the extent of the polygon along pAxis, when clipped to the
308// limits pVoxelLimit. If the polygon exists after clippin, set pMin to
309// the polygon's minimum extent along the axis if <pMin, and set pMax to
310// the polygon's maximum extent along the axis if >pMax.
311//
312// The polygon is described by a set of vectors, where each vector represents
313// a vertex, so that the polygon is described by the vertex sequence:
314// 0th->1st 1st->2nd 2nd->... nth->0th
315//
316// Modifications to the polygon are made
317//
318// NOTE: Execessive copying during clipping
319
320void G4VSolid::ClipPolygon( G4ThreeVectorList& pPolygon,
321 const G4VoxelLimits& pVoxelLimit,
322 const EAxis ) const
323{
324 G4ThreeVectorList outputPolygon;
325
326 if ( pVoxelLimit.IsLimited() )
327 {
328 if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
329 {
330 G4VoxelLimits simpleLimit1;
331 simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
332 // G4cout<<"MinXExtent()"<<G4endl;
333 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
334
335 pPolygon.clear();
336
337 if ( !outputPolygon.size() ) return;
338
339 G4VoxelLimits simpleLimit2;
340 // G4cout<<"MaxXExtent()"<<G4endl;
341 simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
342 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
343
344 if ( !pPolygon.size() ) return;
345 else outputPolygon.clear();
346 }
347 if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
348 {
349 G4VoxelLimits simpleLimit1;
350 simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
351 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
352
353 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
354 // early exit
355
356 pPolygon.clear();
357
358 if ( !outputPolygon.size() ) return;
359
360 G4VoxelLimits simpleLimit2;
361 simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
362 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
363
364 if ( !pPolygon.size() ) return;
365 else outputPolygon.clear();
366 }
367 if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
368 {
369 G4VoxelLimits simpleLimit1;
370 simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
371 ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
372
373 // Must always clear pPolygon - for clip to simpleLimit2 and in case of
374 // early exit
375
376 pPolygon.clear();
377
378 if ( !outputPolygon.size() ) return;
379
380 G4VoxelLimits simpleLimit2;
381 simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
382 ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
383
384 // Return after final clip - no cleanup
385 }
386 }
387}
388
389////////////////////////////////////////////////////////////////////////////
390//
391// pVoxelLimits must be only limited along one axis, and either the maximum
392// along the axis must be +kInfinity, or the minimum -kInfinity
393
394void
395G4VSolid::ClipPolygonToSimpleLimits( G4ThreeVectorList& pPolygon,
396 G4ThreeVectorList& outputPolygon,
397 const G4VoxelLimits& pVoxelLimit ) const
398{
399 G4int i;
400 G4int noVertices=pPolygon.size();
401 G4ThreeVector vEnd,vStart;
402
403 for (i = 0 ; i < noVertices ; i++ )
404 {
405 vStart = pPolygon[i];
406 // G4cout << "i = " << i << G4endl;
407 if ( i == noVertices-1 ) vEnd = pPolygon[0];
408 else vEnd = pPolygon[i+1];
409
410 if ( pVoxelLimit.Inside(vStart) )
411 {
412 if (pVoxelLimit.Inside(vEnd))
413 {
414 // vStart and vEnd inside -> output end point
415 //
416 outputPolygon.push_back(vEnd);
417 }
418 else
419 {
420 // vStart inside, vEnd outside -> output crossing point
421 //
422 // G4cout << "vStart inside, vEnd outside" << G4endl;
423 pVoxelLimit.ClipToLimits(vStart,vEnd);
424 outputPolygon.push_back(vEnd);
425 }
426 }
427 else
428 {
429 if (pVoxelLimit.Inside(vEnd))
430 {
431 // vStart outside, vEnd inside -> output inside section
432 //
433 // G4cout << "vStart outside, vEnd inside" << G4endl;
434 pVoxelLimit.ClipToLimits(vStart,vEnd);
435 outputPolygon.push_back(vStart);
436 outputPolygon.push_back(vEnd);
437 }
438 else // Both point outside -> no output
439 {
440 // outputPolygon.push_back(vStart);
441 // outputPolygon.push_back(vEnd);
442 }
443 }
444 }
445}
446
447const G4VSolid* G4VSolid::GetConstituentSolid(G4int) const
448{ return 0; }
449
450G4VSolid* G4VSolid::GetConstituentSolid(G4int)
451{ return 0; }
452
453const G4DisplacedSolid* G4VSolid::GetDisplacedSolidPtr() const
454{ return 0; }
455
456G4DisplacedSolid* G4VSolid::GetDisplacedSolidPtr()
457{ return 0; }
458
459G4VisExtent G4VSolid::GetExtent () const
460{
461 G4VisExtent extent;
462 G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
463 G4AffineTransform affineTransform;
464 G4double vmin, vmax;
465 CalculateExtent(kXAxis,voxelLimits,affineTransform,vmin,vmax);
466 extent.SetXmin (vmin);
467 extent.SetXmax (vmax);
468 CalculateExtent(kYAxis,voxelLimits,affineTransform,vmin,vmax);
469 extent.SetYmin (vmin);
470 extent.SetYmax (vmax);
471 CalculateExtent(kZAxis,voxelLimits,affineTransform,vmin,vmax);
472 extent.SetZmin (vmin);
473 extent.SetZmax (vmax);
474 return extent;
475}
476
477G4Polyhedron* G4VSolid::CreatePolyhedron () const
478{
479 return 0;
480}
481
482G4NURBS* G4VSolid::CreateNURBS () const
483{
484 return 0;
485}
486
487G4Polyhedron* G4VSolid::GetPolyhedron () const
488{
489 return 0;
490}
491
492////////////////////////////////////////////////////////////////
493//
494// Returns an estimation of the solid volume in internal units.
495// The number of statistics and error accuracy is fixed.
496// This method may be overloaded by derived classes to compute the
497// exact geometrical quantity for solids where this is possible.
498// or anyway to cache the computed value.
499// This implementation does NOT cache the computed value.
500
501G4double G4VSolid::GetCubicVolume()
502{
503 G4int cubVolStatistics = 1000000;
504 G4double cubVolEpsilon = 0.001;
505 return EstimateCubicVolume(cubVolStatistics, cubVolEpsilon);
506}
507
508////////////////////////////////////////////////////////////////
509//
510// Calculate cubic volume based on Inside() method.
511// Accuracy is limited by the second argument or the statistics
512// expressed by the first argument.
513// Implementation is courtesy of Vasiliki Despoina Mitsou,
514// University of Athens.
515
516G4double G4VSolid::EstimateCubicVolume(G4int nStat, G4double epsilon) const
517{
518 G4int iInside=0;
519 G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume;
520 G4bool yesno;
521 G4ThreeVector p;
522 EInside in;
523
524 // values needed for CalculateExtent signature
525
526 G4VoxelLimits limit; // Unlimited
527 G4AffineTransform origin;
528
529 // min max extents of pSolid along X,Y,Z
530
531 yesno = this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
532 yesno = this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
533 yesno = this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
534
535 // limits
536
537 if(nStat < 100) nStat = 100;
538 if(epsilon > 0.01) epsilon = 0.01;
539
540 for(G4int i = 0; i < nStat; i++ )
541 {
542 px = minX+(maxX-minX)*G4UniformRand();
543 py = minY+(maxY-minY)*G4UniformRand();
544 pz = minZ+(maxZ-minZ)*G4UniformRand();
545 p = G4ThreeVector(px,py,pz);
546 in = this->Inside(p);
547 if(in != kOutside) iInside++;
548 }
549 volume = (maxX-minX)*(maxY-minY)*(maxZ-minZ)*iInside/nStat;
550 return volume;
551}
552
553////////////////////////////////////////////////////////////////
554//
555// Returns an estimation of the solid surface area in internal units.
556// The number of statistics and error accuracy is fixed.
557// This method may be overloaded by derived classes to compute the
558// exact geometrical quantity for solids where this is possible.
559// or anyway to cache the computed value.
560// This implementation does NOT cache the computed value.
561
562G4double G4VSolid::GetSurfaceArea()
563{
564 G4int stat = 1000000;
565 G4double ell = -1.;
566 return EstimateSurfaceArea(stat,ell);
567}
568
569////////////////////////////////////////////////////////////////
570//
571// Estimate surface area based on Inside(), DistanceToIn(), and
572// DistanceToOut() methods. Accuracy is limited by the statistics
573// defined by the first argument. Implemented by Mikhail Kosov.
574
575G4double G4VSolid::EstimateSurfaceArea(G4int nStat, G4double ell) const
576{
577 G4int inside=0;
578 G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,surf;
579 G4bool yesno;
580 G4ThreeVector p;
581 EInside in;
582
583 // values needed for CalculateExtent signature
584
585 G4VoxelLimits limit; // Unlimited
586 G4AffineTransform origin;
587
588 // min max extents of pSolid along X,Y,Z
589
590 yesno = this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
591 yesno = this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
592 yesno = this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
593
594 // limits
595
596 if(nStat < 100) { nStat = 100; }
597
598 G4double dX=maxX-minX;
599 G4double dY=maxY-minY;
600 G4double dZ=maxZ-minZ;
601 if(ell<=0.) // Automatic definition of skin thickness
602 {
603 G4double minval=dX;
604 if(dY<dX) { minval=dY; }
605 if(dZ<minval) { minval=dZ; }
606 ell=.01*minval;
607 }
608
609 G4double dd=2*ell;
610 minX-=ell; minY-=ell; minZ-=ell; dX+=dd; dY+=dd; dZ+=dd;
611
612 for(G4int i = 0; i < nStat; i++ )
613 {
614 px = minX+dX*G4UniformRand();
615 py = minY+dY*G4UniformRand();
616 pz = minZ+dZ*G4UniformRand();
617 p = G4ThreeVector(px,py,pz);
618 in = this->Inside(p);
619 if(in != kOutside)
620 {
621 if (DistanceToOut(p)<ell) { inside++; }
622 }
623 else if(DistanceToIn(p)<ell) { inside++; }
624 }
625 // @@ The conformal correction can be upgraded
626 surf = dX*dY*dZ*inside/dd/nStat;
627 return surf;
628}
Note: See TracBrowser for help on using the repository browser.