Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4DrawVoxels Class Reference

#include <G4DrawVoxels.hh>

Public Member Functions

G4PlacedPolyhedronListCreatePlacedPolyhedra (const G4LogicalVolume *) const
 
void DrawVoxels (const G4LogicalVolume *lv) const
 
 G4DrawVoxels ()
 
void SetBoundingBoxVisAttributes (G4VisAttributes &)
 
void SetVoxelsVisAttributes (G4VisAttributes &, G4VisAttributes &, G4VisAttributes &)
 
 ~G4DrawVoxels ()
 

Private Member Functions

void ComputeVoxelPolyhedra (const G4LogicalVolume *, const G4SmartVoxelHeader *, G4VoxelLimits &, G4PlacedPolyhedronList *) const
 
 G4DrawVoxels (const G4DrawVoxels &)=delete
 
G4DrawVoxels operator= (const G4DrawVoxels &)=delete
 

Private Attributes

G4VisAttributes fBoundingBoxVisAttributes
 
G4VisAttributes fVoxelsVisAttributes [3]
 

Detailed Description

Definition at line 48 of file G4DrawVoxels.hh.

Constructor & Destructor Documentation

◆ G4DrawVoxels() [1/2]

G4DrawVoxels::G4DrawVoxels ( )

Definition at line 47 of file G4DrawVoxels.cc.

48{
53}
G4VisAttributes fVoxelsVisAttributes[3]
Definition: G4DrawVoxels.hh:83
G4VisAttributes fBoundingBoxVisAttributes
Definition: G4DrawVoxels.hh:84
void SetColour(const G4Colour &)

References fBoundingBoxVisAttributes, fVoxelsVisAttributes, and G4VisAttributes::SetColour().

◆ ~G4DrawVoxels()

G4DrawVoxels::~G4DrawVoxels ( )

Definition at line 57 of file G4DrawVoxels.cc.

58{
59}

◆ G4DrawVoxels() [2/2]

G4DrawVoxels::G4DrawVoxels ( const G4DrawVoxels )
privatedelete

Member Function Documentation

◆ ComputeVoxelPolyhedra()

void G4DrawVoxels::ComputeVoxelPolyhedra ( const G4LogicalVolume lv,
const G4SmartVoxelHeader header,
G4VoxelLimits limit,
G4PlacedPolyhedronList ppl 
) const
private

Definition at line 80 of file G4DrawVoxels.cc.

84{
85 // Let's draw the selected voxelisation now !
86
87 G4VSolid* solid = lv->GetSolid();
88
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 G4int 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}
#define voxel_width
Definition: G4DrawVoxels.cc:43
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Translate3D G4Translate3D
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void ComputeVoxelPolyhedra(const G4LogicalVolume *, const G4SmartVoxelHeader *, G4VoxelLimits &, G4PlacedPolyhedronList *) const
Definition: G4DrawVoxels.cc:80
G4VSolid * GetSolid() const
size_t GetNoDaughters() const
G4int GetMaxEquivalentSliceNo() const
G4double GetMaxExtent() const
G4double GetMinExtent() const
EAxis GetAxis() const
G4SmartVoxelProxy * GetSlice(G4int n) const
size_t GetNoSlices() const
G4int GetMaxEquivalentSliceNo() const
G4SmartVoxelNode * GetNode() const
G4SmartVoxelHeader * GetHeader() const
G4bool IsHeader() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
static const G4double kInfinity
Definition: geomdefs.hh:41

References G4VoxelLimits::AddLimit(), G4VSolid::CalculateExtent(), ComputeVoxelPolyhedra(), fBoundingBoxVisAttributes, fVoxelsVisAttributes, G4SmartVoxelHeader::GetAxis(), G4SmartVoxelProxy::GetHeader(), G4SmartVoxelHeader::GetMaxEquivalentSliceNo(), G4SmartVoxelNode::GetMaxEquivalentSliceNo(), G4SmartVoxelHeader::GetMaxExtent(), G4SmartVoxelHeader::GetMinExtent(), G4LogicalVolume::GetNoDaughters(), G4SmartVoxelProxy::GetNode(), G4SmartVoxelHeader::GetNoSlices(), G4SmartVoxelHeader::GetSlice(), G4LogicalVolume::GetSolid(), G4SmartVoxelProxy::IsHeader(), kInfinity, kXAxis, kYAxis, kZAxis, G4Visible::SetVisAttributes(), and voxel_width.

Referenced by ComputeVoxelPolyhedra(), and CreatePlacedPolyhedra().

◆ CreatePlacedPolyhedra()

G4PlacedPolyhedronList * G4DrawVoxels::CreatePlacedPolyhedra ( const G4LogicalVolume lv) const

Definition at line 185 of file G4DrawVoxels.cc.

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}
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
G4SmartVoxelHeader * GetVoxelHeader() const

References ComputeVoxelPolyhedra(), and G4LogicalVolume::GetVoxelHeader().

Referenced by G4LogicalVolumeModel::DescribeYourselfTo(), and DrawVoxels().

◆ DrawVoxels()

void G4DrawVoxels::DrawVoxels ( const G4LogicalVolume lv) const

Definition at line 195 of file G4DrawVoxels.cc.

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}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
static G4TransportationManager * GetTransportationManager()
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0

References CreatePlacedPolyhedra(), G4VVisManager::Draw(), G4Exception(), G4VVisManager::GetConcreteInstance(), G4LogicalVolume::GetNoDaughters(), G4TransportationManager::GetTransportationManager(), JustWarning, G4AffineTransform::NetRotation(), and G4AffineTransform::NetTranslation().

◆ operator=()

G4DrawVoxels G4DrawVoxels::operator= ( const G4DrawVoxels )
privatedelete

◆ SetBoundingBoxVisAttributes()

void G4DrawVoxels::SetBoundingBoxVisAttributes ( G4VisAttributes VA_boundingbox)

Definition at line 72 of file G4DrawVoxels.cc.

73{
74 fBoundingBoxVisAttributes = VA_boundingbox;
75}

References fBoundingBoxVisAttributes.

◆ SetVoxelsVisAttributes()

void G4DrawVoxels::SetVoxelsVisAttributes ( G4VisAttributes VA_voxelX,
G4VisAttributes VA_voxelY,
G4VisAttributes VA_voxelZ 
)

Definition at line 63 of file G4DrawVoxels.cc.

66{
67 fVoxelsVisAttributes[0] = VA_voxelX;
68 fVoxelsVisAttributes[1] = VA_voxelY;
69 fVoxelsVisAttributes[2] = VA_voxelZ;
70}

References fVoxelsVisAttributes.

Field Documentation

◆ fBoundingBoxVisAttributes

G4VisAttributes G4DrawVoxels::fBoundingBoxVisAttributes
private

◆ fVoxelsVisAttributes

G4VisAttributes G4DrawVoxels::fVoxelsVisAttributes[3]
private

Definition at line 83 of file G4DrawVoxels.hh.

Referenced by ComputeVoxelPolyhedra(), G4DrawVoxels(), and SetVoxelsVisAttributes().


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