Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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"
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}
60
61// Methods that allow changing colors of the drawing
62//
64 G4VisAttributes& VA_voxelY,
65 G4VisAttributes& VA_voxelZ)
66{
67 fVoxelsVisAttributes[0] = VA_voxelX;
68 fVoxelsVisAttributes[1] = VA_voxelY;
69 fVoxelsVisAttributes[2] = VA_voxelZ;
70}
71
73{
74 fBoundingBoxVisAttributes = VA_boundingbox;
75}
76
77// --------------------------------------------------------------------
78
79void
80G4DrawVoxels::ComputeVoxelPolyhedra(const G4LogicalVolume* lv,
81 const G4SmartVoxelHeader* header,
82 G4VoxelLimits& limit,
83 G4PlacedPolyhedronList* ppl) const
84{
85 // Let's draw the selected voxelisation now !
86
87 G4VSolid* solid = lv->GetSolid();
88
89 G4double dx=kInfinity, dy=kInfinity, dz=kInfinity;
90 G4double xmax=0, xmin=0, ymax=0, ymin=0, zmax=0, zmin=0;
91
92 if (lv->GetNoDaughters()<=0)
93 {
94 return;
95 }
96
97 // Let's get the data for the voxelisation
98
99 solid->CalculateExtent(kXAxis,limit,G4AffineTransform(),xmin,xmax);
100 // G4AffineTransform() is identity
101 solid->CalculateExtent(kYAxis,limit,G4AffineTransform(),ymin,ymax);
102 // extents according to the axis of the local frame
103 solid->CalculateExtent(kZAxis,limit,G4AffineTransform(),zmin,zmax);
104 dx = xmax-xmin;
105 dy = ymax-ymin;
106 dz = zmax-zmin;
107
108 // Preparing the colored bounding polyhedronBox for the pVolume
109 //
110 G4PolyhedronBox bounding_polyhedronBox(dx*0.5,dy*0.5,dz*0.5);
111 bounding_polyhedronBox.SetVisAttributes(&fBoundingBoxVisAttributes);
112 G4ThreeVector t_centerofBoundingBox((xmin+xmax)*0.5,
113 (ymin+ymax)*0.5,
114 (zmin+zmax)*0.5);
115
116 ppl->push_back(G4PlacedPolyhedron(bounding_polyhedronBox,
117 G4Translate3D(t_centerofBoundingBox)));
118
119 G4ThreeVector t_FirstCenterofVoxelPlane;
120 const G4VisAttributes* voxelsVisAttributes = nullptr;
121
122 G4ThreeVector unit_translation_vector;
123 G4ThreeVector current_translation_vector;
124
125 switch(header->GetAxis())
126 {
127 case kXAxis:
128 dx=voxel_width;
129 unit_translation_vector=G4ThreeVector(1,0,0);
130 t_FirstCenterofVoxelPlane=G4ThreeVector(xmin,(ymin+ymax)*0.5,
131 (zmin+zmax)*0.5);
132 voxelsVisAttributes=&fVoxelsVisAttributes[0];
133 break;
134 case kYAxis:
135 dy=voxel_width;
136 t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,ymin,
137 (zmin+zmax)*0.5);
138 unit_translation_vector=G4ThreeVector(0,1,0);
139 voxelsVisAttributes=&fVoxelsVisAttributes[1];
140 break;
141 case kZAxis:
142 dz=voxel_width;
143 t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,
144 (ymin+ymax)*0.5,zmin);
145 unit_translation_vector=G4ThreeVector(0,0,1);
146 voxelsVisAttributes=&fVoxelsVisAttributes[2];
147 break;
148 default:
149 break;
150 };
151
152 G4PolyhedronBox voxel_plane(dx*0.5,dy*0.5,dz*0.5);
153 voxel_plane.SetVisAttributes(voxelsVisAttributes);
154
155 G4SmartVoxelProxy* slice = header->GetSlice(0);
156 std::size_t slice_no = 0, no_slices = header->GetNoSlices();
157 G4double beginning = header->GetMinExtent(),
158 step = (header->GetMaxExtent()-beginning)/no_slices;
159
160 while (slice_no<no_slices)
161 {
162 if (slice->IsHeader())
163 {
164 G4VoxelLimits newlimit(limit);
165 newlimit.AddLimit(header->GetAxis(), beginning+step*slice_no,
166 beginning+step*(slice->GetHeader()->GetMaxEquivalentSliceNo()+1));
167 ComputeVoxelPolyhedra(lv,slice->GetHeader(), newlimit, ppl);
168 }
169 current_translation_vector = unit_translation_vector;
170 current_translation_vector *= step*slice_no;
171
172 ppl->push_back(G4PlacedPolyhedron(voxel_plane,
173 G4Translate3D(current_translation_vector
174 + t_FirstCenterofVoxelPlane)));
175 slice_no = (slice->IsHeader()
176 ? slice->GetHeader()->GetMaxEquivalentSliceNo()+1
177 : slice->GetNode()->GetMaxEquivalentSliceNo()+1);
178 if (slice_no<no_slices) { slice=header->GetSlice(slice_no); }
179 }
180}
181
182// --------------------------------------------------------------------
183
186{
188 G4VoxelLimits limits; // Working object for recursive call.
189 ComputeVoxelPolyhedra(lv,lv->GetVoxelHeader(),limits,pplist);
190 return pplist; //it s up to the calling program to destroy it then!
191}
192
193// --------------------------------------------------------------------
194
196{
198
199 if (lv->GetNoDaughters()<=0)
200 {
201 return;
202 }
203
204 // Computing the transformation according to the world volume
205 // (the drawing is directly in the world volume while the axis
206 // are relative to the mother volume of lv's daughter.)
207
208 G4TouchableHistoryHandle aTouchable =
210 GetNavigatorForTracking()->CreateTouchableHistoryHandle();
211 G4AffineTransform globTransform =
212 aTouchable->GetHistory()->GetTopTransform().Inverse();
213 G4Transform3D transf3D(globTransform.NetRotation(),
214 globTransform.NetTranslation());
215
217 if(pVVisManager != nullptr)
218 {
219 // Drawing the bounding and voxel polyhedra for the pVolume
220 //
221 for (size_t i=0; i<pplist->size(); ++i)
222 {
223 pVVisManager->Draw((*pplist)[i].GetPolyhedron(),
224 (*pplist)[i].GetTransform()*transf3D);
225 }
226 }
227 else
228 {
229 G4Exception("G4DrawVoxels::DrawVoxels()",
230 "GeomNav1002", JustWarning,
231 "Pointer to visualization manager is null!");
232 }
233 delete pplist;
234}
#define voxel_width
Definition: G4DrawVoxels.cc:43
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:83
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
void DrawVoxels(const G4LogicalVolume *lv) const
void SetVoxelsVisAttributes(G4VisAttributes &, G4VisAttributes &, G4VisAttributes &)
Definition: G4DrawVoxels.cc:63
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
void SetBoundingBoxVisAttributes(G4VisAttributes &)
Definition: G4DrawVoxels.cc:72
G4VSolid * GetSolid() const
std::size_t GetNoDaughters() const
G4SmartVoxelHeader * GetVoxelHeader() 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
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