Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DrawVoxels.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// class G4DrawVoxels implementation
27//
28// Define G4DrawVoxelsDebug for debugging information on G4cout
29//
30// 29/07/1999 first comitted version L.G.
31// --------------------------------------------------------------------
32
33#include "G4DrawVoxels.hh"
34#include "G4AffineTransform.hh"
35#include "G4SmartVoxelHeader.hh"
36#include "G4LogicalVolume.hh"
37#include "G4VSolid.hh"
38#include "G4VVisManager.hh"
39#include "G4Colour.hh"
41#include "G4TouchableHandle.hh"
42
43#define voxel_width 0
44
45// Private Constructor
46//
48{
49 fVoxelsVisAttributes[0].SetColour(G4Colour(1.,0.,0.));
50 fVoxelsVisAttributes[1].SetColour(G4Colour(0.,1.,0.));
51 fVoxelsVisAttributes[2].SetColour(G4Colour(0.,0.,1.));
52 fBoundingBoxVisAttributes.SetColour(G4Colour(.3,0.,.2));
53}
54
55// Destructor
56//
58
59// Methods that allow changing colors of the drawing
60//
62 G4VisAttributes& VA_voxelY,
63 G4VisAttributes& VA_voxelZ)
64{
65 fVoxelsVisAttributes[0] = VA_voxelX;
66 fVoxelsVisAttributes[1] = VA_voxelY;
67 fVoxelsVisAttributes[2] = VA_voxelZ;
68}
69
71{
72 fBoundingBoxVisAttributes = VA_boundingbox;
73}
74
75// --------------------------------------------------------------------
76
77void
78G4DrawVoxels::ComputeVoxelPolyhedra(const G4LogicalVolume* lv,
79 const G4SmartVoxelHeader* header,
80 G4VoxelLimits& limit,
81 G4PlacedPolyhedronList* ppl) const
82{
83 // Let's draw the selected voxelisation now !
84
85 G4VSolid* solid = lv->GetSolid();
86
87 G4double dx=kInfinity, dy=kInfinity, dz=kInfinity;
88 G4double xmax=0, xmin=0, ymax=0, ymin=0, zmax=0, zmin=0;
89
90 if (lv->GetNoDaughters()<=0)
91 {
92 return;
93 }
94
95 // Let's get the data for the voxelisation
96
97 solid->CalculateExtent(kXAxis,limit,G4AffineTransform(),xmin,xmax);
98 // G4AffineTransform() is identity
99 solid->CalculateExtent(kYAxis,limit,G4AffineTransform(),ymin,ymax);
100 // extents according to the axis of the local frame
101 solid->CalculateExtent(kZAxis,limit,G4AffineTransform(),zmin,zmax);
102 dx = xmax-xmin;
103 dy = ymax-ymin;
104 dz = zmax-zmin;
105
106 // Preparing the colored bounding polyhedronBox for the pVolume
107 //
108 G4PolyhedronBox bounding_polyhedronBox(dx*0.5,dy*0.5,dz*0.5);
109 bounding_polyhedronBox.SetVisAttributes(&fBoundingBoxVisAttributes);
110 G4ThreeVector t_centerofBoundingBox((xmin+xmax)*0.5,
111 (ymin+ymax)*0.5,
112 (zmin+zmax)*0.5);
113
114 ppl->push_back(G4PlacedPolyhedron(bounding_polyhedronBox,
115 G4Translate3D(t_centerofBoundingBox)));
116
117 G4ThreeVector t_FirstCenterofVoxelPlane;
118 const G4VisAttributes* voxelsVisAttributes = nullptr;
119
120 G4ThreeVector unit_translation_vector;
121 G4ThreeVector current_translation_vector;
122
123 switch(header->GetAxis())
124 {
125 case kXAxis:
126 dx=voxel_width;
127 unit_translation_vector=G4ThreeVector(1,0,0);
128 t_FirstCenterofVoxelPlane=G4ThreeVector(xmin,(ymin+ymax)*0.5,
129 (zmin+zmax)*0.5);
130 voxelsVisAttributes=&fVoxelsVisAttributes[0];
131 break;
132 case kYAxis:
133 dy=voxel_width;
134 t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,ymin,
135 (zmin+zmax)*0.5);
136 unit_translation_vector=G4ThreeVector(0,1,0);
137 voxelsVisAttributes=&fVoxelsVisAttributes[1];
138 break;
139 case kZAxis:
140 dz=voxel_width;
141 t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,
142 (ymin+ymax)*0.5,zmin);
143 unit_translation_vector=G4ThreeVector(0,0,1);
144 voxelsVisAttributes=&fVoxelsVisAttributes[2];
145 break;
146 default:
147 break;
148 };
149
150 G4PolyhedronBox voxel_plane(dx*0.5,dy*0.5,dz*0.5);
151 voxel_plane.SetVisAttributes(voxelsVisAttributes);
152
153 G4SmartVoxelProxy* slice = header->GetSlice(0);
154 std::size_t slice_no = 0, no_slices = header->GetNoSlices();
155 G4double beginning = header->GetMinExtent(),
156 step = (header->GetMaxExtent()-beginning)/no_slices;
157
158 while (slice_no<no_slices)
159 {
160 if (slice->IsHeader())
161 {
162 G4VoxelLimits newlimit(limit);
163 newlimit.AddLimit(header->GetAxis(), beginning+step*slice_no,
164 beginning+step*(slice->GetHeader()->GetMaxEquivalentSliceNo()+1));
165 ComputeVoxelPolyhedra(lv,slice->GetHeader(), newlimit, ppl);
166 }
167 current_translation_vector = unit_translation_vector;
168 current_translation_vector *= step*slice_no;
169
170 ppl->push_back(G4PlacedPolyhedron(voxel_plane,
171 G4Translate3D(current_translation_vector
172 + t_FirstCenterofVoxelPlane)));
173 slice_no = (slice->IsHeader()
174 ? slice->GetHeader()->GetMaxEquivalentSliceNo()+1
175 : slice->GetNode()->GetMaxEquivalentSliceNo()+1);
176 if (slice_no<no_slices) { slice=header->GetSlice(slice_no); }
177 }
178}
179
180// --------------------------------------------------------------------
181
184{
185 auto pplist = new G4PlacedPolyhedronList;
186 G4VoxelLimits limits; // Working object for recursive call.
187 ComputeVoxelPolyhedra(lv,lv->GetVoxelHeader(),limits,pplist);
188 return pplist; //it s up to the calling program to destroy it then!
189}
190
191// --------------------------------------------------------------------
192
194{
196
197 if (lv->GetNoDaughters()<=0)
198 {
199 return;
200 }
201
202 // Computing the transformation according to the world volume
203 // (the drawing is directly in the world volume while the axis
204 // are relative to the mother volume of lv's daughter.)
205
206 G4TouchableHandle aTouchable =
208 GetNavigatorForTracking()->CreateTouchableHistoryHandle();
209 G4AffineTransform globTransform =
210 aTouchable->GetHistory()->GetTopTransform().Inverse();
211 G4Transform3D transf3D(globTransform.NetRotation(),
212 globTransform.NetTranslation());
213
215 if(pVVisManager != nullptr)
216 {
217 // Drawing the bounding and voxel polyhedra for the pVolume
218 //
219 for (const auto & i : *pplist)
220 {
221 pVVisManager->Draw(i.GetPolyhedron(),
222 i.GetTransform()*transf3D);
223 }
224 }
225 else
226 {
227 G4Exception("G4DrawVoxels::DrawVoxels()",
228 "GeomNav1002", JustWarning,
229 "Pointer to visualization manager is null!");
230 }
231 delete pplist;
232}
#define voxel_width
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
G4AffineTransform Inverse() const
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
void DrawVoxels(const G4LogicalVolume *lv) const
void SetVoxelsVisAttributes(G4VisAttributes &, G4VisAttributes &, G4VisAttributes &)
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
void SetBoundingBoxVisAttributes(G4VisAttributes &)
G4VSolid * GetSolid() const
std::size_t GetNoDaughters() const
G4SmartVoxelHeader * GetVoxelHeader() const
const G4AffineTransform & GetTopTransform() const
G4int GetMaxEquivalentSliceNo() const
std::size_t GetNoSlices() const
G4double GetMaxExtent() const
G4double GetMinExtent() const
G4SmartVoxelProxy * GetSlice(std::size_t n) const
EAxis GetAxis() const
G4int GetMaxEquivalentSliceNo() const
G4SmartVoxelNode * GetNode() const
G4SmartVoxelHeader * GetHeader() const
G4bool IsHeader() const
virtual const G4NavigationHistory * GetHistory() const
static G4TransportationManager * GetTransportationManager()
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
void SetColour(const G4Colour &)
@ kYAxis
Definition geomdefs.hh:56
@ kXAxis
Definition geomdefs.hh:55
@ kZAxis
Definition geomdefs.hh:57