Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicalVolumeModel Class Reference

#include <G4PhysicalVolumeModel.hh>

+ Inheritance diagram for G4PhysicalVolumeModel:

Classes

class  G4PhysicalVolumeModelTouchable
 
class  G4PhysicalVolumeNodeID
 
struct  TouchableProperties
 

Public Types

enum  { UNLIMITED = -1 }
 
enum  ClippingMode { subtraction , intersection }
 

Public Member Functions

 G4PhysicalVolumeModel (G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false, const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath=std::vector< G4PhysicalVolumeNodeID >())
 
virtual ~G4PhysicalVolumeModel ()
 
void DescribeYourselfTo (G4VGraphicsScene &)
 
G4String GetCurrentDescription () const
 
G4String GetCurrentTag () const
 
G4VPhysicalVolumeGetTopPhysicalVolume () const
 
G4int GetRequestedDepth () const
 
const G4VSolidGetClippingSolid () const
 
G4int GetCurrentDepth () const
 
const G4Transform3DGetTransformation () const
 
G4VPhysicalVolumeGetCurrentPV () const
 
G4int GetCurrentPVCopyNo () const
 
G4LogicalVolumeGetCurrentLV () const
 
G4MaterialGetCurrentMaterial () const
 
const G4Transform3DGetCurrentTransform () const
 
const std::vector< G4PhysicalVolumeNodeID > & GetBaseFullPVPath () const
 
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath () const
 
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath () const
 
const std::map< G4String, G4AttDef > * GetAttDefs () const
 
std::vector< G4AttValue > * CreateCurrentAttValues () const
 
const std::map< G4int, G4int > & GetNumberOfTouchables () const
 
G4int GetTotalTouchables ()
 
void SetRequestedDepth (G4int requestedDepth)
 
void SetClippingSolid (G4VSolid *pClippingSolid)
 
void SetClippingMode (ClippingMode mode)
 
G4bool Validate (G4bool warn)
 
void Abort () const
 
void CurtailDescent () const
 
void CalculateExtent ()
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () const
 
const G4VisExtentGetExtent () const
 
const G4StringGetGlobalDescription () const
 
const G4StringGetGlobalTag () const
 
void SetModelingParameters (const G4ModelingParameters *)
 
void SetExtent (const G4VisExtent &)
 
void SetType (const G4String &)
 
void SetGlobalDescription (const G4String &)
 
void SetGlobalTag (const G4String &)
 

Static Public Member Functions

static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath (const std::vector< G4PhysicalVolumeNodeID > &)
 
static G4String GetPVNamePathString (const std::vector< G4PhysicalVolumeNodeID > &)
 

Protected Member Functions

void VisitGeometryAndGetVisReps (G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
 
void DescribeAndDescend (G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
 
virtual void DescribeSolid (const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
 

Protected Attributes

G4VPhysicalVolumefpTopPV
 
G4String fTopPVName
 
G4int fTopPVCopyNo
 
G4int fRequestedDepth
 
G4bool fUseFullExtent
 
G4Transform3D fTransform
 
G4int fCurrentDepth
 
G4VPhysicalVolumefpCurrentPV
 
G4int fCurrentPVCopyNo
 
G4LogicalVolumefpCurrentLV
 
G4MaterialfpCurrentMaterial
 
G4Transform3D fCurrentTransform
 
std::vector< G4PhysicalVolumeNodeIDfBaseFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfDrawnPVPath
 
G4bool fAbort
 
G4bool fCurtailDescent
 
G4VSolidfpClippingSolid
 
ClippingMode fClippingMode
 
G4int fNClippers
 
std::map< G4int, G4intfNTouchables
 
G4int fTotalTouchables
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
const G4ModelingParametersfpMP
 

Detailed Description

Definition at line 82 of file G4PhysicalVolumeModel.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
UNLIMITED 

Definition at line 86 of file G4PhysicalVolumeModel.hh.

◆ ClippingMode

Constructor & Destructor Documentation

◆ G4PhysicalVolumeModel()

G4PhysicalVolumeModel::G4PhysicalVolumeModel ( G4VPhysicalVolume * pVPV = 0,
G4int requestedDepth = UNLIMITED,
const G4Transform3D & modelTransformation = G4Transform3D(),
const G4ModelingParameters * pMP = 0,
G4bool useFullExtent = false,
const std::vector< G4PhysicalVolumeNodeID > & baseFullPVPath = std::vector<G4PhysicalVolumeNodeID>() )

Definition at line 60 of file G4PhysicalVolumeModel.cc.

67: G4VModel (pMP)
68, fpTopPV (pVPV)
69, fTopPVCopyNo (pVPV? pVPV->GetCopyNo(): 0)
70, fRequestedDepth (requestedDepth)
71, fUseFullExtent (useFullExtent)
72, fTransform (modelTransform)
73, fCurrentDepth (0)
75, fCurrentPVCopyNo (fpTopPV? fpTopPV->GetCopyNo(): 0)
76, fpCurrentLV (fpTopPV? fpTopPV->GetLogicalVolume(): 0)
77, fpCurrentMaterial (fpCurrentLV? fpCurrentLV->GetMaterial(): 0)
78, fCurrentTransform (modelTransform)
79, fBaseFullPVPath (baseFullPVPath)
81, fAbort (false)
82, fCurtailDescent (false)
85, fNClippers (0)
87{
88 fType = "G4PhysicalVolumeModel";
89
90 if (!fpTopPV) {
91
92 // In some circumstances creating an "empty" G4PhysicalVolumeModel is
93 // allowed, so I have supressed the G4Exception below. If it proves to
94 // be a problem we might have to re-instate it, but it is unlikley to
95 // be used except by visualisation experts. See, for example, /vis/list,
96 // where it is used simply to get a list of G4AttDefs.
97 // G4Exception
98 // ("G4PhysicalVolumeModel::G4PhysicalVolumeModel",
99 // "modeling0010", FatalException, "Null G4PhysicalVolumeModel pointer.");
100
101 fTopPVName = "NULL";
102 fGlobalTag = "Empty";
103 fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
104
105 } else {
106
107 fTopPVName = fpTopPV -> GetName ();
108 std::ostringstream oss;
109 oss << fpTopPV->GetName() << ':' << fpTopPV->GetCopyNo()
110 << " BasePath:" << fBaseFullPVPath;
111 fGlobalTag = oss.str();
112 fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
114 }
115}
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
G4VPhysicalVolume * fpCurrentPV
G4String fGlobalDescription
Definition G4VModel.hh:96
G4String fType
Definition G4VModel.hh:94
G4VModel(const G4ModelingParameters *=0)
Definition G4VModel.cc:37
G4String fGlobalTag
Definition G4VModel.hh:95
virtual G4int GetCopyNo() const =0

Referenced by G4LogicalVolumeModel::DescribeYourselfTo(), and G4LogicalVolumeModel::G4LogicalVolumeModel().

◆ ~G4PhysicalVolumeModel()

G4PhysicalVolumeModel::~G4PhysicalVolumeModel ( )
virtual

Definition at line 117 of file G4PhysicalVolumeModel.cc.

118{
119 delete fpClippingSolid;
120}

Member Function Documentation

◆ Abort()

void G4PhysicalVolumeModel::Abort ( ) const
inline

Definition at line 266 of file G4PhysicalVolumeModel.hh.

266{fAbort = true;}

◆ CalculateExtent()

void G4PhysicalVolumeModel::CalculateExtent ( )

Definition at line 144 of file G4PhysicalVolumeModel.cc.

145{
146 // To handle paramaterisations, set copy number and compute dimensions
147 // to get extent right
148 G4VPVParameterisation* pP = fpTopPV -> GetParameterisation ();
149 if (pP) {
150 fpTopPV -> SetCopyNo (fTopPVCopyNo);
151 G4VSolid* solid = pP -> ComputeSolid (fTopPVCopyNo, fpTopPV);
152 solid -> ComputeDimensions (pP, fTopPVCopyNo, fpTopPV);
153 }
154 if (fUseFullExtent) {
155 fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
156 } else {
157 // Calculate extent of *drawn* volumes, i.e., ignoring culled, e.g.,
158 // invisible volumes, by traversing the whole geometry hierarchy below
159 // this physical volume.
160 G4BoundingExtentScene beScene(this);
161 const G4int tempRequestedDepth = fRequestedDepth;
162 const G4Transform3D tempTransform = fTransform;
163 const G4ModelingParameters* tempMP = fpMP;
164 fRequestedDepth = -1; // Always search to all depths to define extent.
165 fTransform = G4Transform3D(); // Extent is in local cooridinates
166 G4ModelingParameters mParams
167 (0, // No default vis attributes needed.
168 G4ModelingParameters::wf, // wireframe (not relevant for this).
169 true, // Global culling.
170 true, // Cull invisible volumes.
171 false, // Density culling.
172 0., // Density (not relevant if density culling false).
173 true, // Cull daughters of opaque mothers.
174 24); // No of sides (not relevant for this operation).
175 mParams.SetSpecialMeshRendering(true); // Avoids traversing parameterisations
176 fpMP = &mParams;
177 DescribeYourselfTo (beScene);
178 fExtent = beScene.GetBoundingExtent();
179 fpMP = tempMP;
180 fTransform = tempTransform;
181 fRequestedDepth = tempRequestedDepth;
182 }
183 G4double radius = fExtent.GetExtentRadius();
184 if (radius < 0.) { // Nothing in the scene - revert to top extent
185 fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
186 }
187 fExtent.Transform(fTransform);
188}
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void DescribeYourselfTo(G4VGraphicsScene &)
G4VisExtent fExtent
Definition G4VModel.hh:97
const G4VisExtent & GetExtent() const
const G4ModelingParameters * fpMP
Definition G4VModel.hh:98

Referenced by G4PhysicalVolumeModel(), and G4VVisCommandGeometrySet::Set().

◆ CreateCurrentAttValues()

std::vector< G4AttValue > * G4PhysicalVolumeModel::CreateCurrentAttValues ( ) const

Definition at line 987 of file G4PhysicalVolumeModel.cc.

988{
989 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
990
991 if (!fpCurrentLV) {
993 ("G4PhysicalVolumeModel::CreateCurrentAttValues",
994 "modeling0004",
996 "Current logical volume not defined.");
997 return values;
998 }
999
1000 std::ostringstream oss; oss << fFullPVPath;
1001 values->push_back(G4AttValue("PVPath", oss.str(),""));
1002
1003 oss.str(""); oss << fBaseFullPVPath;
1004 values->push_back(G4AttValue("BasePVPath", oss.str(),""));
1005
1006 values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
1007 G4VSolid* pSol = fpCurrentLV->GetSolid();
1008
1009 values->push_back(G4AttValue("Solid", pSol->GetName(),""));
1010
1011 values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
1012
1013 oss.str(""); oss << '\n' << *pSol;
1014 values->push_back(G4AttValue("DmpSol", oss.str(),""));
1015
1016 const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
1017 const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
1018 oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
1019 values->push_back(G4AttValue("LocalTrans", oss.str(),""));
1020
1021 oss.str(""); oss << '\n' << pSol->GetExtent() << std::endl;
1022 values->push_back(G4AttValue("LocalExtent", oss.str(),""));
1023
1024 oss.str(""); oss << '\n' << fCurrentTransform;
1025 values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
1026
1027 oss.str(""); oss << '\n' << (pSol->GetExtent()).Transform(fCurrentTransform) << std::endl;
1028 values->push_back(G4AttValue("GlobalExtent", oss.str(),""));
1029
1030 G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
1031 values->push_back(G4AttValue("Material", matName,""));
1032
1033 G4double matDensity = fpCurrentMaterial? fpCurrentMaterial->GetDensity(): 0.;
1034 values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
1035
1037 oss.str(""); oss << matState;
1038 values->push_back(G4AttValue("State", oss.str(),""));
1039
1040 G4double matRadlen = fpCurrentMaterial? fpCurrentMaterial->GetRadlen(): 0.;
1041 values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
1042
1043 G4Region* region = fpCurrentLV->GetRegion();
1044 G4String regionName = region? region->GetName(): G4String("No region");
1045 values->push_back(G4AttValue("Region", regionName,""));
1046
1047 oss.str(""); oss << fpCurrentLV->IsRootRegion();
1048 values->push_back(G4AttValue("RootRegion", oss.str(),""));
1049
1050 return values;
1051}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4State
@ kStateUndefined
CLHEP::HepRotation G4RotationMatrix
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
const G4String & GetName() const
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0

Referenced by G4VSceneHandler::LoadAtts(), G4ASCIITreeSceneHandler::RequestPrimitives(), and G4VisCommandsTouchable::SetNewValue().

◆ CurtailDescent()

void G4PhysicalVolumeModel::CurtailDescent ( ) const
inline

Definition at line 269 of file G4PhysicalVolumeModel.hh.

269{fCurtailDescent = true;}

Referenced by G4ASCIITreeSceneHandler::RequestPrimitives().

◆ DescribeAndDescend()

void G4PhysicalVolumeModel::DescribeAndDescend ( G4VPhysicalVolume * pVPV,
G4int requestedDepth,
G4LogicalVolume * pLV,
G4VSolid * pSol,
G4Material * pMaterial,
const G4Transform3D & theAT,
G4VGraphicsScene & sceneHandler )
protected

Definition at line 413 of file G4PhysicalVolumeModel.cc.

421{
422 // Maintain useful data members...
423 fpCurrentPV = pVPV;
424 fCurrentPVCopyNo = pVPV->GetCopyNo();
425 fpCurrentLV = pLV;
426 fpCurrentMaterial = pMaterial;
427
428 // Create a nodeID for use below - note the "drawn" flag is true
429 G4int copyNo = fpCurrentPV->GetCopyNo();
430 auto nodeID = G4PhysicalVolumeNodeID
432
433 // Update full path of physical volumes...
434 fFullPVPath.push_back(nodeID);
435
436 const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue ();
437 const G4ThreeVector& translation = pVPV -> GetTranslation ();
438 G4Transform3D theLT (G4Transform3D (objectRotation, translation));
439
440 // Compute the accumulated transformation...
441 // Note that top volume's transformation relative to the world
442 // coordinate system is specified in theAT == startingTransformation
443 // = fTransform (see DescribeYourselfTo), so first time through the
444 // volume's own transformation, which is only relative to its
445 // mother, i.e., not relative to the world coordinate system, should
446 // not be accumulated.
447 G4Transform3D theNewAT (theAT);
448 if (fCurrentDepth != 0) theNewAT = theAT * theLT;
449 fCurrentTransform = theNewAT;
450
451 const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes();
452 // If the volume does not have any vis attributes, create it.
453 G4VisAttributes* tempVisAtts = nullptr;
454 if (!pVisAttribs) {
455 if (fpMP->GetDefaultVisAttributes()) {
456 tempVisAtts = new G4VisAttributes(*fpMP->GetDefaultVisAttributes());
457 } else {
458 tempVisAtts = new G4VisAttributes;
459 }
460 // The user may request /vis/viewer/set/colourByDensity.
461 if (fpMP->GetCBDAlgorithmNumber() == 1) {
462 // Algorithm 1: 3 parameters: Simple rainbow mapping.
463 if (fpMP->GetCBDParameters().size() != 3) {
464 G4Exception("G4PhysicalVolumeModelTouchable::DescribeAndDescend",
465 "modeling0014",
467 "Algorithm-parameter mismatch for Colour By Density");
468 } else {
469 const G4double d = pMaterial? pMaterial->GetDensity(): 0.;
470 const G4double d0 = fpMP->GetCBDParameters()[0]; // Invisible d < d0.
471 const G4double d1 = fpMP->GetCBDParameters()[1]; // Rainbow d0->d1->d2.
472 const G4double d2 = fpMP->GetCBDParameters()[2]; // Blue d > d2.
473 if (d < d0) { // Density < d0 is invisible.
474 tempVisAtts->SetVisibility(false);
475 } else { // Intermediate densities are on a spectrum.
476 G4double red, green, blue;
477 if (d < d1) {
478 red = (d1-d)/(d1-d0); green = (d-d0)/(d1-d0); blue = 0.;
479 } else if (d < d2) {
480 red = 0.; green = (d2-d)/(d2-d1); blue = (d-d1)/(d2-d1);
481 } else { // Density >= d2 is blue.
482 red = 0.; green = 0.; blue = 1.;
483 }
484 tempVisAtts->SetColour(G4Colour(red,green,blue));
485 }
486 }
487 } else if (fpMP->GetCBDAlgorithmNumber() == 2) {
488 // Algorithm 2
489 // ...etc.
490 }
491 pVisAttribs = tempVisAtts;
492 }
493 // From here, can assume pVisAttribs is a valid pointer. This is necessary
494 // because PreAddSolid needs a vis attributes object.
495
496 // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
497 const auto& vams = fpMP->GetVisAttributesModifiers();
498 if (vams.size()) {
499 // OK, we have some VAMs (Vis Attributes Modifiers).
500 for (const auto& vam: vams) {
501 const auto& vamPath = vam.GetPVNameCopyNoPath();
502 if (vamPath.size() == fFullPVPath.size()) {
503 // OK, we have a size match.
504 // Check the volume name/copy number path.
505 auto iVAMNameCopyNo = vamPath.begin();
506 auto iPVNodeId = fFullPVPath.begin();
507 for (; iVAMNameCopyNo != vamPath.end(); ++iVAMNameCopyNo, ++iPVNodeId) {
508 if (!(
509 iVAMNameCopyNo->GetName() ==
510 iPVNodeId->GetPhysicalVolume()->GetName() &&
511 iVAMNameCopyNo->GetCopyNo() ==
512 iPVNodeId->GetPhysicalVolume()->GetCopyNo()
513 )) {
514 // This path element does NOT match.
515 break;
516 }
517 }
518 if (iVAMNameCopyNo == vamPath.end()) {
519 // OK, the paths match (the above loop terminated normally).
520 // Create a vis atts object for the modified vis atts.
521 // It is static so that we may return a reliable pointer to it.
522 static G4VisAttributes modifiedVisAtts;
523 // Initialise it with the current vis atts and reset the pointer.
524 modifiedVisAtts = *pVisAttribs;
525 pVisAttribs = &modifiedVisAtts;
526 const G4VisAttributes& transVisAtts = vam.GetVisAttributes();
527 switch (vam.GetVisAttributesSignifier()) {
529 modifiedVisAtts.SetVisibility(transVisAtts.IsVisible());
530 break;
532 modifiedVisAtts.SetDaughtersInvisible
533 (transVisAtts.IsDaughtersInvisible());
534 break;
536 modifiedVisAtts.SetColour(transVisAtts.GetColour());
537 break;
539 modifiedVisAtts.SetLineStyle(transVisAtts.GetLineStyle());
540 break;
542 modifiedVisAtts.SetLineWidth(transVisAtts.GetLineWidth());
543 break;
545 if (transVisAtts.IsForceDrawingStyle()) {
546 if (transVisAtts.GetForcedDrawingStyle() ==
548 modifiedVisAtts.SetForceWireframe(true);
549 }
550 }
551 break;
553 if (transVisAtts.IsForceDrawingStyle()) {
554 if (transVisAtts.GetForcedDrawingStyle() ==
556 modifiedVisAtts.SetForceSolid(true);
557 }
558 }
559 break;
561 if (transVisAtts.IsForceDrawingStyle()) {
562 if (transVisAtts.GetForcedDrawingStyle() ==
564 modifiedVisAtts.SetForceCloud(true);
565 }
566 }
567 break;
569 modifiedVisAtts.SetForceNumberOfCloudPoints
570 (transVisAtts.GetForcedNumberOfCloudPoints());
571 break;
573 if (transVisAtts.IsForceAuxEdgeVisible()) {
574 modifiedVisAtts.SetForceAuxEdgeVisible
575 (transVisAtts.IsForcedAuxEdgeVisible());
576 }
577 break;
579 modifiedVisAtts.SetForceLineSegmentsPerCircle
580 (transVisAtts.GetForcedLineSegmentsPerCircle());
581 break;
582 }
583 }
584 }
585 }
586 }
587
588 // Check for special mesh rendering
589 if (fpMP->IsSpecialMeshRendering()) {
590 G4bool potentialG4Mesh = false;
591 if (fpMP->GetSpecialMeshVolumes().empty()) {
592 // No volumes specified - all are potentially possible
593 potentialG4Mesh = true;
594 } else {
595 // Name and (optionally) copy number of container volume is specified
596 for (const auto& pvNameCopyNo: fpMP->GetSpecialMeshVolumes()) {
597 if (pVPV->GetName() == pvNameCopyNo.GetName()) {
598 // We have a name match
599 if (pvNameCopyNo.GetCopyNo() < 0) {
600 // Any copy number is OK
601 potentialG4Mesh = true;
602 } else {
603 if (pVPV->GetCopyNo() == pvNameCopyNo.GetCopyNo()) {
604 // We have a name and copy number match
605 potentialG4Mesh = true;
606 }
607 }
608 }
609 }
610 }
611 if (potentialG4Mesh) {
612 // Create - or at least attempt to create - a mesh. If it cannot be created
613 // out of this pVPV the type will be "invalid".
614 G4Mesh mesh(pVPV,theNewAT);
615 if (mesh.GetMeshType() != G4Mesh::invalid) {
616 // Create "artificial" nodeID to represent the replaced volumes
617 G4int artCopyNo = 0;
618 auto artPV = mesh.GetParameterisedVolume();
619 auto artDepth = fCurrentDepth + 1;
620 auto artNodeID = G4PhysicalVolumeNodeID(artPV,artCopyNo,artDepth);
621 fFullPVPath.push_back(artNodeID);
622 fDrawnPVPath.push_back(artNodeID);
623 sceneHandler.AddCompound(mesh);
624 fFullPVPath.pop_back();
625 fDrawnPVPath.pop_back();
626 delete tempVisAtts; // Needs cleaning up (Coverity warning!!)
627 return; // Mesh found and processed - nothing more to do.
628 } // else continue processing
629 }
630 }
631
632 // Make decision to draw...
633 G4bool thisToBeDrawn = true;
634
635 // There are various reasons why this volume
636 // might not be drawn...
637 G4bool culling = fpMP->IsCulling();
638 G4bool cullingInvisible = fpMP->IsCullingInvisible();
639 G4bool markedVisible
640 = pVisAttribs->IsVisible() && pVisAttribs->GetColour().GetAlpha() > 0;
641 G4bool cullingLowDensity = fpMP->IsDensityCulling();
642 G4double density = pMaterial? pMaterial->GetDensity(): 0;
643 G4double densityCut = fpMP -> GetVisibleDensity ();
644
645 // 1) Global culling is on....
646 if (culling) {
647 // 2) Culling of invisible volumes is on...
648 if (cullingInvisible) {
649 // 3) ...and the volume is marked not visible...
650 if (!markedVisible) thisToBeDrawn = false;
651 }
652 // 4) Or culling of low density volumes is on...
653 if (cullingLowDensity) {
654 // 5) ...and density is less than cut value...
655 if (density < densityCut) thisToBeDrawn = false;
656 }
657 }
658 // 6) The user has asked for all further traversing to be aborted...
659 if (fAbort) thisToBeDrawn = false;
660
661 // Set "drawn" flag (it was true by default) - thisToBeDrawn may be false
662 nodeID.SetDrawn(thisToBeDrawn);
663
664 if (thisToBeDrawn) {
665
666 // Update path of drawn physical volumes...
667 fDrawnPVPath.push_back(nodeID);
668
669 if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
670 // For top-level drawn volumes, explode along radius...
671 G4Transform3D centering = G4Translate3D(fpMP->GetExplodeCentre());
672 G4Transform3D centred = centering.inverse() * theNewAT;
673 G4Scale3D oldScale;
674 G4Rotate3D oldRotation;
675 G4Translate3D oldTranslation;
676 centred.getDecomposition(oldScale, oldRotation, oldTranslation);
677 G4double explodeFactor = fpMP->GetExplodeFactor();
678 G4Translate3D newTranslation =
679 G4Translate3D(explodeFactor * oldTranslation.dx(),
680 explodeFactor * oldTranslation.dy(),
681 explodeFactor * oldTranslation.dz());
682 theNewAT = centering * newTranslation * oldRotation * oldScale;
683 }
684
685 auto fullDepth = fCurrentDepth + (G4int)fBaseFullPVPath.size();
686 fNTouchables[fullDepth]++; // Increment for every touchable drawn at each depth
687
688 DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
689
690 }
691
692 // Make decision to draw daughters, if any. There are various
693 // reasons why daughters might not be drawn...
694
695 // First, reasons that do not depend on culling policy...
696 G4int nDaughters = (G4int)pLV->GetNoDaughters();
697 G4bool daughtersToBeDrawn = true;
698 // 1) There are no daughters...
699 if (!nDaughters) daughtersToBeDrawn = false;
700 // 2) We are at the limit if requested depth...
701 else if (requestedDepth == 0) daughtersToBeDrawn = false;
702 // 3) The user has asked for all further traversing to be aborted...
703 else if (fAbort) daughtersToBeDrawn = false;
704 // 4) The user has asked that the descent be curtailed...
705 else if (fCurtailDescent) daughtersToBeDrawn = false;
706
707 // Now, reasons that depend on culling policy...
708 else {
709 G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
710 // Culling of covered daughters request. This is computed in
711 // G4VSceneHandler::CreateModelingParameters() depending on view
712 // parameters...
713 G4bool cullingCovered = fpMP->IsCullingCovered();
714 G4bool surfaceDrawing =
715 fpMP->GetDrawingStyle() == G4ModelingParameters::hsr ||
716 fpMP->GetDrawingStyle() == G4ModelingParameters::hlhsr;
717 if (pVisAttribs->IsForceDrawingStyle()) {
718 switch (pVisAttribs->GetForcedDrawingStyle()) {
719 default:
720 case G4VisAttributes::wireframe: surfaceDrawing = false; break;
721 case G4VisAttributes::solid: surfaceDrawing = true; break;
722 }
723 }
724 G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
725 // 5) Global culling is on....
726 if (culling) {
727 // 6) ..and culling of invisible volumes is on...
728 if (cullingInvisible) {
729 // 7) ...and the mother requests daughters invisible
730 if (daughtersInvisible) daughtersToBeDrawn = false;
731 }
732 // 8) Or culling of covered daughters is requested...
733 if (cullingCovered) {
734 // 9) ...and surface drawing is operating...
735 if (surfaceDrawing) {
736 // 10) ...but only if mother is visible...
737 if (thisToBeDrawn) {
738 // 11) ...and opaque...
739 if (opaque) daughtersToBeDrawn = false;
740 }
741 }
742 }
743 }
744 }
745
746 if (daughtersToBeDrawn) {
747 for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
748 // Store daughter pVPV in local variable ready for recursion...
749 G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
750 // Descend the geometry structure recursively...
753 (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
755 }
756 }
757
758 // Clean up
759 delete tempVisAtts;
760
761 // Reset for normal descending of next volume at this level...
762 fCurtailDescent = false;
763
764 // Pop item from paths physical volumes...
765 fFullPVPath.pop_back();
766 if (thisToBeDrawn) {
767 fDrawnPVPath.pop_back();
768 }
769}
@ FatalErrorInArgument
HepGeom::Translate3D G4Translate3D
HepGeom::Scale3D G4Scale3D
HepGeom::Rotate3D G4Rotate3D
bool G4bool
Definition G4Types.hh:86
G4double GetAlpha() const
Definition G4Colour.hh:173
const G4VisAttributes * GetVisAttributes() const
std::size_t GetNoDaughters() const
G4double GetDensity() const
@ invalid
Definition G4Mesh.hh:52
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
virtual void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
std::map< G4int, G4int > fNTouchables
virtual void AddCompound(const G4VTrajectory &)=0
const G4String & GetName() const
G4int GetForcedNumberOfCloudPoints() const
G4double GetLineWidth() const
G4bool IsDaughtersInvisible() const
void SetColour(const G4Colour &)
void SetVisibility(G4bool=true)
void SetForceAuxEdgeVisible(G4bool=true)
void SetForceCloud(G4bool=true)
G4int GetForcedLineSegmentsPerCircle() const
void SetForceWireframe(G4bool=true)
void SetLineWidth(G4double)
LineStyle GetLineStyle() const
const G4Colour & GetColour() const
G4bool IsVisible() const
G4bool IsForceAuxEdgeVisible() const
G4bool IsForcedAuxEdgeVisible() const
ForcedDrawingStyle GetForcedDrawingStyle() const
void SetForceSolid(G4bool=true)
void SetLineStyle(LineStyle)
void SetForceLineSegmentsPerCircle(G4int nSegments)
void SetDaughtersInvisible(G4bool=true)
void SetForceNumberOfCloudPoints(G4int nPoints)
G4bool IsForceDrawingStyle() const
double dy() const
double dz() const
double dx() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Transform3D inverse() const

Referenced by VisitGeometryAndGetVisReps().

◆ DescribeSolid()

void G4PhysicalVolumeModel::DescribeSolid ( const G4Transform3D & theAT,
G4VSolid * pSol,
const G4VisAttributes * pVisAttribs,
G4VGraphicsScene & sceneHandler )
protectedvirtual

Reimplemented in G4LogicalVolumeModel.

Definition at line 814 of file G4PhysicalVolumeModel.cc.

819{
820 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
821 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
822
824
825 // Normal case - no clipping, etc. - or, illegally, more than one of those
826 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
827 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
828 sceneHandler.PostAddSolid ();
829
830 } else { // fNClippers == 1
831
832 G4VSolid* pResultantSolid = nullptr;
833 G4DisplacedSolid* pDisplacedSolid = nullptr;
834
835 if (fpClippingSolid) {
836 pDisplacedSolid = new G4DisplacedSolid("clipper", fpClippingSolid, theAT.inverse());
837 switch (fClippingMode) {
838 case subtraction:
839 if (SubtractionBoundingLimits(pSol)) {
840 pResultantSolid = new G4SubtractionSolid
841 ("subtracted_clipped_solid", pSol, pDisplacedSolid);
842 }
843 break;
844 case intersection:
845 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
846 pResultantSolid = new G4IntersectionSolid
847 ("intersected_clipped_solid", pSol, pDisplacedSolid);
848 }
849 break;
850 }
851
852 } else if (pSectionSolid) {
853 pDisplacedSolid = new G4DisplacedSolid("intersector", pSectionSolid, theAT.inverse());
854 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
855 pResultantSolid = new G4IntersectionSolid("sectioned_solid", pSol, pDisplacedSolid);
856 }
857
858 } else if (pCutawaySolid) {
859 pDisplacedSolid = new G4DisplacedSolid("cutaway", pCutawaySolid, theAT.inverse());
860 switch (fpMP->GetCutawayMode()) {
862 if (SubtractionBoundingLimits(pSol)) {
863 pResultantSolid = new G4SubtractionSolid("cutaway_solid", pSol, pDisplacedSolid);
864 }
865 break;
867 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
868 pResultantSolid = new G4IntersectionSolid("cutaway_solid", pSol, pDisplacedSolid);
869 }
870 break;
871 }
872 }
873
874 if (pResultantSolid) {
875 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
876 pResultantSolid -> DescribeYourselfTo (sceneHandler);
877 sceneHandler.PostAddSolid ();
878 }
879
880 delete pResultantSolid;
881 delete pDisplacedSolid;
882 }
883}
virtual void PostAddSolid()=0
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0

Referenced by DescribeAndDescend().

◆ DescribeYourselfTo()

void G4PhysicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene & sceneHandler)
virtual

Implements G4VModel.

Definition at line 190 of file G4PhysicalVolumeModel.cc.

192{
193 if (!fpTopPV) {
195 ("G4PhysicalVolumeModel::DescribeYourselfTo",
196 "modeling0012", FatalException, "No model.");
197 return; // Should never reach here, but keeps Coverity happy.
198 }
199
200 if (!fpMP) {
202 ("G4PhysicalVolumeModel::DescribeYourselfTo",
203 "modeling0013", FatalException, "No modeling parameters.");
204 return; // Should never reach here, but keeps Coverity happy.
205 }
206
207 fNClippers = 0;
208 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
209 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
211 if (pSectionSolid) fNClippers++;
212 if (pCutawaySolid) fNClippers++;
213 if (fNClippers > 1) {
215 ed << "More than one solid cutter/clipper:";
216 if (fpClippingSolid) ed << "\nclipper in force";
217 if (pSectionSolid) ed << "\nsectioner in force";
218 if (pCutawaySolid) ed << "\ncutaway in force";
219 G4Exception("G4PhysicalVolumeModel::DescribeSolid", "modeling0016", JustWarning, ed);
220 }
221
222 G4Transform3D startingTransformation = fTransform;
223
224 fNTouchables.clear(); // Keeps count of touchable drawn at each depth
225
227 (fpTopPV,
229 startingTransformation,
230 sceneHandler);
231
233 for (const auto& entry : fNTouchables) {
234 fTotalTouchables += entry.second;
235 }
236
237 // Reset or clear data...
238 fCurrentDepth = 0;
240 fCurrentPVCopyNo = fpTopPV->GetCopyNo();
241 fpCurrentLV = fpTopPV->GetLogicalVolume();
242 fpCurrentMaterial = fpCurrentLV? fpCurrentLV->GetMaterial(): 0;
244 fDrawnPVPath.clear();
245 fAbort = false;
246 fCurtailDescent = false;
247}
@ FatalException
std::ostringstream G4ExceptionDescription

Referenced by CalculateExtent(), DescribeSolid(), G4LogicalVolumeModel::DescribeYourselfTo(), G4VSceneHandler::Draw3DRectMeshAsDots(), G4VSceneHandler::Draw3DRectMeshAsSurfaces(), G4VisManager::DrawGeometry(), G4VSceneHandler::DrawTetMeshAsDots(), G4VSceneHandler::DrawTetMeshAsSurfaces(), G4ASCIITreeSceneHandler::EndModeling(), G4TouchableUtils::FindTouchableProperties(), G4VisCommandSceneAddLocalAxes::SetNewValue(), G4VisCommandSceneAddVolume::SetNewValue(), G4VisCommandSetTouchable::SetNewValue(), G4VisCommandSetVolumeForField::SetNewValue(), G4VisCommandsTouchable::SetNewValue(), G4VisCommandViewerCentreOn::SetNewValue(), and G4VtkUnstructuredGridPipeline::SetUnstructuredGridData().

◆ GetAttDefs()

const std::map< G4String, G4AttDef > * G4PhysicalVolumeModel::GetAttDefs ( ) const

Definition at line 906 of file G4PhysicalVolumeModel.cc.

907{
908 G4bool isNew;
909 std::map<G4String,G4AttDef>* store
910 = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
911 if (isNew) {
912 (*store)["PVPath"] =
913 G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
914 (*store)["BasePVPath"] =
915 G4AttDef("BasePVPath","Base Physical Volume Path","Physics","","G4String");
916 (*store)["LVol"] =
917 G4AttDef("LVol","Logical Volume","Physics","","G4String");
918 (*store)["Solid"] =
919 G4AttDef("Solid","Solid Name","Physics","","G4String");
920 (*store)["EType"] =
921 G4AttDef("EType","Entity Type","Physics","","G4String");
922 (*store)["DmpSol"] =
923 G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
924 (*store)["LocalTrans"] =
925 G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
926 (*store)["LocalExtent"] =
927 G4AttDef("LocalExtent","Local extent of volume","Physics","","G4String");
928 (*store)["GlobalTrans"] =
929 G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
930 (*store)["GlobalExtent"] =
931 G4AttDef("GlobalExtent","Global extent of volume","Physics","","G4String");
932 (*store)["Material"] =
933 G4AttDef("Material","Material Name","Physics","","G4String");
934 (*store)["Density"] =
935 G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
936 (*store)["State"] =
937 G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
938 (*store)["Radlen"] =
939 G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
940 (*store)["Region"] =
941 G4AttDef("Region","Cuts Region","Physics","","G4String");
942 (*store)["RootRegion"] =
943 G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
944 }
945 return store;
946}
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Referenced by G4VSceneHandler::LoadAtts(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4VisCommandList::SetNewValue(), and G4VisCommandsTouchable::SetNewValue().

◆ GetBaseFullPVPath()

const std::vector< G4PhysicalVolumeNodeID > & G4PhysicalVolumeModel::GetBaseFullPVPath ( ) const
inline

Definition at line 203 of file G4PhysicalVolumeModel.hh.

204 {return fBaseFullPVPath;}

◆ GetClippingSolid()

const G4VSolid * G4PhysicalVolumeModel::GetClippingSolid ( ) const
inline

Definition at line 179 of file G4PhysicalVolumeModel.hh.

180 {return fpClippingSolid;}

◆ GetCurrentDepth()

G4int G4PhysicalVolumeModel::GetCurrentDepth ( ) const
inline

Definition at line 182 of file G4PhysicalVolumeModel.hh.

182{return fCurrentDepth;}

Referenced by G4GMocrenFileSceneHandler::AddSolid().

◆ GetCurrentDescription()

G4String G4PhysicalVolumeModel::GetCurrentDescription ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 261 of file G4PhysicalVolumeModel.cc.

262{
263 return "G4PhysicalVolumeModel " + GetCurrentTag ();
264}

◆ GetCurrentLV()

G4LogicalVolume * G4PhysicalVolumeModel::GetCurrentLV ( ) const
inline

◆ GetCurrentMaterial()

G4Material * G4PhysicalVolumeModel::GetCurrentMaterial ( ) const
inline

◆ GetCurrentPV()

G4VPhysicalVolume * G4PhysicalVolumeModel::GetCurrentPV ( ) const
inline

◆ GetCurrentPVCopyNo()

G4int G4PhysicalVolumeModel::GetCurrentPVCopyNo ( ) const
inline

Definition at line 191 of file G4PhysicalVolumeModel.hh.

191{return fCurrentPVCopyNo;}

◆ GetCurrentTag()

G4String G4PhysicalVolumeModel::GetCurrentTag ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 249 of file G4PhysicalVolumeModel.cc.

250{
251 if (fpCurrentPV) {
252 std::ostringstream o;
253 o << fpCurrentPV -> GetCopyNo ();
254 return fpCurrentPV -> GetName () + ":" + o.str();
255 }
256 else {
257 return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
258 }
259}

Referenced by GetCurrentDescription().

◆ GetCurrentTransform()

const G4Transform3D & G4PhysicalVolumeModel::GetCurrentTransform ( ) const
inline

◆ GetDrawnPVPath()

const std::vector< G4PhysicalVolumeNodeID > & G4PhysicalVolumeModel::GetDrawnPVPath ( ) const
inline

◆ GetFullPVPath()

const std::vector< G4PhysicalVolumeNodeID > & G4PhysicalVolumeModel::GetFullPVPath ( ) const
inline

◆ GetNumberOfTouchables()

const std::map< G4int, G4int > & G4PhysicalVolumeModel::GetNumberOfTouchables ( ) const
inline

Definition at line 245 of file G4PhysicalVolumeModel.hh.

245{return fNTouchables;}

◆ GetPVNameCopyNoPath()

G4ModelingParameters::PVNameCopyNoPath G4PhysicalVolumeModel::GetPVNameCopyNoPath ( const std::vector< G4PhysicalVolumeNodeID > & path)
static

Definition at line 122 of file G4PhysicalVolumeModel.cc.

124{
126 for (const auto& node: path) {
127 PVNameCopyNoPath.push_back
128 (G4ModelingParameters::PVNameCopyNo
129 (node.GetPhysicalVolume()->GetName(),node.GetCopyNo()));
130 }
131 return PVNameCopyNoPath;
132}
std::vector< PVNameCopyNo > PVNameCopyNoPath

Referenced by G4VViewer::TouchableSetColour(), G4VViewer::TouchableSetVisibility(), and G4VVisCommand::Twinkle().

◆ GetPVNamePathString()

G4String G4PhysicalVolumeModel::GetPVNamePathString ( const std::vector< G4PhysicalVolumeNodeID > & path)
static

Definition at line 134 of file G4PhysicalVolumeModel.cc.

138{
139 std::ostringstream oss;
140 oss << path;
141 return oss.str();
142}

Referenced by G4VViewer::TouchableSetColour(), and G4VViewer::TouchableSetVisibility().

◆ GetRequestedDepth()

G4int G4PhysicalVolumeModel::GetRequestedDepth ( ) const
inline

Definition at line 177 of file G4PhysicalVolumeModel.hh.

177{return fRequestedDepth;}

Referenced by G4ASCIITreeSceneHandler::EndModeling().

◆ GetTopPhysicalVolume()

G4VPhysicalVolume * G4PhysicalVolumeModel::GetTopPhysicalVolume ( ) const
inline

◆ GetTotalTouchables()

G4int G4PhysicalVolumeModel::GetTotalTouchables ( )
inline

Definition at line 248 of file G4PhysicalVolumeModel.hh.

248{return fTotalTouchables;}

◆ GetTransformation()

const G4Transform3D & G4PhysicalVolumeModel::GetTransformation ( ) const
inline

Definition at line 185 of file G4PhysicalVolumeModel.hh.

185{return fTransform;}

◆ SetClippingMode()

void G4PhysicalVolumeModel::SetClippingMode ( ClippingMode mode)
inline

Definition at line 259 of file G4PhysicalVolumeModel.hh.

259 {
260 fClippingMode = mode;
261 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetClippingSolid()

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid * pClippingSolid)
inline

Definition at line 255 of file G4PhysicalVolumeModel.hh.

255 {
256 fpClippingSolid = pClippingSolid;
257 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetRequestedDepth()

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int requestedDepth)
inline

Definition at line 251 of file G4PhysicalVolumeModel.hh.

251 {
252 fRequestedDepth = requestedDepth;
253 }

◆ Validate()

G4bool G4PhysicalVolumeModel::Validate ( G4bool warn)
virtual

Reimplemented from G4VModel.

Definition at line 885 of file G4PhysicalVolumeModel.cc.

886{
887// Not easy to see how to validate this sort of model. Previously there was
888// a check that a volume of the same name (fTopPVName) existed somewhere in
889// the geometry tree but under some circumstances this consumed lots of CPU
890// time. Instead, let us simply check that the volume (fpTopPV) exists in the
891// physical volume store.
892 const auto& pvStore = G4PhysicalVolumeStore::GetInstance();
893 auto iterator = find(pvStore->begin(),pvStore->end(),fpTopPV);
894 if (iterator == pvStore->end()) {
895 if (warn) {
897 ed << "Attempt to validate a volume that is no longer in the physical volume store.";
898 G4Exception("G4PhysicalVolumeModel::Validate", "modeling0015", JustWarning, ed);
899 }
900 return false;
901 } else {
902 return true;
903 }
904}
static G4PhysicalVolumeStore * GetInstance()

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ VisitGeometryAndGetVisReps()

void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps ( G4VPhysicalVolume * pVPV,
G4int requestedDepth,
const G4Transform3D & theAT,
G4VGraphicsScene & sceneHandler )
protected

Definition at line 266 of file G4PhysicalVolumeModel.cc.

271{
272 // Visits geometry structure to a given depth (requestedDepth), starting
273 // at given physical volume with given starting transformation and
274 // describes volumes to the scene handler.
275 // requestedDepth < 0 (default) implies full visit.
276 // theAT is the Accumulated Transformation.
277
278 // Find corresponding logical volume and (later) solid, storing in
279 // local variables to preserve re-entrancy.
280 G4LogicalVolume* pLV = pVPV -> GetLogicalVolume ();
281 G4VSolid* pSol = nullptr;
282 G4Material* pMaterial = nullptr;
283
284 if (!(pVPV -> IsReplicated ())) {
285 // Non-replicated physical volume.
286 pSol = pLV -> GetSolid ();
287 pMaterial = pLV -> GetMaterial ();
288 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
289 theAT, sceneHandler);
290 }
291 else {
292 // Replicated or parametrised physical volume.
293 EAxis axis;
294 G4int nReplicas;
295 G4double width;
297 G4bool consuming;
298 pVPV -> GetReplicationData (axis, nReplicas, width, offset, consuming);
299 G4int nBegin = 0;
300 G4int nEnd = nReplicas;
301 if (fCurrentDepth == 0) { // i.e., top volume
302 nBegin = fTopPVCopyNo; // Describe only one volume, namely the one
303 nEnd = nBegin + 1; // specified by the given copy number.
304 }
305 G4VPVParameterisation* pP = pVPV -> GetParameterisation ();
306 if (pP) { // Parametrised volume.
307 for (int n = nBegin; n < nEnd; n++) {
308 pSol = pP -> ComputeSolid (n, pVPV);
309 pP -> ComputeTransformation (n, pVPV);
310 pSol -> ComputeDimensions (pP, n, pVPV);
311 pVPV -> SetCopyNo (n);
313 // Create a touchable of current parent for ComputeMaterial.
314 // fFullPVPath has not been updated yet so at this point it
315 // corresponds to the parent.
317 pMaterial = pP -> ComputeMaterial (n, pVPV, &parentTouchable);
318 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
319 theAT, sceneHandler);
320 }
321 }
322 else { // Plain replicated volume. From geometry_guide.txt...
323 // The replica's positions are claculated by means of a linear formula.
324 // Replication may occur along:
325 //
326 // o Cartesian axes (kXAxis,kYAxis,kZAxis)
327 //
328 // The replications, of specified width have coordinates of
329 // form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1
330 // for the case of kXAxis, and are unrotated.
331 //
332 // o Radial axis (cylindrical polar) (kRho)
333 //
334 // The replications are cons/tubs sections, centred on the origin
335 // and are unrotated.
336 // They have radii of width*n+offset to width*(n+1)+offset
337 // where n=0..nReplicas-1
338 //
339 // o Phi axis (cylindrical polar) (kPhi)
340 // The replications are `phi sections' or wedges, and of cons/tubs form
341 // They have phi of offset+n*width to offset+(n+1)*width where
342 // n=0..nReplicas-1
343 //
344 pSol = pLV -> GetSolid ();
345 pMaterial = pLV -> GetMaterial ();
346 G4ThreeVector originalTranslation = pVPV -> GetTranslation ();
347 G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation ();
348 G4double originalRMin = 0., originalRMax = 0.;
349 if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
350 originalRMin = ((G4Tubs*)pSol)->GetInnerRadius();
351 originalRMax = ((G4Tubs*)pSol)->GetOuterRadius();
352 }
353 G4bool visualisable = true;
354 for (int n = nBegin; n < nEnd; n++) {
355 G4ThreeVector translation; // Identity.
356 G4RotationMatrix rotation; // Identity - life enough for visualizing.
357 G4RotationMatrix* pRotation = 0;
358 switch (axis) {
359 default:
360 case kXAxis:
361 translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0);
362 break;
363 case kYAxis:
364 translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0);
365 break;
366 case kZAxis:
367 translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width);
368 break;
369 case kRho:
370 if (pSol->GetEntityType() == "G4Tubs") {
371 ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset);
372 ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset);
373 } else {
374 if (fpMP->IsWarning())
375 G4warn <<
376 "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:"
377 "\n built-in replicated volumes replicated in radius for "
378 << pSol->GetEntityType() <<
379 "-type\n solids (your solid \""
380 << pSol->GetName() <<
381 "\") are not visualisable."
382 << G4endl;
383 visualisable = false;
384 }
385 break;
386 case kPhi:
387 rotation.rotateZ (-(offset+(n+0.5)*width));
388 // Minus Sign because for the physical volume we need the
389 // coordinate system rotation.
390 pRotation = &rotation;
391 break;
392 }
393 pVPV -> SetTranslation (translation);
394 pVPV -> SetRotation (pRotation);
395 pVPV -> SetCopyNo (n);
397 if (visualisable) {
398 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
399 theAT, sceneHandler);
400 }
401 }
402 // Restore originals...
403 pVPV -> SetTranslation (originalTranslation);
404 pVPV -> SetRotation (pOriginalRotation);
405 if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
406 ((G4Tubs*)pSol)->SetInnerRadius(originalRMin);
407 ((G4Tubs*)pSol)->SetOuterRadius(originalRMax);
408 }
409 }
410 }
411}
G4ThreadLocal T * G4GeomSplitter< T >::offset
#define G4warn
Definition G4Scene.cc:41
#define G4endl
Definition G4ios.hh:67
HepRotation & rotateZ(double delta)
Definition Rotation.cc:87
void DescribeAndDescend(G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
EAxis
Definition geomdefs.hh:54
@ kPhi
Definition geomdefs.hh:60
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57
@ kRho
Definition geomdefs.hh:58

Referenced by DescribeAndDescend(), and DescribeYourselfTo().

Member Data Documentation

◆ fAbort

G4bool G4PhysicalVolumeModel::fAbort
mutableprotected

◆ fBaseFullPVPath

std::vector<G4PhysicalVolumeNodeID> G4PhysicalVolumeModel::fBaseFullPVPath
protected

◆ fClippingMode

ClippingMode G4PhysicalVolumeModel::fClippingMode
protected

Definition at line 316 of file G4PhysicalVolumeModel.hh.

Referenced by DescribeSolid(), G4PhysicalVolumeModel(), and SetClippingMode().

◆ fCurrentDepth

G4int G4PhysicalVolumeModel::fCurrentDepth
protected

◆ fCurrentPVCopyNo

G4int G4PhysicalVolumeModel::fCurrentPVCopyNo
protected

◆ fCurrentTransform

G4Transform3D G4PhysicalVolumeModel::fCurrentTransform
protected

◆ fCurtailDescent

G4bool G4PhysicalVolumeModel::fCurtailDescent
mutableprotected

◆ fDrawnPVPath

std::vector<G4PhysicalVolumeNodeID> G4PhysicalVolumeModel::fDrawnPVPath
protected

◆ fFullPVPath

std::vector<G4PhysicalVolumeNodeID> G4PhysicalVolumeModel::fFullPVPath
protected

◆ fNClippers

G4int G4PhysicalVolumeModel::fNClippers
protected

◆ fNTouchables

std::map<G4int,G4int> G4PhysicalVolumeModel::fNTouchables
protected

◆ fpClippingSolid

G4VSolid* G4PhysicalVolumeModel::fpClippingSolid
protected

◆ fpCurrentLV

G4LogicalVolume* G4PhysicalVolumeModel::fpCurrentLV
protected

◆ fpCurrentMaterial

G4Material* G4PhysicalVolumeModel::fpCurrentMaterial
protected

◆ fpCurrentPV

G4VPhysicalVolume* G4PhysicalVolumeModel::fpCurrentPV
protected

◆ fpTopPV

◆ fRequestedDepth

G4int G4PhysicalVolumeModel::fRequestedDepth
protected

◆ fTopPVCopyNo

G4int G4PhysicalVolumeModel::fTopPVCopyNo
protected

◆ fTopPVName

G4String G4PhysicalVolumeModel::fTopPVName
protected

Definition at line 298 of file G4PhysicalVolumeModel.hh.

Referenced by G4PhysicalVolumeModel().

◆ fTotalTouchables

G4int G4PhysicalVolumeModel::fTotalTouchables
protected

◆ fTransform

G4Transform3D G4PhysicalVolumeModel::fTransform
protected

◆ fUseFullExtent

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 302 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent(), and G4PhysicalVolumeModel().


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