Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/visualization/Vtk/src/G4VtkPolydataPipeline.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 26 #include "G4VtkPolydataPipeline.hh"
 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 
 46 std::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;
 68 }
 69 
 70 G4VtkPolydataPipeline::G4VtkPolydataPipeline(G4String nameIn, const G4VtkVisContext& vcIn)
 71   : G4VVtkPipeline(nameIn, "G4VtkPolydataPipeline", vcIn, false, vcIn.fViewer->renderer)
 72 {
 73   // Set pipeline type
 74   SetTypeName(G4String("G4VtkPolydataPipeline"));
 75 
 76   polydataPoints = vtkSmartPointer<vtkPoints>::New();
 77   polydataCells = vtkSmartPointer<vtkCellArray>::New();
 78   polydata = vtkSmartPointer<vtkPolyData>::New();
 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
101   mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
102   mapper->SetInputConnection(GetFinalFilter()->GetOutputPort());
103   mapper->SetColorModeToDirectScalars();
104 
105   // add to actor
106   actor = vtkSmartPointer<vtkActor>::New();
107   actor->SetMapper(mapper);
108   actor->SetVisibility(1);
109 
110   // set actor properties from vis context
111   if (vc.fDrawingStyle == G4ViewParameters::hsr) {
112   }
113   else if (vc.fDrawingStyle == G4ViewParameters::hlr) {
114   }
115   else if (vc.fDrawingStyle == G4ViewParameters::wireframe) {
116     actor->GetProperty()->SetRepresentationToWireframe();
117   }
118 
119   // add to renderer
120   vc.fViewer->renderer->AddActor(GetActor());
121 }
122 
123 void G4VtkPolydataPipeline::Enable()
124 {
125   actor->SetVisibility(1);
126 }
127 
128 void G4VtkPolydataPipeline::Disable()
129 {
130   actor->SetVisibility(0);
131 }
132 
133 void G4VtkPolydataPipeline::Print()
134 {
135   G4cout << "G4VtkPolydataPipeline filters (";
136   for (const auto& f : filters)
137     G4cout << f->GetInformation() << ",";
138   G4cout << ")" << G4endl;
139 
140   G4VVtkPipeline::Print();
141 }
142 
143 void G4VtkPolydataPipeline::Modified()
144 {
145   actor->Modified();
146   polydata->Modified();
147   mapper->Update();
148 
149   G4VVtkPipeline::Modified();
150 }
151 
152 void G4VtkPolydataPipeline::Clear()
153 {
154   renderer->RemoveActor(actor);
155   G4VVtkPipeline::Clear();
156 }
157 
158 void G4VtkPolydataPipeline::SetPolydata(const G4Polyhedron& polyhedron)
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 
169     vtkSmartPointer<vtkIdList> poly = vtkSmartPointer<vtkIdList>::New();
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 
181 void G4VtkPolydataPipeline::SetPolydata(vtkPolyData *polydataIn) {
182   polydata->DeepCopy(polydataIn);
183 }
184 
185 
186 void G4VtkPolydataPipeline::SetPolydata(const G4Polyline& polyline)
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) {
195       vtkSmartPointer<vtkLine> line = vtkSmartPointer<vtkLine>::New();
196       line->GetPointIds()->SetId(0, id);
197       line->GetPointIds()->SetId(1, id + 1);
198       polydataCells->InsertNextCell(line);
199     }
200   }
201 }
202 
203 void G4VtkPolydataPipeline::SetPolydataData(const G4Point3D& p)
204 {
205   SetPolydataData(p.x(), p.y(), p.z());
206 }
207 
208 void G4VtkPolydataPipeline::SetPolydataData(double x, double y, double z)
209 {
210   polydataPoints->InsertNextPoint(x, y, z);
211 }
212 
213 void G4VtkPolydataPipeline::SetActorTransform(G4double dx, G4double dy, G4double dz, G4double r00,
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 
226 void G4VtkPolydataPipeline::SetActorColour(G4double r, G4double g, G4double b, G4double a)
227 {
228   actor->GetProperty()->SetColor(r, g, b);
229   actor->GetProperty()->SetOpacity(a);
230 }
231