Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Static Private Member Functions | Private Attributes
G4PVPlacement Class Reference

#include <G4PVPlacement.hh>

Inheritance diagram for G4PVPlacement:
G4VPhysicalVolume

Public Member Functions

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

Static Public Member Functions

static void Clean ()
 
static const G4PVManagerGetSubInstanceManager ()
 

Protected Member Functions

void InitialiseWorker (G4VPhysicalVolume *pMasterObject, G4RotationMatrix *pRot, const G4ThreeVector &tlate)
 
void TerminateWorker (G4VPhysicalVolume *pMasterObject)
 

Protected Attributes

G4int instanceID
 

Static Protected Attributes

static G4GEOM_DLL G4PVManager subInstanceManager
 

Static Private Member Functions

static G4RotationMatrixNewPtrRotMatrix (const G4RotationMatrix &RotMat)
 

Private Attributes

G4bool fallocatedRotM = false
 
G4int fcopyNo = 0
 
G4LogicalVolumeflmother = nullptr
 
G4LogicalVolumeflogical = nullptr
 
G4bool fmany = false
 
G4String fname
 
G4PVDatapvdata = nullptr
 

Detailed Description

Definition at line 43 of file G4PVPlacement.hh.

Constructor & Destructor Documentation

◆ G4PVPlacement() [1/6]

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

Definition at line 96 of file G4PVPlacement.cc.

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

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

◆ G4PVPlacement() [2/6]

G4PVPlacement::G4PVPlacement ( const G4Transform3D Transform3D,
G4LogicalVolume pCurrentLogical,
const G4String pName,
G4LogicalVolume pMotherLogical,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 121 of file G4PVPlacement.cc.

128 : G4VPhysicalVolume(nullptr, Transform3D.getTranslation(),
129 pName, pCurrentLogical, nullptr),
130 fmany(pMany), fcopyNo(pCopyNo)
131{
132 if (pCurrentLogical == pMotherLogical)
133 {
134 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
135 FatalException, "Cannot place a volume inside itself!");
136 }
138 fallocatedRotM = (GetRotation() != nullptr);
139 SetMotherLogical(pMotherLogical);
140 if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
141 if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
142}
HepRotation inverse() const
static G4RotationMatrix * NewPtrRotMatrix(const G4RotationMatrix &RotMat)
G4bool fallocatedRotM
const G4RotationMatrix * GetRotation() const
void SetRotation(G4RotationMatrix *)
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const

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

◆ G4PVPlacement() [3/6]

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

Definition at line 39 of file G4PVPlacement.cc.

47 : G4VPhysicalVolume(pRot, tlate, pName, pLogical, pMother),
48 fmany(pMany), fcopyNo(pCopyNo)
49{
50 if (pMother)
51 {
52 G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
53 if (pLogical == motherLogical)
54 {
55 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
56 FatalException, "Cannot place a volume inside itself!");
57 }
58 SetMotherLogical(motherLogical);
59 motherLogical->AddDaughter(this);
60 if (pSurfChk) { CheckOverlaps(); }
61 }
62}
G4LogicalVolume * GetLogicalVolume() const

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

◆ G4PVPlacement() [4/6]

G4PVPlacement::G4PVPlacement ( const G4Transform3D Transform3D,
const G4String pName,
G4LogicalVolume pLogical,
G4VPhysicalVolume pMother,
G4bool  pMany,
G4int  pCopyNo,
G4bool  pSurfChk = false 
)

Definition at line 67 of file G4PVPlacement.cc.

75 Transform3D.getTranslation(), pName, pLogical, pMother),
76 fmany(pMany), fcopyNo(pCopyNo)
77{
78 fallocatedRotM = (GetRotation() != 0);
79 if (pMother)
80 {
81 G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
82 if (pLogical == motherLogical)
83 G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
84 FatalException, "Cannot place a volume inside itself!");
85 SetMotherLogical(motherLogical);
86 motherLogical->AddDaughter(this);
87 if (pSurfChk) { CheckOverlaps(); }
88 }
89}

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

◆ ~G4PVPlacement()

G4PVPlacement::~G4PVPlacement ( )
virtual

Definition at line 156 of file G4PVPlacement.cc.

157{
158 if( fallocatedRotM ){ delete this->GetRotation() ; }
159}

References fallocatedRotM, and G4VPhysicalVolume::GetRotation().

◆ G4PVPlacement() [5/6]

G4PVPlacement::G4PVPlacement ( __void__ &  a)

Definition at line 148 of file G4PVPlacement.cc.

150{
151}

◆ G4PVPlacement() [6/6]

G4PVPlacement::G4PVPlacement ( const G4PVPlacement )
delete

Member Function Documentation

◆ CheckOverlaps()

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

Reimplemented from G4VPhysicalVolume.

Definition at line 243 of file G4PVPlacement.cc.

245{
246 if (res <= 0) { return false; }
247
248 G4VSolid* solid = GetLogicalVolume()->GetSolid();
249 G4LogicalVolume* motherLog = GetMotherLogical();
250 if (motherLog == nullptr) { return false; }
251
252 G4int trials = 0;
253 G4bool retval = false;
254
255 if (verbose)
256 {
257 G4cout << "Checking overlaps for volume "
258 << GetName() << ':' << GetCopyNo()
259 << " (" << solid->GetEntityType() << ") ... ";
260 }
261
262 // Check that random points are gererated correctly
263 //
264 G4ThreeVector ptmp = solid->GetPointOnSurface();
265 if (solid->Inside(ptmp) != kSurface)
266 {
267 G4String position[3] = { "outside", "surface", "inside" };
268 std::ostringstream message;
269 message << "Sample point is not on the surface !" << G4endl
270 << " The issue is detected for volume "
271 << GetName() << ':' << GetCopyNo()
272 << " (" << solid->GetEntityType() << ")" << G4endl
273 << " generated point " << ptmp
274 << " is " << position[solid->Inside(ptmp)];
275 G4Exception("G4PVPlacement::CheckOverlaps()",
276 "GeomVol1002", JustWarning, message);
277 return false;
278 }
279
280 // Generate random points on the surface of the solid,
281 // transform them into the mother volume coordinate system
282 // and find the bonding box
283 //
284 std::vector<G4ThreeVector> points(res);
285 G4double xmin = kInfinity, ymin = kInfinity, zmin = kInfinity;
286 G4double xmax = -kInfinity, ymax = -kInfinity, zmax = -kInfinity;
288 for (G4int i = 0; i < res; ++i)
289 {
290 points[i] = Tm.TransformPoint(solid->GetPointOnSurface());
291 xmin = std::min(xmin, points[i].x());
292 ymin = std::min(ymin, points[i].y());
293 zmin = std::min(zmin, points[i].z());
294 xmax = std::max(xmax, points[i].x());
295 ymax = std::max(ymax, points[i].y());
296 zmax = std::max(zmax, points[i].z());
297 }
298 G4ThreeVector scenter(0.5*(xmax+xmin), 0.5*(ymax+ymin), 0.5*(zmax+zmin));
299 G4double sradius = 0.5*G4ThreeVector(xmax-xmin, ymax-ymin, zmax-zmin).mag();
300
301 // Check overlap with the mother volume
302 //
303 G4int overlapCount = 0;
304 G4double overlapSize = -kInfinity;
305 G4ThreeVector overlapPoint;
306 G4VSolid* motherSolid = motherLog->GetSolid();
307 for (G4int i = 0; i < res; ++i)
308 {
309 G4ThreeVector mp = points[i];
310 if (motherSolid->Inside(mp) != kOutside) continue;
311 G4double distin = motherSolid->DistanceToIn(mp);
312 if (distin < tol) continue; // too small overlap
313 ++overlapCount;
314 if (distin <= overlapSize) continue;
315 overlapSize = distin;
316 overlapPoint = mp;
317 }
318
319 // Print information on overlap
320 //
321 if (overlapCount > 0)
322 {
323 ++trials;
324 retval = true;
325 std::ostringstream message;
326 message << "Overlap with mother volume !" << G4endl
327 << " Overlap is detected for volume "
328 << GetName() << ':' << GetCopyNo()
329 << " (" << solid->GetEntityType() << ")"
330 << " with its mother volume " << motherLog->GetName()
331 << " (" << motherSolid->GetEntityType() << ")" << G4endl
332 << " protrusion at mother local point " << overlapPoint
333 << " by " << G4BestUnit(overlapSize, "Length")
334 << " (max of " << overlapCount << " cases)";
335 if (trials >= maxErr)
336 {
337 message << G4endl
338 << "NOTE: Reached maximum fixed number -" << maxErr
339 << "- of overlaps reports for this volume !";
340 }
341 G4Exception("G4PVPlacement::CheckOverlaps()",
342 "GeomVol1002", JustWarning, message);
343 if (trials >= maxErr) { return true; }
344 }
345
346 // Checking overlaps with each 'sister' volumes
347 //
348 G4VSolid* previous = nullptr;
349 G4ThreeVector pmin_local(0.,0.,0.), pmax_local(0.,0.,0.);
350
351 for (size_t k = 0; k < motherLog->GetNoDaughters(); ++k)
352 {
353 G4VPhysicalVolume* daughter = motherLog->GetDaughter(k);
354 if (daughter == this) continue;
355 G4bool check_encapsulation = true;
356
357 G4AffineTransform Td(daughter->GetRotation(), daughter->GetTranslation());
358 G4VSolid* daughterSolid = daughter->GetLogicalVolume()->GetSolid();
359 if (previous != daughterSolid)
360 {
361 daughterSolid->BoundingLimits(pmin_local, pmax_local);
362 previous = daughterSolid;
363 }
364 overlapCount = 0;
365 overlapSize = -kInfinity;
366 if (!Td.IsRotated()) { // no rotation, only translation
367 G4ThreeVector offset = Td.NetTranslation();
368 G4ThreeVector pmin(pmin_local + offset);
369 G4ThreeVector pmax(pmax_local + offset);
370 if (pmin.x() >= xmax) continue;
371 if (pmin.y() >= ymax) continue;
372 if (pmin.z() >= zmax) continue;
373 if (pmax.x() <= xmin) continue;
374 if (pmax.y() <= ymin) continue;
375 if (pmax.z() <= zmin) continue;
376 for (G4int i = 0; i < res; ++i)
377 {
378 G4ThreeVector p = points[i];
379 if (p.x() <= pmin.x()) continue;
380 if (p.x() >= pmax.x()) continue;
381 if (p.y() <= pmin.y()) continue;
382 if (p.y() >= pmax.y()) continue;
383 if (p.z() <= pmin.z()) continue;
384 if (p.z() >= pmax.z()) continue;
385 G4ThreeVector md = p - offset;
386 if (daughterSolid->Inside(md) == kInside)
387 {
388 check_encapsulation = false;
389 G4double distout = daughterSolid->DistanceToOut(md);
390 if (distout < tol) continue; // too small overlap
391 ++overlapCount;
392 if (distout <= overlapSize) continue;
393 overlapSize = distout;
394 overlapPoint = md;
395 }
396 }
397 }
398 else // transformation with rotation
399 {
400 G4ThreeVector pmin(pmin_local), pmax(pmax_local);
401 G4ThreeVector dcenter = Td.TransformPoint(0.5*(pmin + pmax));
402 G4double dradius = 0.5*((pmax - pmin).mag());
403 if ((scenter - dcenter).mag2() >= (sradius + dradius)*(sradius + dradius)) continue;
404 if (dcenter.x() - dradius >= xmax) continue;
405 if (dcenter.y() - dradius >= ymax) continue;
406 if (dcenter.z() - dradius >= zmax) continue;
407 if (dcenter.x() + dradius <= xmin) continue;
408 if (dcenter.y() + dradius <= ymin) continue;
409 if (dcenter.z() + dradius <= zmin) continue;
410
411 G4ThreeVector pbox[8] = {
412 G4ThreeVector(pmin.x(), pmin.y(), pmin.z()),
413 G4ThreeVector(pmax.x(), pmin.y(), pmin.z()),
414 G4ThreeVector(pmin.x(), pmax.y(), pmin.z()),
415 G4ThreeVector(pmax.x(), pmax.y(), pmin.z()),
416 G4ThreeVector(pmin.x(), pmin.y(), pmax.z()),
417 G4ThreeVector(pmax.x(), pmin.y(), pmax.z()),
418 G4ThreeVector(pmin.x(), pmax.y(), pmax.z()),
419 G4ThreeVector(pmax.x(), pmax.y(), pmax.z())
420 };
421 G4double dxmin = kInfinity, dymin = kInfinity, dzmin = kInfinity;
422 G4double dxmax = -kInfinity, dymax = -kInfinity, dzmax = -kInfinity;
423 for (G4int i = 0; i < 8; ++i)
424 {
425 G4ThreeVector p = Td.TransformPoint(pbox[i]);
426 dxmin = std::min(dxmin, p.x());
427 dymin = std::min(dymin, p.y());
428 dzmin = std::min(dzmin, p.z());
429 dxmax = std::max(dxmax, p.x());
430 dymax = std::max(dymax, p.y());
431 dzmax = std::max(dzmax, p.z());
432 }
433 if (dxmin >= xmax) continue;
434 if (dymin >= ymax) continue;
435 if (dzmin >= zmax) continue;
436 if (dxmax <= xmin) continue;
437 if (dymax <= ymin) continue;
438 if (dzmax <= zmin) continue;
439 for (G4int i = 0; i < res; ++i)
440 {
441 G4ThreeVector p = points[i];
442 if (p.x() >= dxmax) continue;
443 if (p.x() <= dxmin) continue;
444 if (p.y() >= dymax) continue;
445 if (p.y() <= dymin) continue;
446 if (p.z() >= dzmax) continue;
447 if (p.z() <= dzmin) continue;
448 G4ThreeVector md = Td.InverseTransformPoint(p);
449 if (daughterSolid->Inside(md) == kInside)
450 {
451 check_encapsulation = false;
452 G4double distout = daughterSolid->DistanceToOut(md);
453 if (distout < tol) continue; // too small overlap
454 ++overlapCount;
455 if (distout <= overlapSize) continue;
456 overlapSize = distout;
457 overlapPoint = md;
458 }
459 }
460 }
461
462 // Print information on overlap
463 //
464 if (overlapCount > 0)
465 {
466 ++trials;
467 retval = true;
468 std::ostringstream message;
469 message << "Overlap with volume already placed !" << G4endl
470 << " Overlap is detected for volume "
471 << GetName() << ':' << GetCopyNo()
472 << " (" << solid->GetEntityType() << ") with "
473 << daughter->GetName() << ':' << daughter->GetCopyNo()
474 << " (" << daughterSolid->GetEntityType() << ")" << G4endl
475 << " overlap at local point " << overlapPoint
476 << " by " << G4BestUnit(overlapSize, "Length")
477 << " (max of " << overlapCount << " cases)";
478 if (trials >= maxErr)
479 {
480 message << G4endl
481 << "NOTE: Reached maximum fixed number -" << maxErr
482 << "- of overlaps reports for this volume !";
483 }
484 G4Exception("G4PVPlacement::CheckOverlaps()",
485 "GeomVol1002", JustWarning, message);
486 if (trials >= maxErr) { return true; }
487 }
488 else if (check_encapsulation)
489 {
490 // Now checking that 'sister' volume is not totally included
491 // and overlapping. Generate a single point inside of
492 // the 'sister' volume and verify that the point in NOT inside
493 // the current volume
494 //
495 G4ThreeVector pSurface = daughterSolid->GetPointOnSurface();
496 G4ThreeVector normal = daughterSolid->SurfaceNormal(pSurface);
497 G4ThreeVector pInside = pSurface - normal*1.e-4; // move point to inside
498 G4ThreeVector dPoint = (daughterSolid->Inside(pInside) == kInside) ?
499 pInside : pSurface;
500
501 // Transform the generated point to the mother's coordinate system
502 // and then to current volume's coordinate system
503 //
504 G4ThreeVector mp2 = Td.TransformPoint(dPoint);
505 G4ThreeVector msi = Tm.InverseTransformPoint(mp2);
506
507 if (solid->Inside(msi) == kInside)
508 {
509 ++trials;
510 retval = true;
511 std::ostringstream message;
512 message << "Overlap with volume already placed !" << G4endl
513 << " Overlap is detected for volume "
514 << GetName() << ':' << GetCopyNo()
515 << " (" << solid->GetEntityType() << ")" << G4endl
516 << " apparently fully encapsulating volume "
517 << daughter->GetName() << ':' << daughter->GetCopyNo()
518 << " (" << daughterSolid->GetEntityType() << ")"
519 << " at the same level!";
520 if (trials >= maxErr)
521 {
522 message << G4endl
523 << "NOTE: Reached maximum fixed number -" << maxErr
524 << "- of overlaps reports for this volume !";
525 }
526 G4Exception("G4PVPlacement::CheckOverlaps()",
527 "GeomVol1002", JustWarning, message);
528 if (trials >= maxErr) { return true; }
529 }
530 }
531 }
532
533 if (verbose && trials == 0) { G4cout << "OK! " << G4endl; }
534 return retval;
535}
@ JustWarning
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
double mag() const
G4VSolid * GetSolid() const
size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
const G4String & GetName() const
G4int GetCopyNo() const
G4LogicalVolume * GetMotherLogical() const
const G4ThreeVector GetTranslation() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4GeometryType GetEntityType() const =0
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68
@ kSurface
Definition: geomdefs.hh:69
static const G4double kInfinity
Definition: geomdefs.hh:41
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References G4VSolid::BoundingLimits(), G4VSolid::DistanceToIn(), G4VSolid::DistanceToOut(), G4BestUnit, G4cout, G4endl, G4Exception(), GetCopyNo(), G4VPhysicalVolume::GetCopyNo(), G4LogicalVolume::GetDaughter(), G4VSolid::GetEntityType(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetMotherLogical(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4LogicalVolume::GetNoDaughters(), G4VSolid::GetPointOnSurface(), G4VPhysicalVolume::GetRotation(), G4LogicalVolume::GetSolid(), G4VPhysicalVolume::GetTranslation(), G4VSolid::Inside(), G4AffineTransform::InverseTransformPoint(), G4AffineTransform::IsRotated(), JustWarning, kInfinity, kInside, kOutside, kSurface, CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4INCL::Math::min(), G4AffineTransform::NetTranslation(), CLHEP::normal(), G4VSolid::SurfaceNormal(), G4AffineTransform::TransformPoint(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by export_G4PVPlacement(), and G4PVPlacement().

◆ Clean()

void G4VPhysicalVolume::Clean ( )
staticinherited

◆ DeduceVolumeType()

EVolume G4VPhysicalVolume::DeduceVolumeType ( ) const
inlineinherited

◆ GetCopyNo()

G4int G4PVPlacement::GetCopyNo ( ) const
inlinevirtual

Implements G4VPhysicalVolume.

Definition at line 117 of file G4PVPlacement.hh.

117{ return fcopyNo; }

References fcopyNo.

Referenced by CheckOverlaps().

◆ GetFrameRotation()

const G4RotationMatrix * G4VPhysicalVolume::GetFrameRotation ( ) const
inherited

◆ GetFrameTranslation()

G4ThreeVector G4VPhysicalVolume::GetFrameTranslation ( ) const
inherited

Definition at line 213 of file G4VPhysicalVolume.cc.

214{
216}
#define G4MT_ty
#define G4MT_tz
#define G4MT_tx

References G4MT_tx, G4MT_ty, and G4MT_tz.

◆ GetInstanceID()

G4int G4VPhysicalVolume::GetInstanceID ( ) const
inlineinherited

◆ GetLogicalVolume()

G4LogicalVolume * G4VPhysicalVolume::GetLogicalVolume ( ) const
inlineinherited

Referenced by G4AdjointCrossSurfChecker::AddanExtSurfaceOfAvolume(), G4VSceneHandler::AddCompound(), G4VtkSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4OpenGLStoredSceneHandler::AddPrimitivePreambleInternal(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4ParallelWorldProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength(), G4ReplicaNavigation::BackLocate(), G4Region::BelongsTo(), G4PhantomParameterisation::BuildContainerSolid(), G4SmartVoxelHeader::BuildNodes(), G4GeometryManager::BuildOptimisations(), G4Track::CalculateVelocityForOpticalPhoton(), G4PVDivision::CheckAndSetParameters(), CheckOverlaps(), G4GeometryWorkspace::CloneReplicaSolid(), G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), G4VPVParameterisation::ComputeMaterial(), G4Navigator::ComputeSafety(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4NormalNavigation::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4ReplicaNavigation::ComputeSafety(), G4ParameterisedNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4VNestedParameterisation::ComputeSolid(), G4VPVParameterisation::ComputeSolid(), G4PhantomParameterisation::ComputeSolid(), G4ParameterisedNavigation::ComputeStep(), G4RegularNavigation::ComputeStep(), G4VoxelNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4NormalNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4PropagatorInField::ComputeStep(), G4RegularNavigation::ComputeStepSkippingEqualMaterials(), G4tgbVolume::ConstructG4Volumes(), G4TheRayTracer::CreateBitMap(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4Radioactivation::DecayIt(), G4RadioactiveDecay::DecayIt(), G4AdjointPosOnPhysVolGenerator::DefinePhysicalVolume(), G4RunManagerKernel::DefineWorldVolume(), G4GeometryManager::DeleteOptimisations(), G4LogicalVolumeModel::DescribeYourselfTo(), G4PhysicalVolumeModel::DescribeYourselfTo(), G4VFieldModel::DescribeYourselfTo(), G4GeometryWorkspace::DestroyWorkspace(), G4GlobalFastSimulationManager::DisplayRegion(), G4GDMLWriteStructure::DivisionvolWrite(), G4VisManager::Draw(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbGeometryDumper::DumpPhysVol(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4TrajectoryOriginVolumeFilter::Evaluate(), export_G4VPhysicalVolume(), G4PropagatorInField::FindAndSetFieldManager(), G4VReadOutGeometry::FindROTouchable(), G4FastTrack::FRecordsAffineTransformation(), G4Mesh::G4Mesh(), G4PVDivision::G4PVDivision(), G4PVParameterised::G4PVParameterised(), G4PVPlacement(), G4PVReplica::G4PVReplica(), G4ReplicatedSlice::G4ReplicatedSlice(), G4VExternalPhysicalVolume::G4VExternalPhysicalVolume(), G4GDMLRead::GeneratePhysvolName(), G4RTPrimaryGeneratorAction::GeneratePrimaries(), G4Navigator::GetGlobalExitNormal(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4VIntersectionLocator::GetLocalSurfaceNormal(), G4LogicalVolume::GetMass(), G4Channeling::GetMatData(), G4Channeling::GetMeanFreePath(), G4VXTRenergyLoss::GetMeanFreePath(), G4NeutrinoElectronProcess::GetMeanFreePath(), G4ElNeutrinoNucleusProcess::GetMeanFreePath(), G4MuNeutrinoNucleusProcess::GetMeanFreePath(), G4VTransitionRadiation::GetMeanFreePath(), G4Navigator::GetMotherToDaughterTransform(), G4ITNavigator1::GetMotherToDaughterTransform(), G4ITNavigator2::GetMotherToDaughterTransform(), G4TransportationManager::GetParallelWorld(), G4ITTransportationManager::GetParallelWorld(), G4tgbGeometryDumper::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4GeometryWorkspace::InitialisePhysicalVolumes(), G4BOptnForceCommonTruncatedExp::Initialize(), G4IStore::IsInWorld(), G4WeightWindowStore::IsInWorld(), G4VEnergyLossProcess::IsRegionForCubcutProcessor(), G4ParameterisedNavigation::LevelLocate(), G4RegularNavigation::LevelLocate(), G4LatticeManager::LoadLattice(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointWithinVolume(), G4ITNavigator1::LocateGlobalPointWithinVolume(), G4ITNavigator2::LocateGlobalPointWithinVolume(), G4FastSimHitMaker::make(), GFlashHitMaker::make(), G4GDMLWriteParamvol::ParametersWrite(), G4GDMLWriteParamvol::ParamvolAlgorithmWrite(), G4GDMLWriteParamvol::ParamvolWrite(), G4GDMLWriteStructure::PhysvolWrite(), G4ParallelWorldProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4Channeling::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4ForwardXrayTR::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4ITTransportation::PostStepDoIt(), G4VTransitionRadiation::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4Transportation::PostStepDoIt(), G4MaxTimeCuts::PostStepGetPhysicalInteractionLength(), G4MinEkineCuts::PostStepGetPhysicalInteractionLength(), G4BiasingProcessInterface::PostStepGetPhysicalInteractionLength(), G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength(), G4StepLimiter::PostStepGetPhysicalInteractionLength(), G4UserSpecialCuts::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4LowECapture::PostStepGetPhysicalInteractionLength(), G4NavigationLogger::PreComputeStepLog(), G4PSFlatSurfaceCurrent::ProcessHits(), G4PSFlatSurfaceFlux::ProcessHits(), G4PSSphereSurfaceFlux::ProcessHits(), G4PSVolumeFlux::ProcessHits(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4ITNavigator2::RecheckDistanceToCurrentBoundary(), G4ReflectionFactory::ReflectPVDivision(), G4ReflectionFactory::ReflectPVPlacement(), G4ReflectionFactory::ReflectPVReplica(), G4LatticeManager::RegisterLattice(), G4GDMLWriteStructure::ReplicavolWrite(), G4PropagatorInField::ReportLoopingParticle(), G4NavigationLogger::ReportOutsideMother(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4VoxelSafety::SafetyForVoxelNode(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4PolarizedIonisationModel::SampleSecondaries(), G4DNAIRT_geometries::Sampling(), G4ProductionCutsTable::ScanAndSetCouple(), G4Region::ScanVolumeTree(), G4LogicalVolume::SetFieldManager(), G4ITStepProcessor::SetInitialStep(), G4SteppingManager::SetInitialStep(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandsTouchable::SetNewValue(), G4RTPrimaryGeneratorAction::SetUp(), G4ScoringProbe::SetupGeometry(), G4ScoringBox::SetupGeometry(), G4ScoringCylinder::SetupGeometry(), G4Navigator::SetupHierarchy(), G4ITNavigator1::SetupHierarchy(), G4ITNavigator2::SetupHierarchy(), G4GlobalFastSimulationManager::ShowSetup(), G4SteppingManager::Stepping(), G4ParallelWorldProcess::SwitchMaterial(), G4GeomTestVolume::TestOverlapInTree(), G4GeomTestVolume::TestRecursiveOverlap(), G4LogicalVolume::TotalVolumeEntities(), G4GDMLWriteStructure::TraverseVolumeTree(), G4Channeling::UpdateParameters(), and G4MSSteppingAction::UserSteppingAction().

◆ GetMotherLogical()

G4LogicalVolume * G4VPhysicalVolume::GetMotherLogical ( ) const
inlineinherited

◆ GetMultiplicity()

G4int G4VPhysicalVolume::GetMultiplicity ( ) const
virtualinherited

◆ GetName()

const G4String & G4VPhysicalVolume::GetName ( ) const
inlineinherited

Referenced by G4TransportationManager::ActivateNavigator(), G4VSceneHandler::AddCompound(), G4LogicalVolume::AddDaughter(), G4HepRepFileSceneHandler::AddHepRepInstance(), G4GDMLWrite::AddModule(), G4GMocrenFileSceneHandler::AddPrimitive(), G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength(), G4CoupledTransportation::AlongStepGetPhysicalInteractionLength(), G4GDMLWriteStructure::BorderSurfaceCache(), G4SmartVoxelHeader::BuildNodes(), G4ReplicatedSlice::CheckAndSetParameters(), G4PVDivision::CheckAndSetParameters(), G4PVReplica::CheckOnlyDaughter(), G4PVParameterised::CheckOverlaps(), CheckOverlaps(), G4Navigator::CheckOverlapsIterative(), G4Navigator::ComputeSafety(), G4ITNavigator1::ComputeSafety(), G4ITNavigator2::ComputeSafety(), G4VoxelNavigation::ComputeSafety(), G4VoxelSafety::ComputeSafety(), G4PolarizedAnnihilation::ComputeSaturationFactor(), G4PolarizedCompton::ComputeSaturationFactor(), G4PolarizedIonisation::ComputeSaturationFactor(), G4ParameterisedNavigation::ComputeStep(), G4ReplicaNavigation::ComputeStep(), G4Navigator::ComputeStep(), G4ITNavigator1::ComputeStep(), G4ITNavigator2::ComputeStep(), G4PropagatorInField::ComputeStep(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), G4ImportanceConfigurator::Configure(), G4WeightCutOffConfigurator::Configure(), G4WeightWindowConfigurator::Configure(), G4tgbDetectorConstruction::Construct(), G4tgbDetectorBuilder::ConstructDetector(), G4tgbVolume::ConstructG4PhysVol(), G4FastSimulationPhysics::ConstructProcess(), G4Qt3DSceneHandler::CreateNewNode(), G4PathFinder::CreateTouchableHandle(), G4ITPathFinder::CreateTouchableHandle(), G4AdjointCrossSurfChecker::CrossingAnInterfaceBetweenTwoVolumes(), G4TransportationManager::DeActivateNavigator(), G4DNAMolecularDissociation::DecayIt(), G4PhysicalVolumeStore::DeRegister(), G4TransportationManager::DeRegisterNavigator(), G4TransportationManager::DeRegisterWorld(), G4ITTransportationManager::DeRegisterWorld(), G4PhysicalVolumeModel::DescribeAndDescend(), G4LogicalVolumeModel::DescribeYourselfTo(), G4GDMLWriteStructure::DivisionvolWrite(), G4TrajectoryDrawByOriginVolume::Draw(), G4tgbVolumeMgr::DumpG4PhysVolLeaf(), G4LogicalBorderSurface::DumpInfo(), G4tgbGeometryDumper::DumpPVParameterised(), G4tgbGeometryDumper::DumpPVPlacement(), G4tgbGeometryDumper::DumpPVReplica(), G4RunManagerKernel::DumpRegion(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4tgbVolumeMgr::DumpSummary(), G4ExceptionHandler::DumpTrackInfo(), G4ASCIITreeSceneHandler::EndModeling(), G4TrajectoryOriginVolumeFilter::Evaluate(), export_G4VPhysicalVolume(), G4FastTrack::FRecordsAffineTransformation(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4IStore::G4IStore(), G4PhysicalVolumeModel::G4PhysicalVolumeModel(), G4PVParameterised::G4PVParameterised(), G4PVReplica::G4PVReplica(), G4Navigator::GetGlobalExitNormal(), G4PSDoseDeposit3D::GetIndex(), G4PSEnergyDeposit3D::GetIndex(), G4LatticeManager::GetLattice(), G4Navigator::GetLocalExitNormal(), G4ITNavigator1::GetLocalExitNormal(), G4ModelingParameters::PVPointerCopyNo::GetName(), G4TransportationManager::GetNavigator(), G4ITTransportationManager::GetNavigator(), G4tgbVolumeMgr::GetTopPhysVol(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolume(), G4AdjointCrossSurfChecker::GoingInOrOutOfaVolumeByExtSurface(), G4SPSPosDistribution::IsSourceConfined(), G4LatticeManager::LoadLattice(), G4ITMultiNavigator::LocateGlobalPointAndSetup(), G4ITNavigator1::LocateGlobalPointAndSetup(), G4ITNavigator2::LocateGlobalPointAndSetup(), G4MultiNavigator::LocateGlobalPointAndSetup(), G4Navigator::LocateGlobalPointAndSetup(), operator<<(), G4GDMLWriteParamvol::ParametersWrite(), Path(), G4GDMLWriteStructure::PhysvolWrite(), G4UCNAbsorption::PostStepDoIt(), G4UCNBoundaryProcess::PostStepDoIt(), G4UCNMultiScattering::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4CoupledTransportation::PostStepDoIt(), G4ITSteppingVerbose::PostStepVerbose(), G4NavigationLogger::PreComputeStepLog(), G4MultiNavigator::PrepareNavigators(), G4ITMultiNavigator::PrepareNavigators(), G4ITSteppingVerbose::PreStepVerbose(), G4MultiNavigator::PrintLimited(), G4PathFinder::PrintLimited(), G4ITMultiNavigator::PrintLimited(), G4ITPathFinder::PrintLimited(), G4Navigator::PrintState(), G4ITNavigator1::PrintState(), G4PropagatorInField::printStatus(), G4PhysicalVolumeMassScene::ProcessVolume(), G4PhysicalVolumesSearchScene::ProcessVolume(), G4ReflectionFactory::ReflectPVDivision(), G4ReflectionFactory::ReflectPVParameterised(), G4ReflectionFactory::ReflectPVPlacement(), G4ReflectionFactory::ReflectPVReplica(), G4PhysicalVolumeStore::Register(), G4tgbVolumeMgr::RegisterMe(), G4PropagatorInField::ReportLoopingParticle(), G4PropagatorInField::ReportStuckParticle(), G4NavigationLogger::ReportVolumeAndIntersection(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VoxelSafety::SafetyForVoxelHeader(), G4Region::ScanVolumeTree(), G4VisCommandsTouchable::SetNewValue(), G4WeightCutOffProcess::SetParallelWorld(), G4WeightWindowProcess::SetParallelWorld(), G4ParallelWorldProcess::SetParallelWorld(), G4ParallelWorldScoringProcess::SetParallelWorld(), G4IStore::SetParallelWorldVolume(), G4IStore::SetWorldVolume(), G4WeightWindowStore::SetWorldVolume(), G4FastSimulationManagerProcess::SetWorldVolume(), G4GlobalFastSimulationManager::ShowSetup(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4GDMLRead::StripNames(), G4ErrorGeomVolumeTarget::TargetReached(), G4GeomTestVolume::TestOverlapInTree(), G4ITSteppingVerbose::TrackingEnded(), G4SteppingVerbose::TrackingStarted(), G4SteppingVerboseWithUnits::TrackingStarted(), G4ITSteppingVerbose::TrackingStarted(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), G4SteppingVerboseWithUnits::VerboseTrack(), G4VScoringMesh::WorkerConstruct(), and G4RunManagerKernel::WorkerUpdateWorldVolume().

◆ GetObjectRotation()

G4RotationMatrix * G4VPhysicalVolume::GetObjectRotation ( ) const
inherited

Definition at line 175 of file G4VPhysicalVolume.cc.

176{
177 static G4RotationMatrix aRotM;
178 static G4RotationMatrix IdentityRM;
179
180 G4RotationMatrix* retval = &IdentityRM;
181
182 // Insure against frot being a null pointer
183 if(this->GetRotation())
184 {
185 aRotM = GetRotation()->inverse();
186 retval= &aRotM;
187 }
188 return retval;
189}

References G4VPhysicalVolume::GetRotation(), and CLHEP::HepRotation::inverse().

◆ GetObjectRotationValue()

G4RotationMatrix G4VPhysicalVolume::GetObjectRotationValue ( ) const
inherited

Definition at line 191 of file G4VPhysicalVolume.cc.

192{
193 G4RotationMatrix aRotM; // Initialised to identity
194
195 // Insure against G4MT_rot being a null pointer
196 if(G4MT_rot)
197 {
198 aRotM = G4MT_rot->inverse();
199 }
200 return aRotM;
201}

References G4MT_rot.

Referenced by G4PhysicalVolumeModel::CreateCurrentAttValues(), export_G4VPhysicalVolume(), G4GDMLWriteParamvol::ParametersWrite(), and G4ReflectionFactory::ReflectPVPlacement().

◆ GetObjectTranslation()

G4ThreeVector G4VPhysicalVolume::GetObjectTranslation ( ) const
inherited

◆ GetParameterisation()

G4VPVParameterisation * G4PVPlacement::GetParameterisation ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 196 of file G4PVPlacement.cc.

197{
198 return nullptr;
199}

◆ GetRegularStructureId()

G4int G4PVPlacement::GetRegularStructureId ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 225 of file G4PVPlacement.cc.

226{
227 return 0;
228}

◆ GetReplicationData()

void G4PVPlacement::GetReplicationData ( EAxis axis,
G4int nReplicas,
G4double width,
G4double offset,
G4bool consuming 
) const
virtual

Implements G4VPhysicalVolume.

Definition at line 204 of file G4PVPlacement.cc.

206{
207 // No-operations
208}

◆ GetRotation() [1/2]

G4RotationMatrix * G4VPhysicalVolume::GetRotation ( )
inherited

Definition at line 165 of file G4VPhysicalVolume.cc.

166{
167 return G4MT_rot;
168}

References G4MT_rot.

◆ GetRotation() [2/2]

const G4RotationMatrix * G4VPhysicalVolume::GetRotation ( ) const
inherited

◆ GetSubInstanceManager()

const G4PVManager & G4VPhysicalVolume::GetSubInstanceManager ( )
staticinherited

Definition at line 140 of file G4VPhysicalVolume.cc.

141{
142 return subInstanceManager;
143}

References G4VPhysicalVolume::subInstanceManager.

Referenced by G4GeometryWorkspace::G4GeometryWorkspace().

◆ GetTranslation()

const G4ThreeVector G4VPhysicalVolume::GetTranslation ( ) const
inherited

◆ InitialiseWorker()

void G4VPhysicalVolume::InitialiseWorker ( G4VPhysicalVolume pMasterObject,
G4RotationMatrix pRot,
const G4ThreeVector tlate 
)
protectedinherited

Definition at line 111 of file G4VPhysicalVolume.cc.

115{
117
118 this->SetRotation( pRot ); // G4MT_rot = pRot;
119 this->SetTranslation( tlate ); // G4MT_trans = tlate;
120 // G4PhysicalVolumeStore::Register(this);
121}
void SlaveCopySubInstanceArray()
void SetTranslation(const G4ThreeVector &v)

References G4VPhysicalVolume::SetRotation(), G4VPhysicalVolume::SetTranslation(), G4GeomSplitter< T >::SlaveCopySubInstanceArray(), and G4VPhysicalVolume::subInstanceManager.

Referenced by G4PVReplica::InitialiseWorker().

◆ IsMany()

G4bool G4PVPlacement::IsMany ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 164 of file G4PVPlacement.cc.

165{
166 return fmany;
167}

References fmany.

◆ IsParameterised()

G4bool G4PVPlacement::IsParameterised ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 188 of file G4PVPlacement.cc.

189{
190 return false;
191}

◆ IsRegularStructure()

G4bool G4PVPlacement::IsRegularStructure ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 215 of file G4PVPlacement.cc.

216{
217 return false;
218}

◆ IsReplicated()

G4bool G4PVPlacement::IsReplicated ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 180 of file G4PVPlacement.cc.

181{
182 return false;
183}

◆ NewPtrRotMatrix()

G4RotationMatrix * G4PVPlacement::NewPtrRotMatrix ( const G4RotationMatrix RotMat)
staticprivate

Definition at line 548 of file G4PVPlacement.cc.

549{
550 G4RotationMatrix* pRotMatrix;
551 if ( RotMat.isIdentity() )
552 {
553 pRotMatrix = nullptr;
554 }
555 else
556 {
557 pRotMatrix = new G4RotationMatrix(RotMat);
558 }
559 return pRotMatrix;
560}
CLHEP::HepRotation G4RotationMatrix
bool isIdentity() const
Definition: Rotation.cc:167

References CLHEP::HepRotation::isIdentity().

Referenced by G4PVPlacement().

◆ operator=()

G4PVPlacement & G4PVPlacement::operator= ( const G4PVPlacement )
delete

◆ operator==()

G4bool G4VPhysicalVolume::operator== ( const G4VPhysicalVolume p) const
inlineinherited

◆ SetCopyNo()

void G4PVPlacement::SetCopyNo ( G4int  CopyNo)
virtual

Implements G4VPhysicalVolume.

Definition at line 172 of file G4PVPlacement.cc.

173{
174 fcopyNo = newCopyNo;
175}

References fcopyNo.

◆ SetLogicalVolume()

void G4VPhysicalVolume::SetLogicalVolume ( G4LogicalVolume pLogical)
inlineinherited

◆ SetMotherLogical()

void G4VPhysicalVolume::SetMotherLogical ( G4LogicalVolume pMother)
inlineinherited

◆ SetName()

void G4VPhysicalVolume::SetName ( const G4String pName)
inherited

◆ SetRotation()

void G4VPhysicalVolume::SetRotation ( G4RotationMatrix pRot)
inherited

◆ SetTranslation()

void G4VPhysicalVolume::SetTranslation ( const G4ThreeVector v)
inherited

Definition at line 155 of file G4VPhysicalVolume.cc.

156{
157 G4MT_tx=vec.x(); G4MT_ty=vec.y(); G4MT_tz=vec.z();
158}

References G4MT_tx, G4MT_ty, G4MT_tz, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4ParameterisationBoxX::ComputeTransformation(), G4ParameterisationBoxY::ComputeTransformation(), G4ParameterisationBoxZ::ComputeTransformation(), G4ParameterisationConsRho::ComputeTransformation(), G4ParameterisationConsPhi::ComputeTransformation(), G4ParameterisationConsZ::ComputeTransformation(), G4ParameterisationParaX::ComputeTransformation(), G4ParameterisationParaY::ComputeTransformation(), G4ParameterisationParaZ::ComputeTransformation(), G4ParameterisationPolyconeRho::ComputeTransformation(), G4ParameterisationPolyconePhi::ComputeTransformation(), G4ParameterisationPolyconeZ::ComputeTransformation(), G4ParameterisationPolyhedraRho::ComputeTransformation(), G4ParameterisationPolyhedraPhi::ComputeTransformation(), G4ParameterisationPolyhedraZ::ComputeTransformation(), G4ParameterisationTrdX::ComputeTransformation(), G4ParameterisationTrdY::ComputeTransformation(), G4ParameterisationTrdZ::ComputeTransformation(), G4ParameterisationTubsRho::ComputeTransformation(), G4ParameterisationTubsPhi::ComputeTransformation(), G4ParameterisationTubsZ::ComputeTransformation(), G4tgbPlaceParamCircle::ComputeTransformation(), G4tgbPlaceParamLinear::ComputeTransformation(), G4tgbPlaceParamSquare::ComputeTransformation(), G4ReplicaNavigation::ComputeTransformation(), G4PartialPhantomParameterisation::ComputeTransformation(), G4PhantomParameterisation::ComputeTransformation(), G4GDMLParameterisation::ComputeTransformation(), export_G4VPhysicalVolume(), G4VPhysicalVolume::G4VPhysicalVolume(), and G4VPhysicalVolume::InitialiseWorker().

◆ TerminateWorker()

void G4VPhysicalVolume::TerminateWorker ( G4VPhysicalVolume pMasterObject)
protectedinherited

Definition at line 134 of file G4VPhysicalVolume.cc.

135{
136}

◆ VolumeType()

EVolume G4PVPlacement::VolumeType ( ) const
virtual

Implements G4VPhysicalVolume.

Definition at line 235 of file G4PVPlacement.cc.

236{
237 return kNormal;
238}
@ kNormal
Definition: geomdefs.hh:84

References kNormal.

Field Documentation

◆ fallocatedRotM

G4bool G4PVPlacement::fallocatedRotM = false
private

Definition at line 167 of file G4PVPlacement.hh.

Referenced by G4PVPlacement(), and ~G4PVPlacement().

◆ fcopyNo

G4int G4PVPlacement::fcopyNo = 0
private

Definition at line 168 of file G4PVPlacement.hh.

Referenced by GetCopyNo(), and SetCopyNo().

◆ flmother

G4LogicalVolume* G4VPhysicalVolume::flmother = nullptr
privateinherited

Definition at line 246 of file G4VPhysicalVolume.hh.

◆ flogical

G4LogicalVolume* G4VPhysicalVolume::flogical = nullptr
privateinherited

Definition at line 242 of file G4VPhysicalVolume.hh.

◆ fmany

G4bool G4PVPlacement::fmany = false
private

Definition at line 166 of file G4PVPlacement.hh.

Referenced by IsMany().

◆ fname

G4String G4VPhysicalVolume::fname
privateinherited

Definition at line 245 of file G4VPhysicalVolume.hh.

Referenced by G4VPhysicalVolume::SetName().

◆ instanceID

G4int G4VPhysicalVolume::instanceID
protectedinherited

Definition at line 234 of file G4VPhysicalVolume.hh.

Referenced by G4VPhysicalVolume::G4VPhysicalVolume().

◆ pvdata

G4PVData* G4VPhysicalVolume::pvdata = nullptr
privateinherited

◆ subInstanceManager

G4PVManager G4VPhysicalVolume::subInstanceManager
staticprotectedinherited

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