Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4ScaledSolid.cc
이 파일의 문서화 페이지로 가기
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 //
27 // $Id:$
28 //
29 // Implementation for G4ScaledSolid class
30 //
31 // History:
32 //
33 // 27.10.15 G.Cosmo: created, based on implementation also provided in Root
34 //
35 // --------------------------------------------------------------------
36 
37 #include "G4ScaledSolid.hh"
38 #include "G4BoundingEnvelope.hh"
39 
40 #include "G4VPVParameterisation.hh"
41 
42 #include "G4ScaleTransform.hh"
43 
44 #include "G4VGraphicsScene.hh"
45 #include "G4Polyhedron.hh"
46 
48 //
49 // Constructor
50 //
52  G4VSolid* pSolid ,
53  const G4Scale3D& pScale )
54  : G4VSolid(pName), fPtrSolid(pSolid),
55  fRebuildPolyhedron(false), fpPolyhedron(0)
56 {
57  fScale = new G4ScaleTransform(pScale);
58 }
59 
61 //
62 // Fake default constructor - sets only member data and allocates memory
63 // for usage restricted to object persistency.
64 //
66  : G4VSolid(a), fPtrSolid(0), fScale(0),
67  fRebuildPolyhedron(false), fpPolyhedron(0)
68 {
69 }
70 
72 //
73 // Destructor
74 //
76 {
77  delete fpPolyhedron; fpPolyhedron= 0;
78  delete fScale; fScale= 0;
79 }
80 
82 //
83 // Copy constructor
84 //
86  : G4VSolid (rhs), fPtrSolid(rhs.fPtrSolid),
87  fRebuildPolyhedron(false), fpPolyhedron(0)
88 {
89  fScale = new G4ScaleTransform(*(rhs.fScale));
90 }
91 
93 //
94 // Assignment operator
95 //
97 {
98  // Check assignment to self
99  //
100  if (this == &rhs) { return *this; }
101 
102  // Copy base class data
103  //
104  G4VSolid::operator=(rhs);
105 
106  // Copy data
107  //
108  fPtrSolid = rhs.fPtrSolid;
109  delete fScale;
110  fScale = new G4ScaleTransform(*(rhs.fScale));
111  fRebuildPolyhedron = false;
112  delete fpPolyhedron; fpPolyhedron= 0;
113 
114  return *this;
115 }
116 
118 //
119 // Return original solid not scaled
120 //
122 {
123  return fPtrSolid;
124 }
125 
127 //
128 // Get bounding box
129 
131  G4ThreeVector& pMax) const
132 {
133  G4ThreeVector bmin,bmax;
135 
136  fPtrSolid->BoundingLimits(bmin,bmax);
137  pMin.set(bmin.x()*scale.x(),bmin.y()*scale.y(),bmin.z()*scale.z());
138  pMax.set(bmax.x()*scale.x(),bmax.y()*scale.y(),bmax.z()*scale.z());
139 
140  // Check correctness of the bounding box
141  //
142  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
143  {
144  std::ostringstream message;
145  message << "Bad bounding box (min >= max) for solid: "
146  << GetName() << " !"
147  << "\npMin = " << pMin
148  << "\npMax = " << pMax;
149  G4Exception("G4ScaledSolid::BoundingLimits()", "GeomMgt0001",
150  JustWarning, message);
151  DumpInfo();
152  }
153 }
154 
156 //
157 // Calculate extent under transform and specified limit
158 //
159 G4bool
161  const G4VoxelLimits& pVoxelLimit,
162  const G4AffineTransform& pTransform,
163  G4double& pMin,
164  G4double& pMax ) const
165 {
166  // Find bounding box of unscaled solid
167  G4ThreeVector bmin,bmax;
168  fPtrSolid->BoundingLimits(bmin,bmax);
169 
170  // Set combined transformation
171  G4Transform3D transform3D =
172  G4Transform3D(pTransform.NetRotation().inverse(),
173  pTransform.NetTranslation())*GetScaleTransform();
174 
175  // Find extent
176  G4BoundingEnvelope bbox(bmin,bmax);
177  return bbox.CalculateExtent(pAxis,pVoxelLimit,transform3D,pMin,pMax);
178 }
179 
181 //
182 // Inside
183 //
185 {
186  return fPtrSolid->Inside(fScale->Transform(p));
187 }
188 
190 //
191 // SurfaceNormal
192 //
195 {
196  // Transform point to unscaled shape frame
197  G4ThreeVector newPoint;
198  fScale->Transform(p, newPoint);
199 
200  // Compute normal in unscaled frame
201  G4ThreeVector newNormal = fPtrSolid->SurfaceNormal(newPoint);
203 
204  // Convert normal to scaled frame
205  fScale->InverseTransformNormal(newNormal, normal);
206  return normal/normal.mag();
207 }
208 
210 //
211 // The same algorithm as in DistanceToIn(p)
212 //
213 G4double
215  const G4ThreeVector& v ) const
216 {
217  // Transform point and direction to unscaled shape frame
218  G4ThreeVector newPoint;
219  fScale->Transform(p, newPoint);
220 
221  // Direction is un-normalized after scale transformation
222  G4ThreeVector newDirection;
223  fScale->Transform(v, newDirection);
224  newDirection = newDirection/newDirection.mag();
225 
226  // Compute distance in unscaled system
227  G4double dist = fPtrSolid->DistanceToIn(newPoint,newDirection);
228 
229  // Return converted distance to global
230  return fScale->InverseTransformDistance(dist, newDirection);
231 }
232 
234 //
235 // Approximate nearest distance from the point p to the solid from outside
236 //
237 G4double
239 {
240  // Transform point to unscaled shape frame
241  G4ThreeVector newPoint;
242  fScale->Transform(p, newPoint);
243 
244  // Compute unscaled safety, then scale it.
245  G4double dist = fPtrSolid->DistanceToIn(newPoint);
246  return fScale->InverseTransformDistance(dist);
247 }
248 
250 //
251 // The same algorithm as DistanceToOut(p)
252 //
253 G4double
255  const G4ThreeVector& v,
256  const G4bool calcNorm,
257  G4bool *validNorm,
258  G4ThreeVector *n ) const
259 {
260  // Transform point and direction to unscaled shape frame
261  G4ThreeVector newPoint;
262  fScale->Transform(p, newPoint);
263 
264  // Direction is un-normalized after scale transformation
265  G4ThreeVector newDirection;
266  fScale->Transform(v, newDirection);
267  newDirection = newDirection/newDirection.mag();
268 
269  // Compute distance in unscaled system
270  G4ThreeVector solNorm;
271  G4double dist = fPtrSolid->DistanceToOut(newPoint,newDirection,
272  calcNorm,validNorm,&solNorm);
273  if(calcNorm)
274  {
276  fScale->TransformNormal(solNorm, normal);
277  *n = normal/normal.mag();
278  }
279 
280  // Return distance converted to global
281  return fScale->InverseTransformDistance(dist, newDirection);
282 }
283 
285 //
286 // Approximate nearest distance from the point p to the solid from inside
287 //
288 G4double
290 {
291  // Transform point to unscaled shape frame
292  G4ThreeVector newPoint;
293  fScale->Transform(p, newPoint);
294 
295  // Compute unscaled safety, then scale it.
296  G4double dist = fPtrSolid->DistanceToOut(newPoint);
297  return fScale->InverseTransformDistance(dist);
298 }
299 
301 //
302 // ComputeDimensions
303 //
304 void
306  const G4int,
307  const G4VPhysicalVolume* )
308 {
309  DumpInfo();
310  G4Exception("G4ScaledSolid::ComputeDimensions()",
311  "GeomSolids0001", FatalException,
312  "Method not applicable in this context!");
313 }
314 
316 //
317 // Returns a point (G4ThreeVector) randomly and uniformly selected
318 // on the solid surface
319 //
321 {
323 }
324 
326 //
327 // Return object type name
328 //
330 {
331  return G4String("G4ScaledSolid");
332 }
333 
335 //
336 // Make a clone of the object
337 //
339 {
340  return new G4ScaledSolid(*this);
341 }
342 
344 //
345 // Returning the scaling transformation
346 //
348 {
349  return G4Scale3D(fScale->GetScale().x(),
350  fScale->GetScale().y(),
351  fScale->GetScale().z());
352 }
353 
355 //
356 // Setting the scaling transformation
357 //
359 {
360  if (fScale) { delete fScale; }
361  fScale = new G4ScaleTransform(scale);
362  fRebuildPolyhedron = true;
363 }
364 
366 //
367 // Stream object contents to an output stream
368 //
369 std::ostream& G4ScaledSolid::StreamInfo(std::ostream& os) const
370 {
371  os << "-----------------------------------------------------------\n"
372  << " *** Dump for Scaled solid - " << GetName() << " ***\n"
373  << " ===================================================\n"
374  << " Solid type: " << GetEntityType() << "\n"
375  << " Parameters of constituent solid: \n"
376  << "===========================================================\n";
377  fPtrSolid->StreamInfo(os);
378  os << "===========================================================\n"
379  << " Scaling: \n"
380  << " Scale transformation : \n"
381  << " " << fScale->GetScale().x() << ", "
382  << fScale->GetScale().y() << ", "
383  << fScale->GetScale().z() << "\n"
384  << "===========================================================\n";
385 
386  return os;
387 }
388 
390 //
391 // DescribeYourselfTo
392 //
393 void
395 {
396  scene.AddSolid (*this);
397 }
398 
400 //
401 // CreatePolyhedron
402 //
403 G4Polyhedron*
405 {
406  G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
407  if (polyhedron)
408  {
409  polyhedron->Transform(GetScaleTransform());
410  }
411  else
412  {
413  DumpInfo();
414  G4Exception("G4ScaledSolid::CreatePolyhedron()",
415  "GeomSolids2003", JustWarning,
416  "No G4Polyhedron for scaled solid");
417  }
418  return polyhedron;
419 }
420 
422 //
423 // GetPolyhedron
424 //
426 {
427  if (!fpPolyhedron ||
431  {
433  fRebuildPolyhedron = false;
434  }
435  return fpPolyhedron;
436 }
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
void set(double x, double y, double z)
G4bool fRebuildPolyhedron
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
const G4ThreeVector & GetScale() const
EInside Inside(const G4ThreeVector &p) const
virtual void AddSolid(const G4Box &)=0
HepGeom::Transform3D G4Transform3D
void TransformNormal(const G4ThreeVector &global, G4ThreeVector &local) const
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:625
const char * p
Definition: xmltok.h:285
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
G4Scale3D GetScaleTransform() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
double z() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4ThreeVector GetPointOnSurface() const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
G4VSolid * Clone() const
void SetScaleTransform(const G4Scale3D &scale)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
virtual EInside Inside(const G4ThreeVector &p) const =0
std::ostream & StreamInfo(std::ostream &os) const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Double_t scale
HepPolyhedron & Transform(const G4Transform3D &t)
void InverseTransform(const G4ThreeVector &local, G4ThreeVector &global) const
G4ScaleTransform * fScale
virtual std::ostream & StreamInfo(std::ostream &os) const =0
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
G4VSolid * fPtrSolid
static G4int GetNumberOfRotationSteps()
HepGeom::Scale3D G4Scale3D
G4ThreeVector NetTranslation() const
G4GeometryType GetEntityType() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
double mag() const
G4ScaledSolid(const G4String &pName, G4VSolid *pSolid, const G4Scale3D &pScale)
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
virtual ~G4ScaledSolid()
EAxis
Definition: geomdefs.hh:54
G4String GetName() const
void DumpInfo() const
G4Polyhedron * fpPolyhedron
double x() const
G4VSolid * GetUnscaledSolid() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
virtual G4Polyhedron * CreatePolyhedron() const
Definition: G4VSolid.cc:660
G4ScaledSolid & operator=(const G4ScaledSolid &rhs)
G4double InverseTransformDistance(G4double dist, const G4ThreeVector &dir) const
Char_t n[5]
void InverseTransformNormal(const G4ThreeVector &local, G4ThreeVector &global) const
double y() const
void Transform(const G4ThreeVector &global, G4ThreeVector &local) const
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
G4Polyhedron * CreatePolyhedron() const
HepRotation inverse() const
G4Polyhedron * GetPolyhedron() const
G4RotationMatrix NetRotation() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const