Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkPolydataPipeline.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25
27
28#include "G4Normal3D.hh"
29#include "G4Point3D.hh"
30#include "G4Polyhedron.hh"
31#include "G4Polyline.hh"
32#include "G4ViewParameters.hh"
33#include "G4VtkViewer.hh"
34#include "G4VtkVisContext.hh"
35
36#include <vtkActor.h>
37#include <vtkCleanPolyData.h>
38#include <vtkLine.h>
39#include <vtkMatrix4x4.h>
40#include <vtkPolyDataAlgorithm.h>
41#include <vtkPolyDataMapper.h>
42#include <vtkPolyDataNormals.h>
43#include <vtkProperty.h>
44#include <vtkTriangleFilter.h>
45
46std::size_t G4VtkPolydataPipeline::MakeHash(const G4Polyhedron& polyhedron,
47 const G4VtkVisContext& vc)
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;
69
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}
122
124{
125 actor->SetVisibility(1);
126}
127
129{
130 actor->SetVisibility(0);
131}
132
134{
135 G4cout << "G4VtkPolydataPipeline filters (";
136 for (const auto& f : filters)
137 G4cout << f->GetInformation() << ",";
138 G4cout << ")" << G4endl;
139
141}
142
144{
145 actor->Modified();
146 polydata->Modified();
147 mapper->Update();
148
150}
151
153{
154 renderer->RemoveActor(actor);
156}
157
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}
180
181void G4VtkPolydataPipeline::SetPolydata(vtkPolyData *polydataIn) {
182 polydata->DeepCopy(polydataIn);
183}
184
185
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}
202
204{
205 SetPolydataData(p.x(), p.y(), p.z());
206}
207
208void G4VtkPolydataPipeline::SetPolydataData(double x, double y, double z)
209{
210 polydataPoints->InsertNextPoint(x, y, z);
211}
212
214 G4double r01, G4double r02, G4double r10,
215 G4double r11, G4double r12, G4double r20,
216 G4double r21, G4double r22)
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}
225
227{
228 actor->GetProperty()->SetColor(r, g, b);
229 actor->GetProperty()->SetOpacity(a);
230}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void Clear()
virtual void Print()
void SetTypeName(G4String typeNameIn)
virtual void Modified()
vtkSmartPointer< vtkRenderer > renderer
G4VtkVisContext vc
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)
vtkSmartPointer< vtkPolyDataAlgorithm > GetFinalFilter()
static std::size_t MakeHash(const G4Polyhedron &p, const G4VtkVisContext &vc)
virtual void SetPolydata(const G4Polyhedron &polyhedron)
vtkSmartPointer< vtkCellArray > polydataCells
vtkSmartPointer< vtkPolyData > polydata
vtkSmartPointer< vtkPoints > polydataPoints
vtkSmartPointer< vtkActor > actor
vtkSmartPointer< vtkPolyDataMapper > mapper
virtual vtkSmartPointer< vtkActor > GetActor()
std::vector< vtkSmartPointer< vtkPolyDataAlgorithm > > filters
virtual void SetActorColour(G4double r, G4double g, G4double b, G4double a)
void AddFilter(vtkSmartPointer< vtkPolyDataAlgorithm > f)
virtual void SetPolydataData(const G4Point3D &p)
G4VtkPolydataPipeline(G4String name, const G4VtkVisContext &vc)
vtkNew< vtkRenderer > renderer
G4ViewParameters::DrawingStyle fDrawingStyle
const G4Transform3D & fTransform
const G4VtkViewer * fViewer
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
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const
void hash_combine(std::size_t)