Geant4 11.2.2
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
 
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
 
- 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)
78, fCurrentTransform (modelTransform)
79, fBaseFullPVPath (baseFullPVPath)
81, fAbort (false)
82, fCurtailDescent (false)
85, fNClippers (0)
86{
87 fType = "G4PhysicalVolumeModel";
88
89 if (!fpTopPV) {
90
91 // In some circumstances creating an "empty" G4PhysicalVolumeModel is
92 // allowed, so I have supressed the G4Exception below. If it proves to
93 // be a problem we might have to re-instate it, but it is unlikley to
94 // be used except by visualisation experts. See, for example, /vis/list,
95 // where it is used simply to get a list of G4AttDefs.
96 // G4Exception
97 // ("G4PhysicalVolumeModel::G4PhysicalVolumeModel",
98 // "modeling0010", FatalException, "Null G4PhysicalVolumeModel pointer.");
99
100 fTopPVName = "NULL";
101 fGlobalTag = "Empty";
102 fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
103
104 } else {
105
106 fTopPVName = fpTopPV -> GetName ();
107 std::ostringstream oss;
108 oss << fpTopPV->GetName() << ':' << fpTopPV->GetCopyNo()
109 << " BasePath:" << fBaseFullPVPath;
110 fGlobalTag = oss.str();
111 fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
113 }
114}
G4Material * GetMaterial() const
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
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetCopyNo() const =0
const G4String & GetName() const

◆ ~G4PhysicalVolumeModel()

G4PhysicalVolumeModel::~G4PhysicalVolumeModel ( )
virtual

Definition at line 116 of file G4PhysicalVolumeModel.cc.

117{
118 delete fpClippingSolid;
119}

Member Function Documentation

◆ Abort()

void G4PhysicalVolumeModel::Abort ( ) const
inline

Definition at line 263 of file G4PhysicalVolumeModel.hh.

263{fAbort = true;}

◆ CalculateExtent()

void G4PhysicalVolumeModel::CalculateExtent ( )

Definition at line 143 of file G4PhysicalVolumeModel.cc.

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

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

◆ CreateCurrentAttValues()

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

Definition at line 981 of file G4PhysicalVolumeModel.cc.

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

266{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 407 of file G4PhysicalVolumeModel.cc.

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

813{
814 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
815 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
816
818
819 // Normal case - no clipping, etc. - or, illegally, more than one of those
820 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
821 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
822 sceneHandler.PostAddSolid ();
823
824 } else { // fNClippers == 1
825
826 G4VSolid* pResultantSolid = nullptr;
827 G4DisplacedSolid* pDisplacedSolid = nullptr;
828
829 if (fpClippingSolid) {
830 pDisplacedSolid = new G4DisplacedSolid("clipper", fpClippingSolid, theAT.inverse());
831 switch (fClippingMode) {
832 case subtraction:
833 if (SubtractionBoundingLimits(pSol)) {
834 pResultantSolid = new G4SubtractionSolid
835 ("subtracted_clipped_solid", pSol, pDisplacedSolid);
836 }
837 break;
838 case intersection:
839 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
840 pResultantSolid = new G4IntersectionSolid
841 ("intersected_clipped_solid", pSol, pDisplacedSolid);
842 }
843 break;
844 }
845
846 } else if (pSectionSolid) {
847 pDisplacedSolid = new G4DisplacedSolid("intersector", pSectionSolid, theAT.inverse());
848 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
849 pResultantSolid = new G4IntersectionSolid("sectioned_solid", pSol, pDisplacedSolid);
850 }
851
852 } else if (pCutawaySolid) {
853 pDisplacedSolid = new G4DisplacedSolid("cutaway", pCutawaySolid, theAT.inverse());
854 switch (fpMP->GetCutawayMode()) {
856 if (SubtractionBoundingLimits(pSol)) {
857 pResultantSolid = new G4SubtractionSolid("cutaway_solid", pSol, pDisplacedSolid);
858 }
859 break;
861 if (IntersectionBoundingLimits(pSol, pDisplacedSolid)) {
862 pResultantSolid = new G4IntersectionSolid("cutaway_solid", pSol, pDisplacedSolid);
863 }
864 break;
865 }
866 }
867
868 if (pResultantSolid) {
869 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
870 pResultantSolid -> DescribeYourselfTo (sceneHandler);
871 sceneHandler.PostAddSolid ();
872 }
873
874 delete pResultantSolid;
875 delete pDisplacedSolid;
876 }
877}
CutawayMode GetCutawayMode() const
G4DisplacedSolid * GetSectionSolid() const
G4DisplacedSolid * GetCutawaySolid() const
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 189 of file G4PhysicalVolumeModel.cc.

191{
192 if (!fpTopPV) {
194 ("G4PhysicalVolumeModel::DescribeYourselfTo",
195 "modeling0012", FatalException, "No model.");
196 return; // Should never reach here, but keeps Coverity happy.
197 }
198
199 if (!fpMP) {
201 ("G4PhysicalVolumeModel::DescribeYourselfTo",
202 "modeling0013", FatalException, "No modeling parameters.");
203 return; // Should never reach here, but keeps Coverity happy.
204 }
205
206 fNClippers = 0;
207 G4DisplacedSolid* pSectionSolid = fpMP->GetSectionSolid();
208 G4DisplacedSolid* pCutawaySolid = fpMP->GetCutawaySolid();
210 if (pSectionSolid) fNClippers++;
211 if (pCutawaySolid) fNClippers++;
212 if (fNClippers > 1) {
214 ed << "More than one solid cutter/clipper:";
215 if (fpClippingSolid) ed << "\nclipper in force";
216 if (pSectionSolid) ed << "\nsectioner in force";
217 if (pCutawaySolid) ed << "\ncutaway in force";
218 G4Exception("G4PhysicalVolumeModel::DescribeSolid", "modeling0016", JustWarning, ed);
219 }
220
221 G4Transform3D startingTransformation = fTransform;
222
223 fNTouchables.clear(); // Keeps count of touchable drawn at each depth
224
226 (fpTopPV,
228 startingTransformation,
229 sceneHandler);
230
231 // Reset or clear data...
232 fCurrentDepth = 0;
238 fDrawnPVPath.clear();
239 fAbort = false;
240 fCurtailDescent = false;
241}
@ 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 900 of file G4PhysicalVolumeModel.cc.

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

256{
257 return "G4PhysicalVolumeModel " + GetCurrentTag ();
258}

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

244{
245 if (fpCurrentPV) {
246 std::ostringstream o;
247 o << fpCurrentPV -> GetCopyNo ();
248 return fpCurrentPV -> GetName () + ":" + o.str();
249 }
250 else {
251 return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
252 }
253}

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

123{
125 for (const auto& node: path) {
126 PVNameCopyNoPath.push_back
128 (node.GetPhysicalVolume()->GetName(),node.GetCopyNo()));
129 }
130 return PVNameCopyNoPath;
131}
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 133 of file G4PhysicalVolumeModel.cc.

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

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

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

256 {
257 fClippingMode = mode;
258 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetClippingSolid()

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid * pClippingSolid)
inline

Definition at line 252 of file G4PhysicalVolumeModel.hh.

252 {
253 fpClippingSolid = pClippingSolid;
254 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetRequestedDepth()

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int requestedDepth)
inline

Definition at line 248 of file G4PhysicalVolumeModel.hh.

248 {
249 fRequestedDepth = requestedDepth;
250 }

◆ Validate()

G4bool G4PhysicalVolumeModel::Validate ( G4bool warn)
virtual

Reimplemented from G4VModel.

Definition at line 879 of file G4PhysicalVolumeModel.cc.

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

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ VisitGeometryAndGetVisReps()

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

Definition at line 260 of file G4PhysicalVolumeModel.cc.

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

Definition at line 310 of file G4PhysicalVolumeModel.hh.

Referenced by Abort(), DescribeAndDescend(), and DescribeYourselfTo().

◆ fBaseFullPVPath

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

◆ fClippingMode

ClippingMode G4PhysicalVolumeModel::fClippingMode
protected

Definition at line 313 of file G4PhysicalVolumeModel.hh.

Referenced by DescribeSolid(), 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

◆ fNClippers

G4int G4PhysicalVolumeModel::fNClippers
protected

Definition at line 314 of file G4PhysicalVolumeModel.hh.

Referenced by DescribeSolid(), and DescribeYourselfTo().

◆ 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

Definition at line 296 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent(), and VisitGeometryAndGetVisReps().

◆ fTopPVName

G4String G4PhysicalVolumeModel::fTopPVName
protected

Definition at line 295 of file G4PhysicalVolumeModel.hh.

Referenced by G4PhysicalVolumeModel().

◆ fTransform

G4Transform3D G4PhysicalVolumeModel::fTransform
protected

◆ fUseFullExtent

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 299 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent().


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