Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PVPlacement Class Reference

#include <G4PVPlacement.hh>

Inheritance diagram for G4PVPlacement:
Collaboration diagram for G4PVPlacement:

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.

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)

Here is the call graph for this function:

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.

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)
const G4RotationMatrix * GetRotation() const
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)
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
CLHEP::Hep3Vector getTranslation() const

Here is the call graph for this function:

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.

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)

Here is the call graph for this function:

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.

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)
const G4RotationMatrix * GetRotation() const
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
void SetMotherLogical(G4LogicalVolume *pMother)
void AddDaughter(G4VPhysicalVolume *p)
CLHEP::Hep3Vector getTranslation() const

Here is the call graph for this function:

G4PVPlacement::~G4PVPlacement ( )
virtual

Definition at line 159 of file G4PVPlacement.cc.

160 {
161  if( fallocatedRotM ){ delete this->GetRotation() ; }
162 }
const G4RotationMatrix * GetRotation() const

Here is the call graph for this function:

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 236 of file G4PVPlacement.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4int G4PVPlacement::GetCopyNo ( ) const
inlinevirtual

Implements G4VPhysicalVolume.

Definition at line 125 of file G4PVPlacement.hh.

125 { return fcopyNo; }
G4VPVParameterisation * G4PVPlacement::GetParameterisation ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 199 of file G4PVPlacement.cc.

200 {
201  return 0;
202 }
G4int G4PVPlacement::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 228 of file G4PVPlacement.cc.

229 {
230  return 0;
231 }
void G4PVPlacement::GetReplicationData ( EAxis axis,
G4int nReplicas,
G4double width,
G4double offset,
G4bool consuming 
) const
virtual

Implements G4VPhysicalVolume.

Definition at line 208 of file G4PVPlacement.cc.

209 {
210  // No-operations
211 }
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 191 of file G4PVPlacement.cc.

192 {
193  return false;
194 }
G4bool G4PVPlacement::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 218 of file G4PVPlacement.cc.

219 {
220  return false;
221 }
G4bool G4PVPlacement::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 183 of file G4PVPlacement.cc.

184 {
185  return false;
186 }
void G4PVPlacement::SetCopyNo ( G4int  CopyNo)
virtual

Implements G4VPhysicalVolume.

Definition at line 175 of file G4PVPlacement.cc.

176 {
177  fcopyNo= newCopyNo;
178 }

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