Geant4 11.1.1
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 ()
 
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}
void SetColour(const G4Colour &)

◆ ~G4DrawVoxels()

G4DrawVoxels::~G4DrawVoxels ( )

Definition at line 57 of file G4DrawVoxels.cc.

58{
59}

Member Function Documentation

◆ 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

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:59
G4ThreeVector NetTranslation() const
G4RotationMatrix NetRotation() const
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
std::size_t GetNoDaughters() 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 72 of file G4DrawVoxels.cc.

73{
74 fBoundingBoxVisAttributes = VA_boundingbox;
75}

◆ 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}

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