Geant4 11.1.1
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, G4bool checkOverlaps=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, 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
 
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 ()
 
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
 
void SetModelingParameters (const G4ModelingParameters *)
 
void SetExtent (const G4VisExtent &)
 
void SetType (const G4String &)
 
void SetGlobalDescription (const G4String &)
 
void SetGlobalTag (const G4String &)
 
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)
 

Protected Attributes

G4LogicalVolumefpLV
 
G4bool fBooleans
 
G4bool fVoxels
 
G4bool fReadout
 
G4bool fCheckOverlaps
 
G4bool fOverlapsPrinted
 
- Protected Attributes inherited from G4PhysicalVolumeModel
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
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
const G4ModelingParametersfpMP
 

Additional Inherited Members

- Public Types inherited from G4PhysicalVolumeModel
enum  { UNLIMITED = -1 }
 
enum  ClippingMode { subtraction , intersection }
 
- Static Public Member Functions inherited from G4PhysicalVolumeModel
static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath (const std::vector< G4PhysicalVolumeNodeID > &)
 

Detailed Description

Definition at line 50 of file G4LogicalVolumeModel.hh.

Constructor & Destructor Documentation

◆ G4LogicalVolumeModel()

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

Definition at line 49 of file G4LogicalVolumeModel.cc.

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

◆ ~G4LogicalVolumeModel()

G4LogicalVolumeModel::~G4LogicalVolumeModel ( )
virtual

Definition at line 88 of file G4LogicalVolumeModel.cc.

88{}

Member Function Documentation

◆ DescribeSolid()

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

Reimplemented from G4PhysicalVolumeModel.

Definition at line 292 of file G4LogicalVolumeModel.cc.

296 {
297
298 if (fBooleans) {
299 // Look for "constituents". Could be a Boolean solid.
300 G4VSolid* pSol0 = pSol -> GetConstituentSolid (0);
301 if (pSol0) { // Composite solid...
302 G4VSolid* pSol1 = pSol -> GetConstituentSolid (1);
303 if (!pSol1) {
305 ("G4PhysicalVolumeModel::DescribeSolid",
306 "modeling0001", FatalException,
307 "2nd component solid in Boolean is missing.");
308 }
309 // Draw these constituents white and "forced wireframe"...
310 G4VisAttributes constituentAttributes;
311 constituentAttributes.SetForceWireframe(true);
312 DescribeSolid (theAT, pSol0, &constituentAttributes, sceneHandler);
313 DescribeSolid (theAT, pSol1, &constituentAttributes, sceneHandler);
314 }
315 }
316
317 // In any case draw the original/resultant solid...
318 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
319 pSol -> DescribeYourselfTo (sceneHandler);
320 sceneHandler.PostAddSolid ();
321}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
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=true)

Referenced by DescribeSolid().

◆ DescribeYourselfTo()

void G4LogicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene sceneHandler)
virtual

Implements G4VModel.

Definition at line 132 of file G4LogicalVolumeModel.cc.

133 {
134
135 // Store current modeling parameters and ensure nothing is culled.
136 const G4ModelingParameters* tmpMP = fpMP;
137 G4ModelingParameters nonCulledMP;
138 if (fpMP) nonCulledMP = *fpMP;
139 nonCulledMP.SetCulling (false);
140 fpMP = &nonCulledMP;
142 fpMP = tmpMP;
143
144 if (fVoxels) {
146 // Add Voxels.
147 G4DrawVoxels dv;
149 dv.CreatePlacedPolyhedra (fpTopPV -> GetLogicalVolume ());
150 for (size_t i = 0; i < pPPL -> size (); i++) {
151 const G4Transform3D& transform = (*pPPL)[i].GetTransform ();
152 const G4Polyhedron& polyhedron = (*pPPL)[i].GetPolyhedron ();
153 sceneHandler.BeginPrimitives (transform);
154 sceneHandler.AddPrimitive (polyhedron);
155 sceneHandler.EndPrimitives ();
156 }
157 delete pPPL;
158 }
159 }
160
161 if (fReadout) {
162 // Draw readout geometry...
164 if (sd) {
165 G4VReadOutGeometry* roGeom = sd->GetROgeometry();
166 if (roGeom) {
167 G4VPhysicalVolume* roWorld = roGeom->GetROWorld();
168// G4cout << "Readout geometry \"" << roGeom->GetName()
169// << "\" with top physical volume \""
170// << roWorld->GetName()
171// << "\"" << G4endl;
172 G4PhysicalVolumeModel pvModel(roWorld);
173 pvModel.SetModelingParameters(fpMP);
174 pvModel.DescribeYourselfTo(sceneHandler);
175 }
176 }
177 }
178
179 if (fCheckOverlaps) {
181 G4VSolid* motherSolid = motherLog->GetSolid();
182 G4int nDaughters = (G4int)motherLog->GetNoDaughters();
183
184 // Models are called repeatedly by the scene handler so be careful...
185 // Print overlaps - but only the first time for a given instantiation of G4LogicalVolume
186 if (!fOverlapsPrinted) {
187 for (G4int iDaughter = 0; iDaughter < nDaughters; ++iDaughter) {
188 G4VPhysicalVolume* daughterPhys = motherLog->GetDaughter(iDaughter);
189 daughterPhys->CheckOverlaps();
190 }
191 fOverlapsPrinted = true;
192 }
193
194 // Draw overlaps
195 solidCopyNoVector.clear();
196 for (G4int iDaughter = 0; iDaughter < nDaughters; ++iDaughter) {
197 G4VPhysicalVolume* daughterPhys = motherLog->GetDaughter(iDaughter);
198 G4PVPlacement* daughterPVPlace = dynamic_cast<G4PVPlacement*>(daughterPhys);
199 G4PVParameterised* daughterPVParam = dynamic_cast<G4PVParameterised*>(daughterPhys);
200 const G4int nPoints = 1000;
201
202 if (daughterPVPlace) {
203
204 // This algorithm is based on G4PVPlacement::CheckOverlaps.
205 G4AffineTransform tDaughter(daughterPhys->GetRotation(),daughterPhys->GetTranslation());
206 G4VSolid* daughterSolid = daughterPhys->GetLogicalVolume()->GetSolid();
207 for (G4int i = 0; i < nPoints; ++i) {
208 G4ThreeVector point = daughterSolid->GetPointOnSurface();
209 // Transform to mother's coordinate system
210 G4ThreeVector motherPoint = tDaughter.TransformPoint(point);
211 // Check overlaps with the mother volume
212 if (motherSolid->Inside(motherPoint)==kOutside) {
213 // Draw mother and daughter and point
214 DrawSolid(sceneHandler,motherSolid,0,G4Transform3D());
215 DrawSolid(sceneHandler,daughterSolid,daughterPhys->GetCopyNo(),tDaughter);
216 DrawPoint(sceneHandler,motherPoint);
217 }
218 // Check other daughters
219 for (G4int iSister = 0; iSister < nDaughters; ++iSister) {
220 if (iSister == iDaughter) continue;
221 G4VPhysicalVolume* sisterPhys = motherLog->GetDaughter(iSister);
222 G4AffineTransform tSister(sisterPhys->GetRotation(),sisterPhys->GetTranslation());
223 // Transform to sister's coordinate system
224 G4ThreeVector sisterPoint = tSister.InverseTransformPoint(motherPoint);
225 G4LogicalVolume* sisterLog = sisterPhys->GetLogicalVolume();
226 G4VSolid* sisterSolid = sisterLog->GetSolid();
227 if (sisterSolid->Inside(sisterPoint)==kInside) {
228 // Draw daughter and sister and point
229 DrawSolid(sceneHandler,daughterSolid,daughterPhys->GetCopyNo(),tDaughter);
230 DrawSolid(sceneHandler,sisterSolid,sisterPhys->GetCopyNo(),tSister);
231 DrawPoint(sceneHandler,motherPoint);
232 }
233 }
234 }
235
236 } else if (daughterPVParam) {
237
238 // This algorithm is based on G4PVParameterised::CheckOverlaps
239 const G4int multiplicity = daughterPVParam->GetMultiplicity();
240 auto* param = daughterPVParam->GetParameterisation();
241 // Cache points for later checking against other parameterisations
242 std::vector<G4ThreeVector> motherPoints;
243 for (G4int iP = 0; iP < multiplicity; iP++) {
244 G4VSolid* daughterSolid = param->ComputeSolid(iP, daughterPhys);
245 daughterSolid->ComputeDimensions(param, iP, daughterPhys);
246 param->ComputeTransformation(iP, daughterPhys);
247 G4AffineTransform tDaughter(daughterPVParam->GetRotation(),daughterPVParam->GetTranslation());
248 for (G4int i = 0; i < nPoints; ++i) {
249 G4ThreeVector point = daughterSolid->GetPointOnSurface();
250 // Transform to mother's coordinate system
251 G4ThreeVector motherPoint = tDaughter.TransformPoint(point);
252 // Check overlaps with the mother volume
253 if (motherSolid->Inside(motherPoint)==kOutside) {
254 // Draw mother and daughter and point
255 DrawSolid(sceneHandler,motherSolid,0,G4Transform3D());
256 DrawSolid(sceneHandler,daughterSolid,iP,tDaughter);
257 DrawPoint(sceneHandler,motherPoint);
258 }
259 motherPoints.push_back(motherPoint);
260 }
261 // Check sister parameterisations
262 for (G4int iPP = iP + 1; iPP < multiplicity; iPP++) {
263 G4VSolid* sisterSolid = param->ComputeSolid(iPP, daughterPhys);
264 sisterSolid->ComputeDimensions(param, iPP, daughterPhys);
265 param->ComputeTransformation(iPP, daughterPhys);
266 G4AffineTransform tSister
267 (daughterPVParam->GetRotation(),daughterPVParam->GetTranslation());
268 for (const auto& motherPoint: motherPoints) {
269 // Transform each point into daughter's frame
270 G4ThreeVector sisterPoint = tSister.InverseTransformPoint(motherPoint);
271 if (sisterSolid->Inside(sisterPoint)==kInside) {
272 // Draw sister
273 DrawSolid(sceneHandler,sisterSolid,iPP,tSister);
274 // Recompute daughter parameterisation before drawing
275 daughterSolid->ComputeDimensions(param, iP, daughterPhys);
276 param->ComputeTransformation(iP, daughterPhys);
277 tDaughter = G4AffineTransform
278 (daughterPVParam->GetRotation(),daughterPVParam->GetTranslation());
279 DrawSolid(sceneHandler,daughterSolid,iP,tDaughter);
280 DrawPoint(sceneHandler,motherPoint);
281 }
282 }
283 }
284 }
285 }
286 }
287 }
288}
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
HepGeom::Transform3D G4Transform3D
int G4int
Definition: G4Types.hh:85
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
G4VSolid * GetSolid() const
G4VSensitiveDetector * GetSensitiveDetector() const
std::size_t GetNoDaughters() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
G4SmartVoxelHeader * GetVoxelHeader() const
void SetCulling(G4bool)
G4VPVParameterisation * GetParameterisation() const
virtual G4int GetMultiplicity() const
Definition: G4PVReplica.cc:297
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:102
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
virtual G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int errMax=1)
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetCopyNo() const =0
G4VPhysicalVolume * GetROWorld() const
G4VReadOutGeometry * GetROgeometry() const
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68

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().

◆ fCheckOverlaps

G4bool G4LogicalVolumeModel::fCheckOverlaps
protected

Definition at line 87 of file G4LogicalVolumeModel.hh.

Referenced by DescribeYourselfTo().

◆ fOverlapsPrinted

G4bool G4LogicalVolumeModel::fOverlapsPrinted
protected

Definition at line 88 of file G4LogicalVolumeModel.hh.

Referenced by DescribeYourselfTo().

◆ 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: