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

#include <G4LogicalVolumeModel.hh>

+ Inheritance diagram for G4LogicalVolumeModel:

Public Member Functions

 G4LogicalVolumeModel (G4LogicalVolume *, G4int soughtDepth=1, G4bool booleans=true, G4bool voxels=true, G4bool readout=true, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
 
virtual ~G4LogicalVolumeModel ()
 
void DescribeYourselfTo (G4VGraphicsScene &)
 
G4bool Validate (G4bool)
 
- Public Member Functions inherited from G4PhysicalVolumeModel
 G4PhysicalVolumeModel (G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false)
 
virtual ~G4PhysicalVolumeModel ()
 
void DescribeYourselfTo (G4VGraphicsScene &)
 
G4String GetCurrentDescription () const
 
G4String GetCurrentTag () const
 
G4VPhysicalVolumeGetTopPhysicalVolume () const
 
G4int GetRequestedDepth () const
 
const G4VSolidGetClippingSolid () const
 
G4int GetCurrentDepth () const
 
G4VPhysicalVolumeGetCurrentPV () const
 
G4LogicalVolumeGetCurrentLV () const
 
G4MaterialGetCurrentMaterial () 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
 
void SetBaseFullPVPath (const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath)
 
void SetRequestedDepth (G4int requestedDepth)
 
void SetClippingSolid (G4VSolid *pClippingSolid)
 
void SetClippingMode (ClippingMode mode)
 
G4bool Validate (G4bool warn)
 
void CurtailDescent ()
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
virtual void DescribeYourselfTo (G4VGraphicsScene &)=0
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () const
 
virtual G4String GetCurrentDescription () const
 
virtual G4String GetCurrentTag () const
 
const G4VisExtentGetExtent () const
 
const G4StringGetGlobalDescription () const
 
const G4StringGetGlobalTag () const
 
const G4Transform3DGetTransformation () const
 
void SetModelingParameters (const G4ModelingParameters *)
 
void SetExtent (const G4VisExtent &)
 
void SetType (const G4String &)
 
void SetGlobalDescription (const G4String &)
 
void SetGlobalTag (const G4String &)
 
void SetTransformation (const G4Transform3D &)
 
virtual G4bool Validate (G4bool warn=true)
 

Protected Member Functions

void DescribeSolid (const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
 
- Protected Member Functions inherited from G4PhysicalVolumeModel
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)
 
void CalculateExtent ()
 

Protected Attributes

G4LogicalVolumefpLV
 
G4bool fBooleans
 
G4bool fVoxels
 
G4bool fReadout
 
- Protected Attributes inherited from G4PhysicalVolumeModel
G4VPhysicalVolumefpTopPV
 
G4String fTopPVName
 
G4int fTopPVCopyNo
 
G4int fRequestedDepth
 
G4bool fUseFullExtent
 
G4int fCurrentDepth
 
G4VPhysicalVolumefpCurrentPV
 
G4LogicalVolumefpCurrentLV
 
G4MaterialfpCurrentMaterial
 
G4Transform3DfpCurrentTransform
 
std::vector< G4PhysicalVolumeNodeIDfBaseFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfDrawnPVPath
 
G4bool fCurtailDescent
 
G4VSolidfpClippingSolid
 
ClippingMode fClippingMode
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
G4Transform3D fTransform
 
const G4ModelingParametersfpMP
 

Additional Inherited Members

- Public Types inherited from G4PhysicalVolumeModel
enum  { UNLIMITED = -1 }
 
enum  ClippingMode { subtraction , intersection }
 

Detailed Description

Definition at line 51 of file G4LogicalVolumeModel.hh.

Constructor & Destructor Documentation

◆ G4LogicalVolumeModel()

G4LogicalVolumeModel::G4LogicalVolumeModel ( G4LogicalVolume pLV,
G4int  soughtDepth = 1,
G4bool  booleans = true,
G4bool  voxels = true,
G4bool  readout = true,
const G4Transform3D modelTransformation = G4Transform3D(),
const G4ModelingParameters pMP = 0 
)

Definition at line 44 of file G4LogicalVolumeModel.cc.

51 :
52 // Instantiate a G4PhysicalVolumeModel with a G4PVPlacement to
53 // represent this logical volume. It has no rotation and a null
54 // translation so that the logical volume will be seen in its own
55 // reference system. It will be added to the physical volume store
56 // but it will not be part of the normal geometry heirarchy so it
57 // has no mother.
59(new G4PVPlacement (0, // No rotation.
60 G4ThreeVector(), // Null traslation.
61 "PhysVol representation of LogVol " + pLV -> GetName (),
62 pLV,
63 0, // No mother.
64 false, // Not "MANY".
65 0), // Copy number.
66 soughtDepth,
67 modelTransformation,
68 pMP,
69 true), // Use full extent.
70 fpLV (pLV),
71 fBooleans (booleans),
72 fVoxels (voxels),
73 fReadout (readout)
74{
75 fType = "G4LogicalVolumeModel";
76 fGlobalTag = fpLV -> GetName ();
77 fGlobalDescription = "G4LogicalVolumeModel " + fGlobalTag;
78}
CLHEP::Hep3Vector G4ThreeVector
G4String fGlobalDescription
Definition: G4VModel.hh:110
G4String fType
Definition: G4VModel.hh:108
G4String fGlobalTag
Definition: G4VModel.hh:109

◆ ~G4LogicalVolumeModel()

G4LogicalVolumeModel::~G4LogicalVolumeModel ( )
virtual

Definition at line 80 of file G4LogicalVolumeModel.cc.

80{}

Member Function Documentation

◆ DescribeSolid()

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

Reimplemented from G4PhysicalVolumeModel.

Definition at line 132 of file G4LogicalVolumeModel.cc.

136 {
137
138 if (fBooleans) {
139 // Look for "constituents". Could be a Boolean solid.
140 G4VSolid* pSol0 = pSol -> GetConstituentSolid (0);
141 if (pSol0) { // Composite solid...
142 G4VSolid* pSol1 = pSol -> GetConstituentSolid (1);
143 if (!pSol1) {
145 ("G4PhysicalVolumeModel::DescribeSolid",
146 "modeling0001", FatalException,
147 "2nd component solid in Boolean is missing.");
148 }
149 // Draw these constituents white and "forced wireframe"...
150 G4VisAttributes constituentAttributes;
151 constituentAttributes.SetForceWireframe(true);
152 DescribeSolid (theAT, pSol0, &constituentAttributes, sceneHandler);
153 DescribeSolid (theAT, pSol1, &constituentAttributes, sceneHandler);
154 }
155 }
156
157 // In any case draw the original/resultant solid...
158 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
159 pSol -> DescribeYourselfTo (sceneHandler);
160 sceneHandler.PostAddSolid ();
161}
@ FatalException
void DescribeYourselfTo(G4VGraphicsScene &)
void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
virtual void PostAddSolid()=0
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
void SetForceWireframe(G4bool)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by DescribeSolid().

◆ DescribeYourselfTo()

void G4LogicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene sceneHandler)
virtual

Implements G4VModel.

Definition at line 82 of file G4LogicalVolumeModel.cc.

83 {
84
85 // Store current modeling parameters and ensure nothing is culled.
86 const G4ModelingParameters* tmpMP = fpMP;
87 G4ModelingParameters nonCulledMP;
88 if (fpMP) nonCulledMP = *fpMP;
89 nonCulledMP.SetCulling (false);
90 fpMP = &nonCulledMP;
92 fpMP = tmpMP;
93
94 if (fVoxels) {
96 // Add Voxels.
97 G4DrawVoxels dv;
99 dv.CreatePlacedPolyhedra (fpTopPV -> GetLogicalVolume ());
100 for (size_t i = 0; i < pPPL -> size (); i++) {
101 const G4Transform3D& transform = (*pPPL)[i].GetTransform ();
102 const G4Polyhedron& polyhedron = (*pPPL)[i].GetPolyhedron ();
103 sceneHandler.BeginPrimitives (transform);
104 sceneHandler.AddPrimitive (polyhedron);
105 sceneHandler.EndPrimitives ();
106 }
107 delete pPPL;
108 }
109 }
110
111 if (fReadout) {
112 // Draw readout geometry...
114 if (sd) {
115 G4VReadOutGeometry* roGeom = sd->GetROgeometry();
116 if (roGeom) {
117 G4VPhysicalVolume* roWorld = roGeom->GetROWorld();
118 G4cout << "Readout geometry \"" << roGeom->GetName()
119 << "\" with top physical volume \""
120 << roWorld->GetName()
121 << "\"" << G4endl;
122 G4PhysicalVolumeModel pvModel(roWorld);
123 pvModel.SetModelingParameters(fpMP);
124 pvModel.DescribeYourselfTo(sceneHandler);
125 }
126 }
127 }
128}
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
G4VSensitiveDetector * GetSensitiveDetector() const
G4SmartVoxelHeader * GetVoxelHeader() const
void SetCulling(G4bool)
G4VPhysicalVolume * fpTopPV
void DescribeYourselfTo(G4VGraphicsScene &)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void EndPrimitives()=0
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4VPhysicalVolume * GetROWorld() const
G4String GetName() const
G4VReadOutGeometry * GetROgeometry() const

Referenced by DescribeSolid().

◆ Validate()

G4bool G4LogicalVolumeModel::Validate ( G4bool  )
inlinevirtual

Reimplemented from G4VModel.

Definition at line 68 of file G4LogicalVolumeModel.hh.

68{return true;}

Member Data Documentation

◆ fBooleans

G4bool G4LogicalVolumeModel::fBooleans
protected

Definition at line 84 of file G4LogicalVolumeModel.hh.

Referenced by DescribeSolid().

◆ fpLV

G4LogicalVolume* G4LogicalVolumeModel::fpLV
protected

Definition at line 83 of file G4LogicalVolumeModel.hh.

Referenced by DescribeYourselfTo(), and G4LogicalVolumeModel().

◆ fReadout

G4bool G4LogicalVolumeModel::fReadout
protected

Definition at line 86 of file G4LogicalVolumeModel.hh.

Referenced by DescribeYourselfTo().

◆ fVoxels

G4bool G4LogicalVolumeModel::fVoxels
protected

Definition at line 85 of file G4LogicalVolumeModel.hh.

Referenced by DescribeYourselfTo().


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