Geant4 11.4.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  G4PhysicalVolumeNodeID
class  G4PhysicalVolumeModelTouchable
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 > & GetMapOfDrawnTouchables () const
G4int GetTotalDrawnTouchables ()
G4int GetMaxFullDepth () const
const std::map< G4int, G4int > & GetMapOfAllTouchables () const
G4int GetTotalAllTouchables ()
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 > &)
Static Public Member Functions inherited from G4VModel
static const G4ModelingParametersGetCurrentModelingParameters ()
static void SetCurrentModelingParameters (const G4ModelingParameters *)

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, G4intfMapDrawnTouchables
G4int fTotalDrawnTouchables
std::map< G4int, G4intfMapAllTouchables
G4int fTotalAllTouchables
G4int fMaxFullDepth
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 62 of file G4PhysicalVolumeModel.cc.

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

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

◆ ~G4PhysicalVolumeModel()

G4PhysicalVolumeModel::~G4PhysicalVolumeModel ( )
virtual

Definition at line 121 of file G4PhysicalVolumeModel.cc.

122{
123 delete fpClippingSolid;
124}

Member Function Documentation

◆ Abort()

void G4PhysicalVolumeModel::Abort ( ) const
inline

Definition at line 275 of file G4PhysicalVolumeModel.hh.

275{fAbort = true;}

◆ CalculateExtent()

void G4PhysicalVolumeModel::CalculateExtent ( )

Definition at line 148 of file G4PhysicalVolumeModel.cc.

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

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

◆ CreateCurrentAttValues()

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

Definition at line 1043 of file G4PhysicalVolumeModel.cc.

1044{
1045 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
1046
1047 if (!fpCurrentLV) {
1049 ("G4PhysicalVolumeModel::CreateCurrentAttValues",
1050 "modeling0004",
1052 "Current logical volume not defined.");
1053 return values;
1054 }
1055
1056 std::ostringstream oss; oss << fFullPVPath;
1057 values->push_back(G4AttValue("PVPath", oss.str(),""));
1058
1059 oss.str(""); oss << fBaseFullPVPath;
1060 values->push_back(G4AttValue("BasePVPath", oss.str(),""));
1061
1062 values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
1063 G4VSolid* pSol = fpCurrentLV->GetSolid();
1064
1065 values->push_back(G4AttValue("Solid", pSol->GetName(),""));
1066
1067 values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
1068
1069 oss.str(""); oss << '\n' << *pSol;
1070 values->push_back(G4AttValue("DmpSol", oss.str(),""));
1071
1072 const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
1073 const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
1074 oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
1075 values->push_back(G4AttValue("LocalTrans", oss.str(),""));
1076
1077 oss.str(""); oss << '\n' << pSol->GetExtent() << std::endl;
1078 values->push_back(G4AttValue("LocalExtent", oss.str(),""));
1079
1080 oss.str(""); oss << '\n' << fCurrentTransform;
1081 values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
1082
1083 oss.str(""); oss << '\n' << (pSol->GetExtent()).Transform(fCurrentTransform) << std::endl;
1084 values->push_back(G4AttValue("GlobalExtent", oss.str(),""));
1085
1086 G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
1087 values->push_back(G4AttValue("Material", matName,""));
1088
1089 G4double matDensity = fpCurrentMaterial? fpCurrentMaterial->GetDensity(): 0.;
1090 values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
1091
1093 oss.str(""); oss << matState;
1094 values->push_back(G4AttValue("State", oss.str(),""));
1095
1096 G4double matRadlen = fpCurrentMaterial? fpCurrentMaterial->GetRadlen(): 0.;
1097 values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
1098
1099 G4Region* region = fpCurrentLV->GetRegion();
1100 G4String regionName = region? region->GetName(): G4String("No region");
1101 values->push_back(G4AttValue("Region", regionName,""));
1102
1103 oss.str(""); oss << fpCurrentLV->IsRootRegion();
1104 values->push_back(G4AttValue("RootRegion", oss.str(),""));
1105
1106 return values;
1107}
@ 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 278 of file G4PhysicalVolumeModel.hh.

278{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 423 of file G4PhysicalVolumeModel.cc.

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

875{
876 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
877 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
878
880
881 // Normal case - no clipping, etc. - or, illegally, more than one of those
882 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
883 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
884 sceneHandler.PostAddSolid ();
885
886 } else { // fNClippers == 1
887
888 G4VSolid* pResultantSolid = nullptr;
889 G4DisplacedSolid* pDisplacedSolid = nullptr;
890
891 if (fpClippingSolid) {
892 pDisplacedSolid = new G4DisplacedSolid("clipper", fpClippingSolid, theAT.inverse());
893 switch (fClippingMode) {
894 case subtraction:
895 if (SubtractionBoundingLimits(pSol)) {
896 pResultantSolid = new G4SubtractionSolid
897 ("subtracted_clipped_solid", pSol, pDisplacedSolid);
898 }
899 break;
900 case intersection:
901 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
902 pResultantSolid = new G4IntersectionSolid
903 ("intersected_clipped_solid", pSol, pDisplacedSolid);
904 }
905 break;
906 }
907
908 } else if (pSectionSolid) {
909 pDisplacedSolid = new G4DisplacedSolid("intersector", pSectionSolid, theAT.inverse());
910 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
911 pResultantSolid = new G4IntersectionSolid("sectioned_solid", pSol, pDisplacedSolid);
912 }
913
914 } else if (pCutawaySolid) {
915 pDisplacedSolid = new G4DisplacedSolid("cutaway", pCutawaySolid, theAT.inverse());
916 switch (fpMP->GetCutawayMode()) {
918 if (SubtractionBoundingLimits(pSol)) {
919 pResultantSolid = new G4SubtractionSolid("cutaway_solid", pSol, pDisplacedSolid);
920 }
921 break;
923 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
924 pResultantSolid = new G4IntersectionSolid("cutaway_solid", pSol, pDisplacedSolid);
925 }
926 break;
927 }
928 }
929
930 if (pResultantSolid) {
931 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
932 pResultantSolid -> DescribeYourselfTo (sceneHandler);
933 sceneHandler.PostAddSolid ();
934 }
935
936 delete pResultantSolid;
937 delete pDisplacedSolid;
938 }
939}
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 194 of file G4PhysicalVolumeModel.cc.

196{
197 if (!fpTopPV) {
199 ("G4PhysicalVolumeModel::DescribeYourselfTo",
200 "modeling0012", FatalException, "No model.");
201 return; // Should never reach here, but keeps Coverity happy.
202 }
203
204 if (!fpMP) {
206 ("G4PhysicalVolumeModel::DescribeYourselfTo",
207 "modeling0013", FatalException, "No modeling parameters.");
208 return; // Should never reach here, but keeps Coverity happy.
209 }
210
211 fNClippers = 0;
212 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
213 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
215 if (pSectionSolid) fNClippers++;
216 if (pCutawaySolid) fNClippers++;
217 if (fNClippers > 1) {
219 ed << "More than one solid cutter/clipper:";
220 if (fpClippingSolid) ed << "\nclipper in force";
221 if (pSectionSolid) ed << "\nsectioner in force";
222 if (pCutawaySolid) ed << "\ncutaway in force";
223 G4Exception("G4PhysicalVolumeModel::DescribeSolid", "modeling0016", JustWarning, ed);
224 }
225
226 G4Transform3D startingTransformation = fTransform;
227
228 fMapDrawnTouchables.clear(); // Keeps count of touchables drawn at each depth
229 fMapAllTouchables.clear(); // Keeps count of all touchables at each depth
230
232 (fpTopPV,
234 startingTransformation,
235 sceneHandler);
236
238 for (const auto& entry : fMapDrawnTouchables) {
239 fTotalDrawnTouchables += entry.second;
240 }
241
243 for (const auto& entry : fMapAllTouchables) {
244 fTotalAllTouchables += entry.second;
245 }
246
247 // Reset or clear data...
248 fCurrentDepth = 0;
250 fCurrentPVCopyNo = fpTopPV->GetCopyNo();
251 fpCurrentLV = fpTopPV->GetLogicalVolume();
252 fpCurrentMaterial = fpCurrentLV? fpCurrentLV->GetMaterial(): 0;
254 fDrawnPVPath.clear();
255 fAbort = false;
256 fCurtailDescent = false;
257}
@ 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 962 of file G4PhysicalVolumeModel.cc.

963{
964 G4bool isNew;
965 std::map<G4String,G4AttDef>* store
966 = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
967 if (isNew) {
968 (*store)["PVPath"] =
969 G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
970 (*store)["BasePVPath"] =
971 G4AttDef("BasePVPath","Base Physical Volume Path","Physics","","G4String");
972 (*store)["LVol"] =
973 G4AttDef("LVol","Logical Volume","Physics","","G4String");
974 (*store)["Solid"] =
975 G4AttDef("Solid","Solid Name","Physics","","G4String");
976 (*store)["EType"] =
977 G4AttDef("EType","Entity Type","Physics","","G4String");
978 (*store)["DmpSol"] =
979 G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
980 (*store)["LocalTrans"] =
981 G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
982 (*store)["LocalExtent"] =
983 G4AttDef("LocalExtent","Local extent of volume","Physics","","G4String");
984 (*store)["GlobalTrans"] =
985 G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
986 (*store)["GlobalExtent"] =
987 G4AttDef("GlobalExtent","Global extent of volume","Physics","","G4String");
988 (*store)["Material"] =
989 G4AttDef("Material","Material Name","Physics","","G4String");
990 (*store)["Density"] =
991 G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
992 (*store)["State"] =
993 G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
994 (*store)["Radlen"] =
995 G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
996 (*store)["Region"] =
997 G4AttDef("Region","Cuts Region","Physics","","G4String");
998 (*store)["RootRegion"] =
999 G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
1000 }
1001 return store;
1002}
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 271 of file G4PhysicalVolumeModel.cc.

272{
273 return "G4PhysicalVolumeModel " + GetCurrentTag ();
274}

◆ 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 259 of file G4PhysicalVolumeModel.cc.

260{
261 if (fpCurrentPV) {
262 std::ostringstream o;
263 o << fpCurrentPV -> GetCopyNo ();
264 return fpCurrentPV -> GetName () + ":" + o.str();
265 }
266 else {
267 return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
268 }
269}
G4int GetCopyNo() const override

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

◆ GetMapOfAllTouchables()

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

Definition at line 254 of file G4PhysicalVolumeModel.hh.

254{return fMapAllTouchables;}

◆ GetMapOfDrawnTouchables()

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

Definition at line 245 of file G4PhysicalVolumeModel.hh.

245{return fMapDrawnTouchables;}

◆ GetMaxFullDepth()

G4int G4PhysicalVolumeModel::GetMaxFullDepth ( ) const
inline

Definition at line 251 of file G4PhysicalVolumeModel.hh.

251{return fMaxFullDepth;}

◆ GetPVNameCopyNoPath()

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

Definition at line 126 of file G4PhysicalVolumeModel.cc.

128{
130 for (const auto& node: path) {
131 PVNameCopyNoPath.push_back
132 (G4ModelingParameters::PVNameCopyNo
133 (node.GetPhysicalVolume()->GetName(),node.GetCopyNo()));
134 }
135 return PVNameCopyNoPath;
136}
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 138 of file G4PhysicalVolumeModel.cc.

142{
143 std::ostringstream oss;
144 oss << path;
145 return oss.str();
146}

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

◆ GetTotalAllTouchables()

G4int G4PhysicalVolumeModel::GetTotalAllTouchables ( )
inline

Definition at line 257 of file G4PhysicalVolumeModel.hh.

257{return fTotalAllTouchables;}

◆ GetTotalDrawnTouchables()

G4int G4PhysicalVolumeModel::GetTotalDrawnTouchables ( )
inline

Definition at line 248 of file G4PhysicalVolumeModel.hh.

248{return fTotalDrawnTouchables;}

◆ 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 268 of file G4PhysicalVolumeModel.hh.

268 {
269 fClippingMode = mode;
270 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetClippingSolid()

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid * pClippingSolid)
inline

Definition at line 264 of file G4PhysicalVolumeModel.hh.

264 {
265 fpClippingSolid = pClippingSolid;
266 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetRequestedDepth()

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int requestedDepth)
inline

Definition at line 260 of file G4PhysicalVolumeModel.hh.

260 {
261 fRequestedDepth = requestedDepth;
262 }

◆ Validate()

G4bool G4PhysicalVolumeModel::Validate ( G4bool warn)
virtual

Reimplemented from G4VModel.

Definition at line 941 of file G4PhysicalVolumeModel.cc.

942{
943// Not easy to see how to validate this sort of model. Previously there was
944// a check that a volume of the same name (fTopPVName) existed somewhere in
945// the geometry tree but under some circumstances this consumed lots of CPU
946// time. Instead, let us simply check that the volume (fpTopPV) exists in the
947// physical volume store.
948 const auto& pvStore = G4PhysicalVolumeStore::GetInstance();
949 auto iterator = find(pvStore->begin(),pvStore->end(),fpTopPV);
950 if (iterator == pvStore->end()) {
951 if (warn) {
953 ed << "Attempt to validate a volume that is no longer in the physical volume store.";
954 G4Exception("G4PhysicalVolumeModel::Validate", "modeling0015", JustWarning, ed);
955 }
956 return false;
957 } else {
958 return true;
959 }
960}
static G4PhysicalVolumeStore * GetInstance()

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ VisitGeometryAndGetVisReps()

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

Definition at line 276 of file G4PhysicalVolumeModel.cc.

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

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 325 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

◆ fMapAllTouchables

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

◆ fMapDrawnTouchables

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

◆ fMaxFullDepth

G4int G4PhysicalVolumeModel::fMaxFullDepth
protected

◆ fNClippers

G4int G4PhysicalVolumeModel::fNClippers
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 307 of file G4PhysicalVolumeModel.hh.

Referenced by G4PhysicalVolumeModel().

◆ fTotalAllTouchables

G4int G4PhysicalVolumeModel::fTotalAllTouchables
protected

◆ fTotalDrawnTouchables

G4int G4PhysicalVolumeModel::fTotalDrawnTouchables
protected

◆ fTransform

G4Transform3D G4PhysicalVolumeModel::fTransform
protected

◆ fUseFullExtent

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 311 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent(), and G4PhysicalVolumeModel().


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