Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VtkUnstructuredGridPipeline.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
26#include <vtkPoints.h>
27#include <vtkCellArray.h>
28#include <vtkCharArray.h>
29#include <vtkUnstructuredGrid.h>
30#include <vtkDataSetMapper.h>
31#include <vtkUnstructuredGridVolumeMapper.h>
32#include <vtkUnstructuredGridVolumeRayCastMapper.h>
33#include <vtkUnstructuredGridVolumeZSweepMapper.h>
34#include <vtkOpenGLProjectedTetrahedraMapper.h>
35#include <vtkActor.h>
36#include <vtkColorTransferFunction.h>
37#include <vtkPiecewiseFunction.h>
38#include <vtkLookupTable.h>
39#include <vtkVolumeProperty.h>
40#include <vtkVolume.h>
41#include <vtkNew.h>
42#include <vtkTetra.h>
43#include <vtkVoxel.h>
44// #include <vtkCleanUnstructuredGrid.h>
45#include <vtkDataSetTriangleFilter.h>
46#include <vtkClipDataSet.h>
47
49#include "G4VtkVisContext.hh"
50#include "G4VtkViewer.hh"
51
52#include "G4Mesh.hh"
53#include "G4LogicalVolume.hh"
57#include "G4PseudoScene.hh"
58#include "G4Transform3D.hh"
59#include "G4VSolid.hh"
60#include "G4Tet.hh"
61#include "G4Material.hh"
62
64 G4VVtkPipeline(nameIn, "G4VtkUnstructuredPipeline", vcIn, false, vcIn.fViewer->renderer) {
65
67
68 pointColourValues = vtkSmartPointer<vtkDoubleArray>::New();
69 cellColourValues = vtkSmartPointer<vtkDoubleArray>::New();
70 pointColourValues->SetNumberOfComponents(4);
71 cellColourValues->SetNumberOfComponents(4);
72
73 pointColourIndices = vtkSmartPointer<vtkDoubleArray>::New();
74 cellColourIndices = vtkSmartPointer<vtkDoubleArray>::New();
75 pointColourIndices->SetNumberOfComponents(1);
76 cellColourIndices->SetNumberOfComponents(1);
77
79 colourLUT->DiscretizeOn();
80
82 unstructuredGrid->SetPoints(points);
83 unstructuredGrid->GetPointData()->SetScalars(pointColourValues);
84 unstructuredGrid->GetCellData()->SetScalars(cellColourValues);
85
86 // clean filter
87#if 0
89 clean->SetInputData(unstructuredGrid);
90 clean->ToleranceIsAbsoluteOff();
91 clean->SetTolerance(1e-6);
92#endif
93
94 // Clip filter
96 vtkNew<vtkPlane> plane;
97 clip->SetClipFunction(plane);
98 clip->SetInputData(unstructuredGrid);
99
100 // Triangle filter
102 tri->SetInputData(unstructuredGrid);
103
104 // create dataset mapper
106 mapper->SetScalarModeToUseCellData();
107 mapper->SetColorModeToDirectScalars();
108 mapper->SetInputData(unstructuredGrid);
109 //mapper->SetInputConnection(clip->GetOutputPort());
110 //mapper->SetInputConnection(tri->GetOutputPort());
111
112 // create volume mapper
114 volumeMapper->SetScalarModeToUseCellData();
115 volumeMapper->SetInputConnection(tri->GetOutputPort());
116
117 // create volume properties
118 auto volumeProp = vtkSmartPointer<vtkVolumeProperty>::New();
119 volumeProp->SetColor(colourLUT);
120
121 // create actor
123 actor->SetMapper(mapper);
124 actor->SetVisibility(1);
125
126 // create volume
128 volume->SetMapper(volumeMapper);
129 //volume->SetProperty(volumeProp);
130 volume->SetVisibility(1);
131
132 // add to renderer
133 vc.fViewer->renderer->AddActor(actor);
134 //vc.fViewer->renderer->AddVolume(volume);
135}
136
137void G4VtkUnstructuredGridPipeline::PseudoSceneForTetCells::AddSolid(const G4VSolid& solid) {
138 if (fpPVModel->GetCurrentDepth() == fDepth) { // Leaf-level cells only
139 // Need to know it's a tet !!!! or implement G4VSceneHandler::AddSolid (const G4Tet&) !!!!
140 try {
141 const G4Tet& tet = dynamic_cast<const G4Tet&>(solid);
142 const auto* pVisAtts = fpPVModel->GetCurrentLV()->GetVisAttributes();
143 const auto& colour = pVisAtts->GetColour();
144
145 G4ThreeVector p0, p1, p2, p3;
146 tet.GetVertices(p0, p1, p2, p3);
147
148 G4int idx[4] = {-1, -1, -1, -1};
149
150 if (iCell > 0 && iCell < 100000000) {
151#if 0
152 G4ThreeVector pts[4] = {p0, p1, p2, p3};
153 for (G4int i = 0; i < 4; i++) {
154 auto h = MakeHash(pts[i]);
155 if (fpPointMap.find(h) == fpPointMap.end()) {
156 fpPointMap.insert(std::make_pair(h, iPoint));
157 fpPoints->InsertNextPoint(pts[i].x(), pts[i].y(), pts[i].z());
158 fpPointVector.push_back(pts[i]);
159 idx[i] = iPoint;
160 iPoint++;
161 }
162 else {
163 idx[i] = fpPointMap[h];
164 }
165 }
166#endif
167
168#if 1
169 fpPoints->InsertNextPoint(p0.x(), p0.y(), p0.z());
170 fpPoints->InsertNextPoint(p1.x(), p1.y(), p1.z());
171 fpPoints->InsertNextPoint(p2.x(), p2.y(), p2.z());
172 fpPoints->InsertNextPoint(p3.x(), p3.y(), p3.z());
173 iPoint += 4;
174
175 idx[0] = 4 * iCellAdd;
176 idx[1] = 4 * iCellAdd + 1;
177 idx[2] = 4 * iCellAdd + 2;
178 idx[3] = 4 * iCellAdd + 3;
179#endif
180
181 vtkNew<vtkTetra> tetra;
182 tetra->GetPointIds()->SetId(0, idx[0]);
183 tetra->GetPointIds()->SetId(1, idx[1]);
184 tetra->GetPointIds()->SetId(2, idx[2]);
185 tetra->GetPointIds()->SetId(3, idx[3]);
186
187 fpGrid->InsertNextCell(tetra->GetCellType(), tetra->GetPointIds());
188
189 auto hash = MakeHash(colour);
190 unsigned short int iColour = 0;
191 if (fpColourMap.find(hash) == fpColourMap.end()) {
192 iColour = fpColourMap.size();
193 fpColourMap.insert(std::make_pair(hash,iColour));
194 double c[4] = {colour.GetRed(), colour.GetGreen(), colour.GetBlue(), colour.GetAlpha()};
195 fpColourLUT->SetIndexedColorRGBA(iColour,c);
196 }
197 else {
198 iColour = fpColourMap[hash];
199 }
200
201 fpCellColourValues->InsertNextTuple4(colour.GetRed(), colour.GetGreen(), colour.GetBlue(), colour.GetAlpha());
202 fpCellColourIndices->InsertNextTuple1(iColour);
203
204 iCellAdd++;
205 }
206 iCell++;
207 }
208 catch (const std::bad_cast&) {
210 ed << "Called for a mesh that is not a tetrahedron mesh: " << solid.GetName();
211 G4Exception("PseudoSceneForTetVertices","visman0108",JustWarning,ed);
212 }
213 }
214}
215
216void G4VtkUnstructuredGridPipeline::PseudoSceneForCubicalCells::AddSolid(const G4Box& box) {
217 if (fpPVModel->GetCurrentDepth() == fDepth) { // Leaf-level cells only
218 try {
219 const auto* pVisAtts = fpPVModel->GetCurrentLV()->GetVisAttributes();
220 const auto& colour = pVisAtts->GetColour();
221 const G4ThreeVector& position = fpCurrentObjectTransformation->getTranslation();
222
223 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()-box.GetZHalfLength());
224 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()-box.GetZHalfLength());
225 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()-box.GetZHalfLength());
226 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()-box.GetZHalfLength());
227
228 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()+box.GetZHalfLength());
229 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()+box.GetZHalfLength());
230 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()+box.GetZHalfLength());
231 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()+box.GetZHalfLength());
232
233 vtkNew<vtkVoxel> voxel;
234 voxel->GetPointIds()->SetId(0, 8*iCell);
235 voxel->GetPointIds()->SetId(1, 8*iCell+1);
236 voxel->GetPointIds()->SetId(2, 8*iCell+2);
237 voxel->GetPointIds()->SetId(3, 8*iCell+3);
238 voxel->GetPointIds()->SetId(4, 8*iCell+4);
239 voxel->GetPointIds()->SetId(5, 8*iCell+5);
240 voxel->GetPointIds()->SetId(6, 8*iCell+6);
241 voxel->GetPointIds()->SetId(7, 8*iCell+7);
242
243 fpGrid->InsertNextCell(voxel->GetCellType(), voxel->GetPointIds());
244
245 double cols[] = {colour.GetRed(),
246 colour.GetGreen(),
247 colour.GetBlue(),
248 colour.GetAlpha()};
249 fpCellColourValues->InsertNextTuple(cols);
250 //fpCellColourIndices->InsertNextTuple1(iColour);
251
252
253 iCell++;
254
255 }
256 catch (const std::bad_cast&) {
258 ed << "Called for a mesh that is not a cubical mesh: " << box.GetName();
259 G4Exception("PseudoSceneForCubicalVertices", "visman0108", JustWarning, ed);
260 }
261 }
262}
263
265
266 // Modelling parameters
268 tmpMP.SetCulling(true); // This avoids drawing transparent...
269 tmpMP.SetCullingInvisible(true); // ... or invisble volumes.
270
271 // Physical volume model
272 const auto& container = mesh.GetContainerVolume();
273 const G4bool useFullExtent = true; // To avoid calculating the extent
274 G4PhysicalVolumeModel tmpPVModel(container,
276 G4Transform3D(), // so that positions are in local coordinates
277 &tmpMP,
278 useFullExtent);
279
280 // Graphics scene
281 if(mesh.GetMeshType() == G4Mesh::tetrahedron) {
282 PseudoSceneForTetCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
283 pointColourValues, cellColourValues,
284 pointColourIndices, cellColourIndices,
285 colourLUT, unstructuredGrid);
286 tmpPVModel.DescribeYourselfTo(pseudoScene);
287 //G4cout << pseudoScene.GetNumberOfPoints() << " " << pseudoScene.GetNumberOfCells() << " " << pseudoScene.GetNumberOfAddedCells() << G4endl;
288 //pseudoScene.DumpColourMap();
289 }
290 else if(mesh.GetMeshType() == G4Mesh::nested3DRectangular) {
291 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
292 pointColourValues, cellColourValues,
293 pointColourIndices, cellColourIndices,
294 colourLUT, unstructuredGrid);
295 tmpPVModel.DescribeYourselfTo(pseudoScene);
296 }
297 else if(mesh.GetMeshType() == G4Mesh::rectangle) {
298 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points,
299 pointColourValues, cellColourValues,
300 pointColourIndices, cellColourIndices,
301 colourLUT, unstructuredGrid);
302 tmpPVModel.DescribeYourselfTo(pseudoScene);
303 }
304
305}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
std::size_t MakeHash(const G4ThreeVector &v)
double z() const
double x() const
double y() const
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
const G4VisAttributes * GetVisAttributes() const
MeshType GetMeshType() const
Definition G4Mesh.hh:75
G4VPhysicalVolume * GetContainerVolume() const
Definition G4Mesh.hh:73
@ rectangle
Definition G4Mesh.hh:53
@ nested3DRectangular
Definition G4Mesh.hh:54
@ tetrahedron
Definition G4Mesh.hh:57
G4int GetMeshDepth() const
Definition G4Mesh.hh:76
void SetCulling(G4bool)
void SetCullingInvisible(G4bool)
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
Definition G4Tet.hh:56
void GetVertices(G4ThreeVector &anchor, G4ThreeVector &p1, G4ThreeVector &p2, G4ThreeVector &p3) const
Definition G4Tet.cc:286
G4String GetName() const
vtkSmartPointer< vtkRenderer > renderer
G4VtkVisContext vc
const G4Colour & GetColour() const
vtkSmartPointer< vtkDiscretizableColorTransferFunction > fpColourLUT
G4VtkUnstructuredGridPipeline(G4String name, const G4VtkVisContext &vc)