Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4PVPlacement Class Reference

#include <G4PVPlacement.hh>

Inheritance diagram for G4PVPlacement:
G4VPhysicalVolume LXeMainVolume LXeWLSFiber LXeWLSSlab

Public Member Functions

 G4PVPlacement (G4RotationMatrix *pRot, const G4ThreeVector &tlate, G4LogicalVolume *pCurrentLogical, const G4String &pName, G4LogicalVolume *pMotherLogical, G4bool pMany, G4int pCopyNo, G4bool pSurfChk=false)
 
 G4PVPlacement (const G4Transform3D &Transform3D, G4LogicalVolume *pCurrentLogical, const G4String &pName, G4LogicalVolume *pMotherLogical, G4bool pMany, G4int pCopyNo, G4bool pSurfChk=false)
 
 G4PVPlacement (G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother, G4bool pMany, G4int pCopyNo, G4bool pSurfChk=false)
 
 G4PVPlacement (const G4Transform3D &Transform3D, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother, G4bool pMany, G4int pCopyNo, G4bool pSurfChk=false)
 
virtual ~G4PVPlacement ()
 
G4int GetCopyNo () const
 
void SetCopyNo (G4int CopyNo)
 
G4bool CheckOverlaps (G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
 
 G4PVPlacement (__void__ &)
 
G4bool IsMany () const
 
G4bool IsReplicated () const
 
G4bool IsParameterised () const
 
G4VPVParameterisationGetParameterisation () const
 
void GetReplicationData (EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
 
G4bool IsRegularStructure () const
 
G4int GetRegularStructureId () const
 
- Public Member Functions inherited from G4VPhysicalVolume
 G4VPhysicalVolume (G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
 
virtual ~G4VPhysicalVolume ()
 
G4bool operator== (const G4VPhysicalVolume &p) const
 
G4RotationMatrixGetObjectRotation () const
 
G4RotationMatrix GetObjectRotationValue () const
 
G4ThreeVector GetObjectTranslation () const
 
const G4RotationMatrixGetFrameRotation () const
 
G4ThreeVector GetFrameTranslation () const
 
const G4ThreeVectorGetTranslation () const
 
const G4RotationMatrixGetRotation () const
 
void SetTranslation (const G4ThreeVector &v)
 
G4RotationMatrixGetRotation ()
 
void SetRotation (G4RotationMatrix *)
 
G4LogicalVolumeGetLogicalVolume () const
 
void SetLogicalVolume (G4LogicalVolume *pLogical)
 
G4LogicalVolumeGetMotherLogical () const
 
void SetMotherLogical (G4LogicalVolume *pMother)
 
const G4StringGetName () const
 
void SetName (const G4String &pName)
 
EVolume VolumeType () const
 
virtual G4int GetMultiplicity () const
 
 G4VPhysicalVolume (__void__ &)
 
G4int GetInstanceID () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicalVolume
static const G4PVManagerGetSubInstanceManager ()
 
- Protected Member Functions inherited from G4VPhysicalVolume
void InitialiseWorker (G4VPhysicalVolume *pMasterObject, G4RotationMatrix *pRot, const G4ThreeVector &tlate)
 
void TerminateWorker (G4VPhysicalVolume *pMasterObject)
 
- Protected Attributes inherited from G4VPhysicalVolume
G4int instanceID
 
- Static Protected Attributes inherited from G4VPhysicalVolume
static G4GEOM_DLL G4PVManager subInstanceManager
 

Detailed Description

Definition at line 51 of file G4PVPlacement.hh.

Constructor & Destructor Documentation

G4PVPlacement::G4PVPlacement ( G4RotationMatrix pRot,
const G4ThreeVector tlate,
G4LogicalVolume pCurrentLogical,
const G4String pName,
G4LogicalVolume pMotherLogical,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 100 of file G4PVPlacement.cc.

References G4LogicalVolume::AddDaughter(), CheckOverlaps(), FatalException, G4Exception(), and G4VPhysicalVolume::SetMotherLogical().

Referenced by LXeMainVolume::LXeMainVolume(), and LXeWLSFiber::LXeWLSFiber().

108  : G4VPhysicalVolume(pRot,tlate,pName,pCurrentLogical,0),
109  fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
110 {
111  if (pCurrentLogical == pMotherLogical)
112  {
113  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
114  FatalException, "Cannot place a volume inside itself!");
115  }
116  SetMotherLogical(pMotherLogical);
117  if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
118  if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
119 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
G4PVPlacement::G4PVPlacement ( const G4Transform3D Transform3D,
G4LogicalVolume pCurrentLogical,
const G4String pName,
G4LogicalVolume pMotherLogical,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 125 of file G4PVPlacement.cc.

References G4LogicalVolume::AddDaughter(), CheckOverlaps(), FatalException, G4Exception(), G4VPhysicalVolume::GetRotation(), HepGeom::Transform3D::getRotation(), CLHEP::HepRotation::inverse(), G4VPhysicalVolume::SetMotherLogical(), and G4VPhysicalVolume::SetRotation().

132  : G4VPhysicalVolume(0,Transform3D.getTranslation(),pName,pCurrentLogical,0),
133  fmany(pMany), fcopyNo(pCopyNo)
134 {
135  if (pCurrentLogical == pMotherLogical)
136  {
137  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
138  FatalException, "Cannot place a volume inside itself!");
139  }
140  SetRotation( NewPtrRotMatrix(Transform3D.getRotation().inverse()) );
141  fallocatedRotM = (GetRotation() != 0);
142  SetMotherLogical(pMotherLogical);
143  if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
144  if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
145 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
HepRotation inverse() const
void SetRotation(G4RotationMatrix *)
CLHEP::HepRotation getRotation() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
const G4RotationMatrix * GetRotation() const
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
CLHEP::Hep3Vector getTranslation() const
G4PVPlacement::G4PVPlacement ( G4RotationMatrix pRot,
const G4ThreeVector tlate,
const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMother,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 43 of file G4PVPlacement.cc.

References G4LogicalVolume::AddDaughter(), CheckOverlaps(), FatalException, G4Exception(), G4VPhysicalVolume::GetLogicalVolume(), and G4VPhysicalVolume::SetMotherLogical().

51  : G4VPhysicalVolume(pRot,tlate,pName,pLogical,pMother),
52  fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
53 {
54  if (pMother)
55  {
56  G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
57  if (pLogical == motherLogical)
58  {
59  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
60  FatalException, "Cannot place a volume inside itself!");
61  }
62  SetMotherLogical(motherLogical);
63  motherLogical->AddDaughter(this);
64  if (pSurfChk) { CheckOverlaps(); }
65  }
66 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
G4LogicalVolume * GetLogicalVolume() const
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
G4PVPlacement::G4PVPlacement ( const G4Transform3D Transform3D,
const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMother,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 71 of file G4PVPlacement.cc.

References G4LogicalVolume::AddDaughter(), CheckOverlaps(), FatalException, G4Exception(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetRotation(), and G4VPhysicalVolume::SetMotherLogical().

78  : G4VPhysicalVolume(NewPtrRotMatrix(Transform3D.getRotation().inverse()),
79  Transform3D.getTranslation(),pName,pLogical,pMother),
80  fmany(pMany), fcopyNo(pCopyNo)
81 {
82  fallocatedRotM = (GetRotation() != 0);
83  if (pMother)
84  {
85  G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
86  if (pLogical == motherLogical)
87  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
88  FatalException, "Cannot place a volume inside itself!");
89  SetMotherLogical(motherLogical);
90  motherLogical->AddDaughter(this);
91  if (pSurfChk) { CheckOverlaps(); }
92  }
93 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)
HepRotation inverse() const
CLHEP::HepRotation getRotation() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
G4LogicalVolume * GetLogicalVolume() const
const G4RotationMatrix * GetRotation() const
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
CLHEP::Hep3Vector getTranslation() const
G4PVPlacement::~G4PVPlacement ( )
virtual

Definition at line 159 of file G4PVPlacement.cc.

References G4VPhysicalVolume::GetRotation().

160 {
161  if( fallocatedRotM ){ delete this->GetRotation() ; }
162 }
const G4RotationMatrix * GetRotation() const
G4PVPlacement::G4PVPlacement ( __void__ &  a)

Definition at line 151 of file G4PVPlacement.cc.

152  : G4VPhysicalVolume(a), fmany(false), fallocatedRotM(0), fcopyNo(0)
153 {
154 }
G4VPhysicalVolume(G4RotationMatrix *pRot, const G4ThreeVector &tlate, const G4String &pName, G4LogicalVolume *pLogical, G4VPhysicalVolume *pMother)

Member Function Documentation

G4bool G4PVPlacement::CheckOverlaps ( G4int  res = 1000,
G4double  tol = 0.,
G4bool  verbose = true,
G4int  maxErr = 1 
)
virtual

Reimplemented from G4VPhysicalVolume.

Definition at line 244 of file G4PVPlacement.cc.

References G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), G4BestUnit, G4cout, G4endl, G4Exception(), G4LogicalVolume::GetDaughter(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetMotherLogical(), G4VPhysicalVolume::GetName(), G4LogicalVolume::GetName(), G4LogicalVolume::GetNoDaughters(), G4VSolid::GetPointOnSurface(), G4VPhysicalVolume::GetRotation(), G4LogicalVolume::GetSolid(), G4VPhysicalVolume::GetTranslation(), G4VSolid::Inside(), G4AffineTransform::Inverse(), JustWarning, kInside, kOutside, n, and G4AffineTransform::TransformPoint().

Referenced by export_G4PVPlacement(), and G4PVPlacement().

246 {
247  if (res<=0) { return false; }
248 
249  G4VSolid* solid = GetLogicalVolume()->GetSolid();
250  G4LogicalVolume* motherLog = GetMotherLogical();
251  if (!motherLog) { return false; }
252 
253  G4int trials = 0;
254  G4bool retval = false;
255 
256  G4VSolid* motherSolid = motherLog->GetSolid();
257 
258  if (verbose)
259  {
260  G4cout << "Checking overlaps for volume " << GetName() << " ... ";
261  }
262 
263  // Create the transformation from daughter to mother
264  //
266 
267  for (G4int n=0; n<res; n++)
268  {
269  // Generate a random point on the solid's surface
270  //
271  G4ThreeVector point = solid->GetPointOnSurface();
272 
273  // Transform the generated point to the mother's coordinate system
274  //
275  G4ThreeVector mp = Tm.TransformPoint(point);
276 
277  // Checking overlaps with the mother volume
278  //
279  if (motherSolid->Inside(mp)==kOutside)
280  {
281  G4double distin = motherSolid->DistanceToIn(mp);
282  if (distin > tol)
283  {
284  trials++; retval = true;
285  std::ostringstream message;
286  message << "Overlap with mother volume !" << G4endl
287  << " Overlap is detected for volume "
288  << GetName() << G4endl
289  << " with its mother volume "
290  << motherLog->GetName() << G4endl
291  << " at mother local point " << mp << ", "
292  << "overlapping by at least: "
293  << G4BestUnit(distin, "Length");
294  if (trials>=maxErr)
295  {
296  message << G4endl
297  << "NOTE: Reached maximum fixed number -" << maxErr
298  << "- of overlaps reports for this volume !";
299  }
300  G4Exception("G4PVPlacement::CheckOverlaps()",
301  "GeomVol1002", JustWarning, message);
302  if (trials>=maxErr) { return true; }
303  }
304  }
305 
306  // Checking overlaps with each 'sister' volume
307  //
308  for (G4int i=0; i<motherLog->GetNoDaughters(); i++)
309  {
310  G4VPhysicalVolume* daughter = motherLog->GetDaughter(i);
311 
312  if (daughter == this) { continue; }
313 
314  // Create the transformation for daughter volume and transform point
315  //
316  G4AffineTransform Td( daughter->GetRotation(),
317  daughter->GetTranslation() );
318  G4ThreeVector md = Td.Inverse().TransformPoint(mp);
319 
320  G4VSolid* daughterSolid = daughter->GetLogicalVolume()->GetSolid();
321  if (daughterSolid->Inside(md)==kInside)
322  {
323  G4double distout = daughterSolid->DistanceToOut(md);
324  if (distout > tol)
325  {
326  trials++; retval = true;
327  std::ostringstream message;
328  message << "Overlap with volume already placed !" << G4endl
329  << " Overlap is detected for volume "
330  << GetName() << G4endl
331  << " with " << daughter->GetName() << " volume's"
332  << G4endl
333  << " local point " << md << ", "
334  << "overlapping by at least: "
335  << G4BestUnit(distout,"Length");
336  if (trials>=maxErr)
337  {
338  message << G4endl
339  << "NOTE: Reached maximum fixed number -" << maxErr
340  << "- of overlaps reports for this volume !";
341  }
342  G4Exception("G4PVPlacement::CheckOverlaps()",
343  "GeomVol1002", JustWarning, message);
344  if (trials>=maxErr) { return true; }
345  }
346  }
347 
348  // Now checking that 'sister' volume is not totally included and
349  // overlapping. Do it only once, for the first point generated
350  //
351  if (n==0)
352  {
353  // Generate a single point on the surface of the 'sister' volume
354  // and verify that the point is NOT inside the current volume
355 
356  G4ThreeVector dPoint = daughterSolid->GetPointOnSurface();
357 
358  // Transform the generated point to the mother's coordinate system
359  // and finally to current volume's coordinate system
360  //
361  G4ThreeVector mp2 = Td.TransformPoint(dPoint);
362  G4ThreeVector msi = Tm.Inverse().TransformPoint(mp2);
363 
364  if (solid->Inside(msi)==kInside)
365  {
366  trials++; retval = true;
367  std::ostringstream message;
368  message << "Overlap with volume already placed !" << G4endl
369  << " Overlap is detected for volume "
370  << GetName() << G4endl
371  << " apparently fully encapsulating volume "
372  << daughter->GetName() << G4endl
373  << " at the same level !";
374  G4Exception("G4PVPlacement::CheckOverlaps()",
375  "GeomVol1002", JustWarning, message);
376  if (trials>=maxErr) { return true; }
377  }
378  }
379  }
380  }
381 
382  if (verbose)
383  {
384  G4cout << "OK! " << G4endl;
385  }
386 
387  return retval;
388 }
const G4ThreeVector & GetTranslation() const
G4String GetName() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
const G4int n
G4LogicalVolume * GetMotherLogical() const
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * GetLogicalVolume() const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
const G4RotationMatrix * GetRotation() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4VSolid * GetSolid() const
G4int G4PVPlacement::GetCopyNo ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 175 of file G4PVPlacement.cc.

176 {
177  return fcopyNo;
178 }
G4VPVParameterisation * G4PVPlacement::GetParameterisation ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 207 of file G4PVPlacement.cc.

208 {
209  return 0;
210 }
G4int G4PVPlacement::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 236 of file G4PVPlacement.cc.

237 {
238  return 0;
239 }
void G4PVPlacement::GetReplicationData ( EAxis axis,
G4int nReplicas,
G4double width,
G4double offset,
G4bool consuming 
) const
virtual

Implements G4VPhysicalVolume.

Definition at line 216 of file G4PVPlacement.cc.

217 {
218  // No-operations
219 }
G4bool G4PVPlacement::IsMany ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 167 of file G4PVPlacement.cc.

168 {
169  return fmany;
170 }
G4bool G4PVPlacement::IsParameterised ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 199 of file G4PVPlacement.cc.

200 {
201  return false;
202 }
G4bool G4PVPlacement::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 226 of file G4PVPlacement.cc.

227 {
228  return false;
229 }
G4bool G4PVPlacement::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 191 of file G4PVPlacement.cc.

192 {
193  return false;
194 }
void G4PVPlacement::SetCopyNo ( G4int  CopyNo)
virtual

Implements G4VPhysicalVolume.

Definition at line 183 of file G4PVPlacement.cc.

184 {
185  fcopyNo= newCopyNo;
186 }

The documentation for this class was generated from the following files: