Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LogicalVolumeModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//
29// John Allison 26th July 1999.
30// Model for logical volumes.
31
33
34#include "G4VSolid.hh"
35#include "G4LogicalVolume.hh"
36#include "G4PVPlacement.hh"
37#include "G4PVParameterised.hh"
40#include "G4VGraphicsScene.hh"
41#include "G4DrawVoxels.hh"
43#include "G4VReadOutGeometry.hh"
44#include "G4Circle.hh"
45
46#include <vector>
47#include <utility>
48
50(G4LogicalVolume* pLV,
51 G4int soughtDepth,
52 G4bool booleans,
53 G4bool voxels,
54 G4bool readout,
55 G4bool checkOverlaps,
56 const G4Transform3D& modelTransformation,
57 const G4ModelingParameters* pMP):
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}
87
89
90namespace {
91 // Vis attributes
92 const G4Colour highlightSolidColour(1.0,0.8,0.8);
93 const G4double highlightSolidLineWidth(10./*pixels*/);
94 const G4Colour highlightPointColour(0.5,0.5,1.0);
95 const G4double highlightPointDiameter(20./*pixels*/);
96 // Keep a vector of solid-copy number pairs to avoid duplication.
97 typedef std::pair<G4VSolid*,G4int> solidCopyNoPair;
98 std::vector<solidCopyNoPair> solidCopyNoVector;
99 void DrawSolid
100 (G4VGraphicsScene& sceneHandler,
101 G4VSolid* sol, G4int copyNo, const G4Transform3D& t) {
102 // Avoid duplication.
103 std::pair<G4VSolid*,G4int> pair(sol,copyNo);
104 auto iter = solidCopyNoVector.begin();
105 for ( ; iter != solidCopyNoVector.end(); ++iter) {
106 if (*iter == pair) break;
107 }
108 if (iter == solidCopyNoVector.end()) {
109 solidCopyNoVector.push_back(pair);
110 G4VisAttributes highlightSolidVisAtts(highlightSolidColour);
111 highlightSolidVisAtts.SetLineWidth(highlightSolidLineWidth);
112 sceneHandler.PreAddSolid(t,highlightSolidVisAtts);
113 sceneHandler.AddSolid(*sol);
114 sceneHandler.PostAddSolid();
115 }
116 }
117 void DrawPoint
118 (G4VGraphicsScene& sceneHandler,
119 const G4ThreeVector& point) {
120 G4VisAttributes highlightPointVisAtts(highlightPointColour);
121 G4Circle overlapPoint;
122 overlapPoint.SetVisAttributes(highlightPointVisAtts);
123 overlapPoint.SetPosition(point);
124 overlapPoint.SetDiameter(G4VMarker::SizeType::screen,highlightPointDiameter);
126 sceneHandler.BeginPrimitives();
127 sceneHandler.AddPrimitive(overlapPoint);
128 sceneHandler.EndPrimitives();
129 }
130}
131
133(G4VGraphicsScene& sceneHandler) {
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);
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}
289
290// This called from G4PhysicalVolumeModel::DescribeAndDescend by the
291// virtual function mechanism.
293(const G4Transform3D& theAT,
294 G4VSolid* pSol,
295 const G4VisAttributes* pVisAttribs,
296 G4VGraphicsScene& sceneHandler) {
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)
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector InverseTransformPoint(const G4ThreeVector &vec) const
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolumeModel(G4LogicalVolume *, G4int soughtDepth=1, G4bool booleans=true, G4bool voxels=true, G4bool readout=true, G4bool checkOverlaps=true, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
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 override
G4int GetMultiplicity() const override
void DescribeYourselfTo(G4VGraphicsScene &)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void PostAddSolid()=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void AddSolid(const G4Box &)=0
virtual void EndPrimitives()=0
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
void SetDiameter(SizeType, G4double)
void SetFillStyle(FillStyle)
void SetPosition(const G4Point3D &)
void SetModelingParameters(const G4ModelingParameters *)
G4String fGlobalDescription
Definition G4VModel.hh:96
G4String fType
Definition G4VModel.hh:94
const G4ModelingParameters * fpMP
Definition G4VModel.hh:98
G4String fGlobalTag
Definition G4VModel.hh:95
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
void SetForceWireframe(G4bool=true)
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68