Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/Boolean/src/G4ScaledSolid.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // Implementation for G4ScaledSolid class
 27 //
 28 // 27.10.15 G.Cosmo: created, based on implementation also provided in Root
 29 // --------------------------------------------------------------------
 30 
 31 #include "G4ScaledSolid.hh"
 32 #include "G4BoundingEnvelope.hh"
 33 
 34 #include "G4VPVParameterisation.hh"
 35 
 36 #include "G4ScaleTransform.hh"
 37 
 38 #include "G4VGraphicsScene.hh"
 39 #include "G4Polyhedron.hh"
 40 
 41 ///////////////////////////////////////////////////////////////////
 42 //
 43 // Constructor
 44 //
 45 G4ScaledSolid::G4ScaledSolid( const G4String& pName,
 46                                     G4VSolid* pSolid,
 47                               const G4Scale3D& pScale )
 48   : G4VSolid(pName), fPtrSolid(pSolid)
 49 {
 50   fScale = new G4ScaleTransform(pScale);
 51 }
 52 
 53 ///////////////////////////////////////////////////////////////////
 54 //
 55 // Fake default constructor - sets only member data and allocates memory
 56 //                            for usage restricted to object persistency
 57 //
 58 G4ScaledSolid::G4ScaledSolid( __void__& a )
 59   : G4VSolid(a)
 60 {
 61 }
 62 
 63 ///////////////////////////////////////////////////////////////////
 64 //
 65 // Destructor
 66 //
 67 G4ScaledSolid::~G4ScaledSolid()
 68 {
 69   delete fpPolyhedron; fpPolyhedron = nullptr;
 70   delete fScale; fScale = nullptr;
 71 }
 72 
 73 ///////////////////////////////////////////////////////////////
 74 //
 75 // Copy constructor
 76 //
 77 G4ScaledSolid::G4ScaledSolid(const G4ScaledSolid& rhs)
 78   : G4VSolid (rhs), fPtrSolid(rhs.fPtrSolid),
 79     fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea)
 80 {
 81   fScale = new G4ScaleTransform(*(rhs.fScale));
 82 }
 83 
 84 ///////////////////////////////////////////////////////////////
 85 //
 86 // Assignment operator
 87 //
 88 G4ScaledSolid& G4ScaledSolid::operator = (const G4ScaledSolid& rhs)
 89 {
 90   // Check assignment to self
 91   //
 92   if (this == &rhs) { return *this; }
 93 
 94   // Copy base class data
 95   //
 96   G4VSolid::operator=(rhs);
 97 
 98   // Copy data
 99   //
100   fPtrSolid = rhs.fPtrSolid;
101   delete fScale;
102   fScale = new G4ScaleTransform(*(rhs.fScale));
103   fCubicVolume = rhs.fCubicVolume;
104   fSurfaceArea = rhs.fSurfaceArea;
105   fRebuildPolyhedron = false;
106   delete fpPolyhedron; fpPolyhedron = nullptr;
107 
108   return *this;
109 }
110 
111 //////////////////////////////////////////////////////////////////////////
112 //
113 // Return original solid not scaled
114 //
115 G4VSolid* G4ScaledSolid::GetUnscaledSolid() const
116 {
117   return fPtrSolid;
118 }
119 
120 //////////////////////////////////////////////////////////////////////////
121 //
122 // Get bounding box
123 //
124 void G4ScaledSolid::BoundingLimits(G4ThreeVector& pMin,
125                                    G4ThreeVector& pMax) const
126 {
127   G4ThreeVector bmin,bmax;
128   G4ThreeVector scale = fScale->GetScale();
129  
130   fPtrSolid->BoundingLimits(bmin,bmax);
131   pMin.set(bmin.x()*scale.x(),bmin.y()*scale.y(),bmin.z()*scale.z());
132   pMax.set(bmax.x()*scale.x(),bmax.y()*scale.y(),bmax.z()*scale.z());
133 
134   // Check correctness of the bounding box
135   //
136   if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
137   {
138     std::ostringstream message;
139     message << "Bad bounding box (min >= max) for solid: "
140             << GetName() << " !"
141             << "\npMin = " << pMin
142            << "\npMax = " << pMax;
143     G4Exception("G4ScaledSolid::BoundingLimits()", "GeomMgt0001",
144                 JustWarning, message);
145     DumpInfo();
146   }
147 }
148 
149 //////////////////////////////////////////////////////////////////////////
150 //
151 // Calculate extent under transform and specified limit
152 //
153 G4bool 
154 G4ScaledSolid::CalculateExtent( const EAxis pAxis,
155                                 const G4VoxelLimits& pVoxelLimit,
156                                 const G4AffineTransform& pTransform,
157                                       G4double& pMin,
158                                       G4double& pMax ) const
159 {
160   // Find bounding box of unscaled solid
161   G4ThreeVector bmin,bmax;
162   fPtrSolid->BoundingLimits(bmin,bmax);
163 
164   // Set combined transformation
165   G4Transform3D transform3D =
166     G4Transform3D(pTransform.NetRotation().inverse(),
167                   pTransform.NetTranslation())*GetScaleTransform();
168 
169   // Find extent
170   G4BoundingEnvelope bbox(bmin,bmax);
171   return bbox.CalculateExtent(pAxis,pVoxelLimit,transform3D,pMin,pMax);
172 }
173  
174 /////////////////////////////////////////////////////
175 //
176 // Inside
177 //
178 EInside G4ScaledSolid::Inside(const G4ThreeVector& p) const
179 {
180   return fPtrSolid->Inside(fScale->Transform(p));
181 }
182 
183 //////////////////////////////////////////////////////////////
184 //
185 // SurfaceNormal
186 //
187 G4ThreeVector
188 G4ScaledSolid::SurfaceNormal( const G4ThreeVector& p ) const
189 {
190   // Transform point to unscaled shape frame
191   G4ThreeVector newPoint;
192   fScale->Transform(p, newPoint);
193 
194   // Compute normal in unscaled frame
195   G4ThreeVector newNormal = fPtrSolid->SurfaceNormal(newPoint);
196   G4ThreeVector normal;
197 
198   // Convert normal to scaled frame
199   fScale->InverseTransformNormal(newNormal, normal);
200   return normal/normal.mag();
201 }
202 
203 /////////////////////////////////////////////////////////////
204 //
205 // The same algorithm as in DistanceToIn(p)
206 //
207 G4double
208 G4ScaledSolid::DistanceToIn( const G4ThreeVector& p,
209                              const G4ThreeVector& v ) const
210 {
211   // Transform point and direction to unscaled shape frame
212   G4ThreeVector newPoint;
213   fScale->Transform(p, newPoint);
214 
215   // Direction is un-normalized after scale transformation
216   G4ThreeVector newDirection;
217   fScale->Transform(v, newDirection);
218   newDirection = newDirection/newDirection.mag();
219 
220   // Compute distance in unscaled system
221   G4double dist = fPtrSolid->DistanceToIn(newPoint,newDirection);
222 
223   // Return converted distance to global
224   return fScale->InverseTransformDistance(dist, newDirection);
225 }
226 
227 ////////////////////////////////////////////////////////
228 //
229 // Approximate nearest distance from the point p to the solid from outside
230 //
231 G4double
232 G4ScaledSolid::DistanceToIn( const G4ThreeVector& p ) const
233 {
234   // Transform point to unscaled shape frame
235   G4ThreeVector newPoint;
236   fScale->Transform(p, newPoint);
237 
238   // Compute unscaled safety, then scale it
239   G4double dist = fPtrSolid->DistanceToIn(newPoint);
240   return fScale->InverseTransformDistance(dist);
241 }
242 
243 //////////////////////////////////////////////////////////
244 //
245 // The same algorithm as DistanceToOut(p)
246 //
247 G4double
248 G4ScaledSolid::DistanceToOut( const G4ThreeVector& p,
249                               const G4ThreeVector& v,
250                               const G4bool calcNorm,
251                                     G4bool *validNorm,
252                                     G4ThreeVector *n ) const
253 {
254   // Transform point and direction to unscaled shape frame
255   G4ThreeVector newPoint;
256   fScale->Transform(p, newPoint);
257 
258   // Direction is un-normalized after scale transformation
259   G4ThreeVector newDirection;
260   fScale->Transform(v, newDirection);
261   newDirection = newDirection/newDirection.mag();
262 
263   // Compute distance in unscaled system
264   G4ThreeVector solNorm;
265   G4double dist = fPtrSolid->DistanceToOut(newPoint,newDirection,
266                                            calcNorm,validNorm,&solNorm);
267   if(calcNorm)
268   {
269     G4ThreeVector normal;
270     fScale->TransformNormal(solNorm, normal);
271     *n = normal.unit();
272   }
273 
274   // Return distance converted to global
275   return fScale->InverseTransformDistance(dist, newDirection);
276 }
277 
278 //////////////////////////////////////////////////////////////
279 //
280 //  Approximate nearest distance from the point p to the solid from inside
281 //
282 G4double
283 G4ScaledSolid::DistanceToOut( const G4ThreeVector& p ) const
284 {
285   // Transform point to unscaled shape frame
286   G4ThreeVector newPoint;
287   fScale->Transform(p, newPoint);
288 
289   // Compute unscaled safety, then scale it
290   G4double dist = fPtrSolid->DistanceToOut(newPoint);
291   return fScale->InverseTransformDistance(dist);
292 }
293 
294 //////////////////////////////////////////////////////////////
295 //
296 // ComputeDimensions
297 //
298 void
299 G4ScaledSolid::ComputeDimensions( G4VPVParameterisation*,
300                                   const G4int,
301                                   const G4VPhysicalVolume* ) 
302 {
303   DumpInfo();
304   G4Exception("G4ScaledSolid::ComputeDimensions()",
305               "GeomSolids0001", FatalException,
306               "Method not applicable in this context!");
307 }
308 
309 //////////////////////////////////////////////////////////////////////////
310 //
311 // Returns a point (G4ThreeVector) randomly and uniformly selected
312 // on the solid surface
313 //
314 G4ThreeVector G4ScaledSolid::GetPointOnSurface() const
315 {
316   return fScale->InverseTransform(fPtrSolid->GetPointOnSurface());
317 }
318 //////////////////////////////////////////////////////////////////////////
319 //
320 // Return the number of constituents used for construction of the solid
321 //
322 G4int G4ScaledSolid::GetNumOfConstituents() const
323 {
324   return fPtrSolid->GetNumOfConstituents();
325 }
326 
327 //////////////////////////////////////////////////////////////////////////
328 //
329 // Return true if the solid has only planar faces
330 //
331 G4bool G4ScaledSolid::IsFaceted() const
332 {
333   return fPtrSolid->IsFaceted();  
334 }
335 
336 //////////////////////////////////////////////////////////////////////////
337 //
338 // Return object type name
339 //
340 G4GeometryType G4ScaledSolid::GetEntityType() const
341 {
342   return {"G4ScaledSolid"};
343 }
344 
345 //////////////////////////////////////////////////////////////////////////
346 //
347 // Make a clone of the object
348 //
349 G4VSolid* G4ScaledSolid::Clone() const
350 {
351   return new G4ScaledSolid(*this);
352 }
353 
354 //////////////////////////////////////////////////////////////////////////
355 //
356 // Returning the scaling transformation
357 //
358 G4Scale3D G4ScaledSolid::GetScaleTransform() const
359 {
360   return { fScale->GetScale().x(),
361            fScale->GetScale().y(),
362            fScale->GetScale().z() };
363 }
364 
365 //////////////////////////////////////////////////////////////////////////
366 //
367 // Setting the scaling transformation
368 //
369 void G4ScaledSolid::SetScaleTransform(const G4Scale3D& scale)
370 {
371   delete fScale; 
372   fScale = new G4ScaleTransform(scale);
373   fRebuildPolyhedron = true;
374 }
375 
376 //////////////////////////////////////////////////////////////////////////
377 //
378 // Get volume of the scaled solid
379 //
380 G4double G4ScaledSolid::GetCubicVolume()
381 {
382   if(fCubicVolume < 0.)
383   {
384     fCubicVolume = fPtrSolid->GetCubicVolume() *
385                    fScale->GetScale().x() *
386                    fScale->GetScale().y() *
387                    fScale->GetScale().z();
388   }
389   return fCubicVolume;
390 }
391 
392 //////////////////////////////////////////////////////////////////////////
393 //
394 // Get estimated surface area of the scaled solid
395 //
396 G4double G4ScaledSolid::GetSurfaceArea()
397 {
398   if(fSurfaceArea < 0.)
399   {
400     fSurfaceArea = G4VSolid::GetSurfaceArea();
401   }
402   return fSurfaceArea;
403 }
404 
405 //////////////////////////////////////////////////////////////////////////
406 //
407 // Stream object contents to an output stream
408 //
409 std::ostream& G4ScaledSolid::StreamInfo(std::ostream& os) const
410 {
411   os << "-----------------------------------------------------------\n"
412      << "    *** Dump for Scaled solid - " << GetName() << " ***\n"
413      << "    ===================================================\n"
414      << " Solid type: " << GetEntityType() << "\n"
415      << " Parameters of constituent solid: \n"
416      << "===========================================================\n";
417   fPtrSolid->StreamInfo(os);
418   os << "===========================================================\n"
419      << " Scaling: \n"
420      << "    Scale transformation : \n"
421      << "           " << fScale->GetScale().x() << ", "
422                       << fScale->GetScale().y() << ", "
423                       << fScale->GetScale().z() << "\n"
424      << "===========================================================\n";
425 
426   return os;
427 }
428 
429 //////////////////////////////////////////////////////////////////////////
430 //
431 // DescribeYourselfTo
432 //
433 void 
434 G4ScaledSolid::DescribeYourselfTo ( G4VGraphicsScene& scene ) const
435 {
436   scene.AddSolid (*this);
437 }
438 
439 //////////////////////////////////////////////////////////////////////////
440 //
441 // CreatePolyhedron
442 //
443 G4Polyhedron* 
444 G4ScaledSolid::CreatePolyhedron () const
445 {
446   G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
447   if (polyhedron != nullptr)
448   {
449     polyhedron->Transform(GetScaleTransform());
450   }
451   else
452   {
453     DumpInfo();
454     G4Exception("G4ScaledSolid::CreatePolyhedron()",
455                 "GeomSolids2003", JustWarning,
456                 "No G4Polyhedron for scaled solid");
457   }
458   return polyhedron;
459 }
460 
461 //////////////////////////////////////////////////////////////////////////
462 //
463 // GetPolyhedron
464 //
465 G4Polyhedron* G4ScaledSolid::GetPolyhedron () const
466 {
467   if (fpPolyhedron == nullptr ||
468       fRebuildPolyhedron ||
469       fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
470       fpPolyhedron->GetNumberOfRotationSteps())
471     {
472       fpPolyhedron = CreatePolyhedron();
473       fRebuildPolyhedron = false;
474     }
475   return fpPolyhedron;
476 }
477