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

#include <G4VtkPolydataPipeline.hh>

+ Inheritance diagram for G4VtkPolydataPipeline:

Public Member Functions

 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)
 

Protected Attributes

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 48 of file G4VtkPolydataPipeline.hh.

Constructor & Destructor Documentation

◆ G4VtkPolydataPipeline()

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

Definition at line 70 of file G4VtkPolydataPipeline.cc.

71 : G4VVtkPipeline(nameIn, "G4VtkPolydataPipeline", vcIn, false, vcIn.fViewer->renderer)
72{
73 // Set pipeline type
74 SetTypeName(G4String("G4VtkPolydataPipeline"));
75
79
80 polydata->SetPoints(polydataPoints);
81 polydata->SetPolys(polydataCells);
82
83 // clean input polydata
84 auto filterClean = vtkSmartPointer<vtkCleanPolyData>::New();
85 filterClean->PointMergingOn();
86 filterClean->AddInputData(polydata);
87 AddFilter(filterClean);
88
89 // ensure triangular mesh
90 auto filterTriangle = vtkSmartPointer<vtkTriangleFilter>::New();
91 filterTriangle->SetInputConnection(filterClean->GetOutputPort());
92 AddFilter(filterTriangle);
93
94 // calculate normals with a feature angle of 45 degrees
95 auto filterNormals = vtkSmartPointer<vtkPolyDataNormals>::New();
96 filterNormals->SetFeatureAngle(45);
97 filterNormals->SetInputConnection(filterTriangle->GetOutputPort());
98 AddFilter(filterNormals);
99
100 // mapper
102 mapper->SetInputConnection(GetFinalFilter()->GetOutputPort());
103 mapper->SetColorModeToDirectScalars();
104
105 // add to actor
107 actor->SetMapper(mapper);
108 actor->SetVisibility(1);
109
110 // set actor properties from vis context
112 }
114 }
116 actor->GetProperty()->SetRepresentationToWireframe();
117 }
118
119 // add to renderer
120 vc.fViewer->renderer->AddActor(GetActor());
121}
void SetTypeName(G4String typeNameIn)
G4VtkVisContext vc
vtkSmartPointer< vtkPolyDataAlgorithm > GetFinalFilter()
vtkSmartPointer< vtkCellArray > polydataCells
vtkSmartPointer< vtkPolyData > polydata
vtkSmartPointer< vtkPoints > polydataPoints
vtkSmartPointer< vtkActor > actor
vtkSmartPointer< vtkPolyDataMapper > mapper
virtual vtkSmartPointer< vtkActor > GetActor()
void AddFilter(vtkSmartPointer< vtkPolyDataAlgorithm > f)
vtkNew< vtkRenderer > renderer
G4ViewParameters::DrawingStyle fDrawingStyle
const G4VtkViewer * fViewer

◆ ~G4VtkPolydataPipeline()

G4VtkPolydataPipeline::~G4VtkPolydataPipeline ( )
overridedefault

Member Function Documentation

◆ AddFilter()

◆ Clear()

void G4VtkPolydataPipeline::Clear ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Reimplemented in G4VtkPolydataPolyline2DPipeline.

Definition at line 152 of file G4VtkPolydataPipeline.cc.

153{
154 renderer->RemoveActor(actor);
156}
virtual void Clear()
vtkSmartPointer< vtkRenderer > renderer

Referenced by G4VtkPolydataPolyline2DPipeline::Clear().

◆ Disable()

void G4VtkPolydataPipeline::Disable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 128 of file G4VtkPolydataPipeline.cc.

129{
130 actor->SetVisibility(0);
131}

◆ Enable()

void G4VtkPolydataPipeline::Enable ( )
overridevirtual

Implements G4VVtkPipeline.

Definition at line 123 of file G4VtkPolydataPipeline.cc.

124{
125 actor->SetVisibility(1);
126}

◆ GetActor()

◆ GetBounds()

virtual G4double * G4VtkPolydataPipeline::GetBounds ( )
inlinevirtual

Definition at line 81 of file G4VtkPolydataPipeline.hh.

81{ return actor->GetBounds(); }

◆ GetFilter()

vtkSmartPointer< vtkPolyDataAlgorithm > G4VtkPolydataPipeline::GetFilter ( G4int iFilter)
inline

Definition at line 55 of file G4VtkPolydataPipeline.hh.

55{ return filters[iFilter]; }

Referenced by G4VtkPolydataInstanceAppendPipeline::addInstance().

◆ GetFinalFilter()

◆ GetNumberOfFilters()

G4int G4VtkPolydataPipeline::GetNumberOfFilters ( )
inline

Definition at line 56 of file G4VtkPolydataPipeline.hh.

56{ return (G4int)filters.size(); }
int G4int
Definition G4Types.hh:85

Referenced by G4VtkPolydataInstanceAppendPipeline::addInstance().

◆ GetPolydata()

virtual vtkSmartPointer< vtkPolyData > G4VtkPolydataPipeline::GetPolydata ( )
inlinevirtual

Definition at line 79 of file G4VtkPolydataPipeline.hh.

79{ return polydata; }

◆ MakeHash()

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

Definition at line 46 of file G4VtkPolydataPipeline.cc.

48{
49 std::size_t hash = std::hash<G4Polyhedron>{}(polyhedron);
50
51 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.dx()));
52 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.dy()));
53 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.dz()));
54
55 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.xx()));
56 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.xy()));
57 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.xz()));
58
59 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.yx()));
60 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.yy()));
61 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.yz()));
62
63 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.zx()));
64 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.zy()));
65 std::hash_combine(hash, std::hash<G4double>{}(vc.fTransform.zz()));
66
67 return hash;
68}
const G4Transform3D & fTransform
double dy() const
double zz() const
double yz() const
double dz() const
double dx() const
double xy() const
double zx() const
double yx() const
double zy() const
double xx() const
double yy() const
double xz() const
void hash_combine(std::size_t)

Referenced by G4VtkStore::AddPrimitiveSeparate().

◆ Modified()

void G4VtkPolydataPipeline::Modified ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Reimplemented in G4VtkPolydataPolyline2DPipeline.

Definition at line 143 of file G4VtkPolydataPipeline.cc.

144{
145 actor->Modified();
146 polydata->Modified();
147 mapper->Update();
148
150}
virtual void Modified()

Referenced by G4VtkPolydataPolyline2DPipeline::Modified().

◆ Print()

void G4VtkPolydataPipeline::Print ( )
overridevirtual

Reimplemented from G4VVtkPipeline.

Reimplemented in G4VtkPolydataPolylinePipeline, and G4VtkPolydataSpherePipeline.

Definition at line 133 of file G4VtkPolydataPipeline.cc.

134{
135 G4cout << "G4VtkPolydataPipeline filters (";
136 for (const auto& f : filters)
137 G4cout << f->GetInformation() << ",";
138 G4cout << ")" << G4endl;
139
141}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void Print()

Referenced by G4VtkPolydataInstanceAppendPipeline::Print(), G4VtkPolydataInstanceBakePipeline::Print(), G4VtkPolydataInstanceTensorPipeline::Print(), G4VtkPolydataPolylinePipeline::Print(), and G4VtkPolydataSpherePipeline::Print().

◆ SetActorColour()

void G4VtkPolydataPipeline::SetActorColour ( G4double r,
G4double g,
G4double b,
G4double a )
virtual

Definition at line 226 of file G4VtkPolydataPipeline.cc.

227{
228 actor->GetProperty()->SetColor(r, g, b);
229 actor->GetProperty()->SetOpacity(a);
230}

◆ SetActorTransform()

void G4VtkPolydataPipeline::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

Definition at line 213 of file G4VtkPolydataPipeline.cc.

217{
218 // create transform
219 auto transform = vtkSmartPointer<vtkMatrix4x4>::New();
220 double transformArray[16] = {r00, r01, r02, dx, r10, r11, r12, dy,
221 r20, r21, r22, dz, 0., 0., 0., 1.};
222 transform->DeepCopy(transformArray);
223 actor->SetUserMatrix(transform);
224}

◆ SetPolydata() [1/3]

void G4VtkPolydataPipeline::SetPolydata ( const G4Polyhedron & polyhedron)
virtual

Reimplemented in G4VtkPolydataInstanceBakePipeline, and G4VtkPolydataPolyline2DPipeline.

Definition at line 158 of file G4VtkPolydataPipeline.cc.

159{
160 G4bool notLastFace;
161 int iVert = 0;
162 do {
163 G4Point3D vertex[4];
164 G4int edgeFlag[4];
165 G4Normal3D normals[4];
166 G4int nEdges;
167 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals);
168
170 // loop over vertices
171 for (int i = 0; i < nEdges; i++) {
172 polydataPoints->InsertNextPoint(vertex[i].x(), vertex[i].y(), vertex[i].z());
173 poly->InsertNextId(iVert);
174 iVert++;
175 }
176 polydataCells->InsertNextCell(poly);
177
178 } while (notLastFace);
179}
bool G4bool
Definition G4Types.hh:86
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const

◆ SetPolydata() [2/3]

void G4VtkPolydataPipeline::SetPolydata ( const G4Polyline & polyline)
virtual

Reimplemented in G4VtkPolydataPolyline2DPipeline.

Definition at line 186 of file G4VtkPolydataPipeline.cc.

187{
188 // Data data
189 const size_t nLines = polyline.size();
190
191 for (size_t i = 0; i < nLines; ++i) {
192 auto id = polydataPoints->InsertNextPoint(polyline[i].x(), polyline[i].y(), polyline[i].z());
193
194 if (i < nLines - 1) {
196 line->GetPointIds()->SetId(0, id);
197 line->GetPointIds()->SetId(1, id + 1);
198 polydataCells->InsertNextCell(line);
199 }
200 }
201}

◆ SetPolydata() [3/3]

void G4VtkPolydataPipeline::SetPolydata ( vtkPolyData * polydata)
virtual

Reimplemented in G4VtkPolydataPolyline2DPipeline.

Definition at line 181 of file G4VtkPolydataPipeline.cc.

181 {
182 polydata->DeepCopy(polydataIn);
183}

◆ SetPolydataData() [1/2]

void G4VtkPolydataPipeline::SetPolydataData ( const G4Point3D & p)
virtual

Definition at line 203 of file G4VtkPolydataPipeline.cc.

204{
205 SetPolydataData(p.x(), p.y(), p.z());
206}
virtual void SetPolydataData(const G4Point3D &p)

Referenced by SetPolydataData().

◆ SetPolydataData() [2/2]

void G4VtkPolydataPipeline::SetPolydataData ( double x,
double y,
double z )
virtual

Definition at line 208 of file G4VtkPolydataPipeline.cc.

209{
210 polydataPoints->InsertNextPoint(x, y, z);
211}

Member Data Documentation

◆ actor

◆ filters

std::vector<vtkSmartPointer<vtkPolyDataAlgorithm> > G4VtkPolydataPipeline::filters
protected

◆ mapper

◆ polydata

◆ polydataCells

◆ polydataPoints


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