Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DrawVoxels Class Reference

#include <G4DrawVoxels.hh>

Public Member Functions

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

Detailed Description

Definition at line 48 of file G4DrawVoxels.hh.

Constructor & Destructor Documentation

◆ G4DrawVoxels()

G4DrawVoxels::G4DrawVoxels ( )

Definition at line 47 of file G4DrawVoxels.cc.

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}

Referenced by ~G4DrawVoxels().

◆ ~G4DrawVoxels()

G4DrawVoxels::~G4DrawVoxels ( )
default

Member Function Documentation

◆ CreatePlacedPolyhedra()

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

Definition at line 179 of file G4DrawVoxels.cc.

180{
181 auto pplist = new G4PlacedPolyhedronList;
182 G4VoxelLimits limits; // Working object for recursive call.
183 ComputeVoxelPolyhedra(lv,lv->GetVoxelHeader(),limits,pplist);
184 return pplist; //it s up to the calling program to destroy it then!
185}
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList
G4SmartVoxelHeader * GetVoxelHeader() const

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

◆ DrawVoxels()

void G4DrawVoxels::DrawVoxels ( const G4LogicalVolume * lv) const

Definition at line 189 of file G4DrawVoxels.cc.

190{
191 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
192
193 if (lv->GetNoDaughters()<=0)
194 {
195 return;
196 }
197
198 // Computing the transformation according to the world volume
199 // (the drawing is directly in the world volume while the axis
200 // are relative to the mother volume of lv's daughter.)
201
202 G4TouchableHandle aTouchable =
204 GetNavigatorForTracking()->CreateTouchableHistoryHandle();
205 G4AffineTransform globTransform =
206 aTouchable->GetHistory()->GetTopTransform().Inverse();
207 G4Transform3D transf3D(globTransform.NetRotation(),
208 globTransform.NetTranslation());
209
211 if(pVVisManager != nullptr)
212 {
213 // Drawing the bounding and voxel polyhedra for the pVolume
214 //
215 for (const auto & i : *pplist)
216 {
217 pVVisManager->Draw(i.GetPolyhedron(),
218 i.GetTransform()*transf3D);
219 }
220 }
221 else
222 {
223 G4Exception("G4DrawVoxels::DrawVoxels()",
224 "GeomNav1002", JustWarning,
225 "Pointer to visualization manager is null!");
226 }
227 delete pplist;
228}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
HepGeom::Transform3D G4Transform3D
G4AffineTransform Inverse() const
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
std::size_t GetNoDaughters() const
const G4AffineTransform & GetTopTransform() const
virtual const G4NavigationHistory * GetHistory() const
static G4TransportationManager * GetTransportationManager()
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0

◆ SetBoundingBoxVisAttributes()

void G4DrawVoxels::SetBoundingBoxVisAttributes ( G4VisAttributes & VA_boundingbox)

Definition at line 66 of file G4DrawVoxels.cc.

67{
68 fBoundingBoxVisAttributes = VA_boundingbox;
69}

◆ SetVoxelsVisAttributes()

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

Definition at line 57 of file G4DrawVoxels.cc.

60{
61 fVoxelsVisAttributes[0] = VA_voxelX;
62 fVoxelsVisAttributes[1] = VA_voxelY;
63 fVoxelsVisAttributes[2] = VA_voxelZ;
64}

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