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

#include <G4VtkPolydataInstanceBakePipeline.hh>

+ Inheritance diagram for G4VtkPolydataInstanceBakePipeline:

Public Member Functions

 G4VtkPolydataInstanceBakePipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkPolydataInstanceBakePipeline () override=default
 
void Print () override
 
void SetPolydata (const G4Polyhedron &polyhedron) override
 
void addInstance (G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22, G4double r, G4double g, G4double b, G4double a, const G4String &name) override
 
void removeInstance (const G4String &name) override
 
virtual void addInstance (G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22, const G4String &name)
 
- Public Member Functions inherited from G4VtkPolydataInstancePipeline
 G4VtkPolydataInstancePipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkPolydataInstancePipeline () override=default
 
- Public Member Functions inherited from G4VtkPolydataPipeline
 G4VtkPolydataPipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkPolydataPipeline () override=default
 
void AddFilter (vtkSmartPointer< vtkPolyDataAlgorithm > f)
 
vtkSmartPointer< vtkPolyDataAlgorithm > GetFilter (G4int iFilter)
 
G4int GetNumberOfFilters ()
 
vtkSmartPointer< vtkPolyDataAlgorithm > GetFinalFilter ()
 
void Enable () override
 
void Disable () override
 
void Print () override
 
void Modified () override
 
void Clear () override
 
virtual vtkSmartPointer< vtkActor > GetActor ()
 
virtual void SetPolydata (const G4Polyline &polyline)
 
virtual void SetPolydata (vtkPolyData *polydata)
 
virtual void SetPolydataData (const G4Point3D &p)
 
virtual void SetPolydataData (double x, double y, double z)
 
virtual void SetActorTransform (G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22)
 
virtual void SetActorColour (G4double r, G4double g, G4double b, G4double a)
 
virtual vtkSmartPointer< vtkPolyData > GetPolydata ()
 
virtual G4doubleGetBounds ()
 
- Public Member Functions inherited from G4VVtkPipeline
 G4VVtkPipeline ()
 
 G4VVtkPipeline (G4String nameIn, G4String typeIn, const G4VtkVisContext &vcIn, G4bool disableParentIn=false, vtkSmartPointer< vtkRenderer > rendererIn=nullptr)
 
virtual ~G4VVtkPipeline ()
 
void SetDisableParent (G4bool disableParentIn)
 
bool GetDisableParent ()
 
void SetName (G4String nameIn)
 
const G4String GetName ()
 
void SetTypeName (G4String typeNameIn)
 
G4String GetTypeName ()
 
G4VtkVisContextGetVtkVisContext ()
 
void AddChildPipeline (G4VVtkPipeline *child)
 
G4VVtkPipelineGetChildPipeline (G4int iPipeline)
 
G4int GetNumberOfChildPipelines ()
 
std::vector< G4VVtkPipeline * > GetChildPipelines ()
 
void ClearChildPipeline ()
 

Static Public Member Functions

static std::size_t MakeHash (const G4Polyhedron &p, const G4VtkVisContext &vc)
 
- Static Public Member Functions inherited from G4VtkPolydataPipeline
static std::size_t MakeHash (const G4Polyhedron &p, const G4VtkVisContext &vc)
 

Protected Attributes

vtkSmartPointer< vtkDoubleArray > polydataPointData
 
const G4PolyhedronpolyhedronPrototype
 
G4int iVert
 
G4int iFace
 
std::map< G4String, std::pair< vtkIdType, vtkIdType > > instanceVertexMap
 
std::map< G4String, std::pair< vtkIdType, vtkIdType > > instanceFaceMap
 
- Protected Attributes inherited from G4VtkPolydataInstancePipeline
vtkSmartPointer< vtkDoubleArray > instanceColour
 
vtkSmartPointer< vtkPoints > instancePosition
 
vtkSmartPointer< vtkDoubleArray > instanceTransform
 
std::map< G4String, vtkIdType > instanceMap
 
- Protected Attributes inherited from G4VtkPolydataPipeline
vtkSmartPointer< vtkPoints > polydataPoints
 
vtkSmartPointer< vtkCellArray > polydataCells
 
vtkSmartPointer< vtkPolyData > polydata
 
std::vector< vtkSmartPointer< vtkPolyDataAlgorithm > > filters
 
vtkSmartPointer< vtkPolyDataMapper > mapper
 
vtkSmartPointer< vtkActor > actor
 
- Protected Attributes inherited from G4VVtkPipeline
G4String name
 
G4String type
 
G4bool disableParent
 
std::vector< G4VVtkPipeline * > childPipelines
 
vtkSmartPointer< vtkRenderer > renderer
 
G4VtkVisContext vc
 

Detailed Description

Definition at line 34 of file G4VtkPolydataInstanceBakePipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkPolydataInstanceBakePipeline()

G4VtkPolydataInstanceBakePipeline::G4VtkPolydataInstanceBakePipeline ( G4String name,
const G4VtkVisContext & vc )

Definition at line 57 of file G4VtkPolydataInstanceBakePipeline.cc.

59 : G4VtkPolydataInstancePipeline(nameIn, vcIn)
60{
61 iVert = 0;
62 iFace = 0;
64 polydataPointData->SetNumberOfComponents(3);
65 polydataPointData->SetName("Colors");
66 polydata->GetPointData()->SetScalars(polydataPointData);
67
68 mapper->SetColorModeToDirectScalars();
69
70 actor->GetProperty()->SetOpacity(vc.alpha);
71}
G4VtkVisContext vc
vtkSmartPointer< vtkDoubleArray > polydataPointData
G4VtkPolydataInstancePipeline(G4String name, const G4VtkVisContext &vc)
vtkSmartPointer< vtkPolyData > polydata
vtkSmartPointer< vtkActor > actor
vtkSmartPointer< vtkPolyDataMapper > mapper

◆ ~G4VtkPolydataInstanceBakePipeline()

G4VtkPolydataInstanceBakePipeline::~G4VtkPolydataInstanceBakePipeline ( )
overridedefault

Member Function Documentation

◆ addInstance() [1/2]

void G4VtkPolydataInstancePipeline::addInstance ( G4double dx,
G4double dy,
G4double dz,
G4double r00,
G4double r01,
G4double r02,
G4double r10,
G4double r11,
G4double r12,
G4double r20,
G4double r21,
G4double r22,
const G4String & name )
virtual

Reimplemented from G4VtkPolydataInstancePipeline.

Definition at line 42 of file G4VtkPolydataInstancePipeline.cc.

53{
54 // add the instance without colour or alpha
55 addInstance(dx, dy, dz, r00, r01, r02, r10, r11, r12, r20, r21, r22, 0, 0, 0, 0, nameIn);
56}
void addInstance(G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22, G4double r, G4double g, G4double b, G4double a, const G4String &name) override

◆ addInstance() [2/2]

void G4VtkPolydataInstanceBakePipeline::addInstance ( G4double dx,
G4double dy,
G4double dz,
G4double r00,
G4double r01,
G4double r02,
G4double r10,
G4double r11,
G4double r12,
G4double r20,
G4double r21,
G4double r22,
G4double r,
G4double g,
G4double b,
G4double a,
const G4String & name )
overridevirtual

Reimplemented from G4VtkPolydataInstancePipeline.

Definition at line 84 of file G4VtkPolydataInstanceBakePipeline.cc.

90{
91 G4VtkPolydataInstancePipeline::addInstance(dx, dy, dz, r00, r01, r02, r10, r11, r12, r20, r21,
92 r22, r, g, b, a, nameIn);
93
94 vtkIdType vStart;
95 vtkIdType vEnd;
96 vtkIdType fStart;
97 vtkIdType fEnd;
98
99 vStart = iVert;
100 fStart = iFace;
101
102 G4bool notLastFace;
103 do {
104 G4Point3D vertex[4];
105 G4int edgeFlag[4];
106 G4Normal3D normals[4];
107 G4int nEdges;
108 notLastFace = polyhedronPrototype->GetNextFacet(nEdges, vertex, edgeFlag, normals);
109
111 // loop over vertices
112 for (int i = 0; i < nEdges; i++) {
113 // note : G4Transform3D does not have a matrix element constructor
114 G4Point3D bakedVertex =
115 G4Point3D(vertex[i].x() * r00 + vertex[i].y() * r01 + vertex[i].z() * r02 + dx,
116 vertex[i].x() * r10 + vertex[i].y() * r11 + vertex[i].z() * r12 + dy,
117 vertex[i].x() * r20 + vertex[i].y() * r21 + vertex[i].z() * r22 + dz);
118 polydataPoints->InsertNextPoint(bakedVertex.x(), bakedVertex.y(), bakedVertex.z());
119 poly->InsertNextId(iVert);
120 iVert++;
121 double cols[] = {r, g, b, 0.5};
122 polydataPointData->InsertNextTuple(cols);
123 }
124
125 polydataCells->InsertNextCell(poly);
126 iFace++;
127
128 } while (notLastFace);
129
130 vEnd = iVert;
131 fEnd = iFace;
132
133 // Add vertex and face to maps to allow for removal later
134 instanceVertexMap[name] = std::pair<vtkIdType, vtkIdType>(vStart, vEnd);
135 instanceFaceMap[name] = std::pair<vtkIdType, vtkIdType>(fStart, fEnd);
136}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
std::map< G4String, std::pair< vtkIdType, vtkIdType > > instanceVertexMap
std::map< G4String, std::pair< vtkIdType, vtkIdType > > instanceFaceMap
virtual void addInstance(G4double dx, G4double dy, G4double dz, G4double r00, G4double r01, G4double r02, G4double r10, G4double r11, G4double r12, G4double r20, G4double r21, G4double r22, const G4String &name)
vtkSmartPointer< vtkCellArray > polydataCells
vtkSmartPointer< vtkPoints > polydataPoints
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const

◆ MakeHash()

std::size_t G4VtkPolydataInstanceBakePipeline::MakeHash ( const G4Polyhedron & p,
const G4VtkVisContext & vc )
static

Definition at line 39 of file G4VtkPolydataInstanceBakePipeline.cc.

41{
42 // Get view parameters that the user can force through the vis attributes, thereby over-riding the
43 // current view parameter.
44 const G4VisAttributes* pVA =
45 vc.fViewer->GetApplicableVisAttributes(polyhedron.GetVisAttributes());
46 G4Colour colour = pVA->GetColour();
47
48 // Hash the vis attributes
49 std::size_t hash = std::hash<G4double>{}(colour.GetAlpha());
50 std::size_t shash = std::hash<G4double>{}(vc.fDrawingStyle);
51
52 std::hash_combine(hash, shash);
53
54 return hash;
55}
G4double GetAlpha() const
Definition G4Colour.hh:155
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4Colour & GetColour() const
G4ViewParameters::DrawingStyle fDrawingStyle
const G4VtkViewer * fViewer
void hash_combine(std::size_t)

Referenced by G4VtkStore::AddPrimitiveTransformBake().

◆ Print()

void G4VtkPolydataInstanceBakePipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 78 of file G4VtkPolydataInstanceBakePipeline.cc.

79{
80 G4cout << "G4VtkPolydataInstanceBakePipeline " << GetName() << G4endl;
82}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String GetName()

◆ removeInstance()

void G4VtkPolydataInstanceBakePipeline::removeInstance ( const G4String & name)
overridevirtual

Reimplemented from G4VtkPolydataInstancePipeline.

Definition at line 138 of file G4VtkPolydataInstanceBakePipeline.cc.

138{}

◆ SetPolydata()

void G4VtkPolydataInstanceBakePipeline::SetPolydata ( const G4Polyhedron & polyhedron)
overridevirtual

Reimplemented from G4VtkPolydataPipeline.

Definition at line 73 of file G4VtkPolydataInstanceBakePipeline.cc.

74{
75 polyhedronPrototype = &polyhedron;
76}

Member Data Documentation

◆ iFace

G4int G4VtkPolydataInstanceBakePipeline::iFace
protected

◆ instanceFaceMap

std::map<G4String, std::pair<vtkIdType, vtkIdType> > G4VtkPolydataInstanceBakePipeline::instanceFaceMap
protected

Definition at line 66 of file G4VtkPolydataInstanceBakePipeline.hh.

Referenced by addInstance().

◆ instanceVertexMap

std::map<G4String, std::pair<vtkIdType, vtkIdType> > G4VtkPolydataInstanceBakePipeline::instanceVertexMap
protected

Definition at line 65 of file G4VtkPolydataInstanceBakePipeline.hh.

Referenced by addInstance().

◆ iVert

G4int G4VtkPolydataInstanceBakePipeline::iVert
protected

◆ polydataPointData

vtkSmartPointer<vtkDoubleArray> G4VtkPolydataInstanceBakePipeline::polydataPointData
protected

◆ polyhedronPrototype

const G4Polyhedron* G4VtkPolydataInstanceBakePipeline::polyhedronPrototype
protected

Definition at line 58 of file G4VtkPolydataInstanceBakePipeline.hh.

Referenced by addInstance(), and SetPolydata().


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