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

#include <G4VtkPolydataInstanceAppendPipeline.hh>

+ Inheritance diagram for G4VtkPolydataInstanceAppendPipeline:

Public Member Functions

 G4VtkPolydataInstanceAppendPipeline (G4String name, const G4VtkVisContext &vc)
 
 ~G4VtkPolydataInstanceAppendPipeline () override=default
 
void Print () 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 G4Polyhedron &polyhedron)
 
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

std::map< G4String, vtkSmartPointer< vtkTransformPolyDataFilter > > transformFilterMap
 
vtkSmartPointer< vtkAppendPolyData > appendFilter
 
- 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 37 of file G4VtkPolydataInstanceAppendPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkPolydataInstanceAppendPipeline()

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

Definition at line 68 of file G4VtkPolydataInstanceAppendPipeline.cc.

70 : G4VtkPolydataInstancePipeline(nameIn, vcIn)
71{
72 // Set pipeline type
73 SetTypeName(G4String("G4VtkPolydataInstanceAppendPipeline"));
74
75 // append filter
78
79 // set polydata mapper
80 mapper->SetInputConnection(GetFinalFilter()->GetOutputPort());
81
82 // set actor
83 actor->SetMapper(mapper);
84 actor->SetVisibility(1);
85
86 // colour parameters
87 actor->GetProperty()->SetOpacity(vc.alpha);
88 actor->GetProperty()->SetColor(vc.red, vc.green, vc.blue);
89
90 // shading parameters
91 actor->GetProperty()->SetAmbient(0.2);
92 actor->GetProperty()->SetDiffuse(0.7);
93 actor->GetProperty()->SetSpecular(0.1);
94 actor->GetProperty()->SetSpecularPower(1);
95
96 // add to renderer
97 vc.fViewer->renderer->AddActor(GetActor());
98}
void SetTypeName(G4String typeNameIn)
G4VtkVisContext vc
vtkSmartPointer< vtkAppendPolyData > appendFilter
G4VtkPolydataInstancePipeline(G4String name, const G4VtkVisContext &vc)
vtkSmartPointer< vtkPolyDataAlgorithm > GetFinalFilter()
vtkSmartPointer< vtkActor > actor
vtkSmartPointer< vtkPolyDataMapper > mapper
virtual vtkSmartPointer< vtkActor > GetActor()
void AddFilter(vtkSmartPointer< vtkPolyDataAlgorithm > f)
vtkNew< vtkRenderer > renderer
const G4VtkViewer * fViewer

◆ ~G4VtkPolydataInstanceAppendPipeline()

G4VtkPolydataInstanceAppendPipeline::~G4VtkPolydataInstanceAppendPipeline ( )
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 G4VtkPolydataInstanceAppendPipeline::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 106 of file G4VtkPolydataInstanceAppendPipeline.cc.

112{
113 // Add to base class
114 G4VtkPolydataInstancePipeline::addInstance(dx, dy, dz, r00, r01, r02, r10, r11, r12, r20, r21,
115 r22, r, g, b, a, nameIn);
116
117 // create transform
119 double transformArray[16] = {r00, r01, r02, dx, r10, r11, r12, dy,
120 r20, r21, r22, dz, 0., 0., 0., 1.};
121 transform->Concatenate(transformArray);
122 transform->Update();
123
124 // Create transform filter and add to local filter map
127 tf->SetTransform(transform);
128 tf->SetInputConnection(GetFilter(GetNumberOfFilters() - 2)->GetOutputPort());
129 tf->Update();
131
132 // Add transform filter to append filter
133 GetFinalFilter()->AddInputConnection(tf->GetOutputPort());
134}
std::map< G4String, vtkSmartPointer< vtkTransformPolyDataFilter > > transformFilterMap
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< vtkPolyDataAlgorithm > GetFilter(G4int iFilter)

◆ MakeHash()

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

Definition at line 42 of file G4VtkPolydataInstanceAppendPipeline.cc.

44{
45 // Get view parameters that the user can force through the vis attributes, thereby over-riding the
46 // current view parameter.
47 const G4VisAttributes* pVA =
48 vc.fViewer->GetApplicableVisAttributes(polyhedron.GetVisAttributes());
49 G4Color colour = pVA->GetColour();
50
51 // Hash the vis attributes
52 std::size_t hash = std::hash<G4double>{}(colour.GetAlpha());
53 std::size_t rhash = std::hash<G4double>{}(colour.GetRed());
54 std::size_t ghash = std::hash<G4double>{}(colour.GetGreen());
55 std::size_t bhash = std::hash<G4double>{}(colour.GetBlue());
56 std::size_t phash = std::hash<G4Polyhedron>{}(polyhedron);
57 std::size_t shash = std::hash<G4double>{}(vc.fDrawingStyle);
58
59 std::hash_combine(hash, phash);
60 std::hash_combine(hash, rhash);
61 std::hash_combine(hash, bhash);
62 std::hash_combine(hash, ghash);
63 std::hash_combine(hash, shash);
64
65 return hash;
66}
G4double GetBlue() const
Definition G4Colour.hh:154
G4double GetAlpha() const
Definition G4Colour.hh:155
G4double GetRed() const
Definition G4Colour.hh:152
G4double GetGreen() const
Definition G4Colour.hh:153
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4Colour & GetColour() const
G4ViewParameters::DrawingStyle fDrawingStyle
void hash_combine(std::size_t)

Referenced by G4VtkStore::AddPrimitiveAppend().

◆ Print()

void G4VtkPolydataInstanceAppendPipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Definition at line 100 of file G4VtkPolydataInstanceAppendPipeline.cc.

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

◆ removeInstance()

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

Reimplemented from G4VtkPolydataInstancePipeline.

Definition at line 136 of file G4VtkPolydataInstanceAppendPipeline.cc.

137{
138 // Remove from base class
140
141 // Remove transform filter from append filter
142 appendFilter->RemoveInputConnection(0, transformFilterMap[name]->GetOutputPort());
143
144 // Remove from local filter map
145 transformFilterMap.erase(nameIn);
146}
virtual void removeInstance(const G4String &name)

Member Data Documentation

◆ appendFilter

vtkSmartPointer<vtkAppendPolyData> G4VtkPolydataInstanceAppendPipeline::appendFilter
protected

◆ transformFilterMap

std::map<G4String, vtkSmartPointer<vtkTransformPolyDataFilter> > G4VtkPolydataInstanceAppendPipeline::transformFilterMap
protected

Definition at line 55 of file G4VtkPolydataInstanceAppendPipeline.hh.

Referenced by addInstance(), and removeInstance().


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