Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VtkSceneHandler.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//
27//
28//
29// John Allison 5th April 2001
30// A template for a simplest possible graphics driver.
31//?? Lines or sections marked like this require specialisation for your driver.
32
33#include "G4VtkSceneHandler.hh"
34
37#include "G4VPhysicalVolume.hh"
38#include "G4LogicalVolume.hh"
39#include "G4Polyline.hh"
40#include "G4Text.hh"
41#include "G4Circle.hh"
42#include "G4Square.hh"
43#include "G4Polyhedron.hh"
44#include "G4Mesh.hh"
45#include "G4PseudoScene.hh"
46#include "G4UnitsTable.hh"
47
48#include "G4SystemOfUnits.hh"
49#include "G4Material.hh"
50#include "G4Box.hh"
51
52#include <stdlib.h>
53
54namespace std
55{
56 inline void hash_combine(std::size_t) {}
57
58 template <typename T, typename... Rest>
59 inline void hash_combine(std::size_t &seed, const T &v, Rest... rest) {
60 std::hash<T> hasher;
61 seed ^= hasher(v) + 0x9e3779b9 + (seed << 6) + (seed >> 2);
62 std::hash_combine(seed, rest...);
63 }
64
65 template<> struct hash<G4VisAttributes> {
66 std::size_t operator()(const G4VisAttributes &va) const {
67 using std::size_t;
68 using std::hash;
69
70 std::size_t h = 0;
71
78 std::hash_combine(h,static_cast<int>(va.GetLineStyle()));
79
80 return h;
81 }
82 };
83
84 template<> struct hash<G4Polyhedron> {
85 std::size_t operator()(const G4Polyhedron &ph) const {
86 using std::size_t;
87 using std::hash;
88
89 G4bool notLastFace;
90 G4Point3D vertex[4];
91 G4int edgeFlag[4];
92 G4Normal3D normals[4];
93 G4int nEdges;
94
95 std::size_t h = 0;
96
97 do {
98 notLastFace = ph.GetNextFacet(nEdges, vertex, edgeFlag, normals);
99
100 for (int i = 0; i < nEdges; i++) {
101 std::size_t hx = std::hash<double>()(vertex[i].x());
102 std::size_t hy = std::hash<double>()(vertex[i].y());
103 std::size_t hz = std::hash<double>()(vertex[i].z());
104 std::hash_combine(h,hx);
105 std::hash_combine(h,hy);
106 std::hash_combine(h,hz);
107 }
108 } while (notLastFace);
109
110 return h;
111 }
112 };
113}
114
116// Counter for XXX scene handlers.
117
119 const G4String& name) :
120 G4VSceneHandler(system, fSceneIdCount++, name)
121{}
122
123#ifdef G4VTKDEBUG
124void G4VtkSceneHandler::PrintThings() {
125 G4cout << " with transformation " << fObjectTransformation.xx() << G4endl;
126 if (fpModel) {
127 G4cout << " from " << fpModel->GetCurrentDescription()
128 << " (tag " << fpModel->GetCurrentTag()
129 << ')';
130 }
131 else {
132 G4cout << "(not from a model)";
133 }
134 G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
135 if (pPVModel) {
136 G4cout << "\n current physical volume: " << pPVModel->GetCurrentPV()->GetName()
137 << "\n current logical volume : " << pPVModel->GetCurrentLV()->GetName() // There might be a problem with the LV pointer if this is a G4LogicalVolumeModel
138 << "\n current depth of geometry tree: " << pPVModel->GetCurrentDepth();
139 }
140 G4cout << G4endl;
141}
142#endif
143
145
147 // GetMarkerSize(text,sizeType);
148 if(fProcessing2D) {sizeType = screen;}
149 else {sizeType = world;}
150
151 // Get vis attributes - pick up defaults if none.
152 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(polyline.GetVisAttributes());
153 G4Color colour = pVA->GetColour();
154 G4double opacity = colour.GetAlpha();
155 G4double lineWidth = pVA->GetLineWidth();
156
157#ifdef G4VTKDEBUG
158 G4cout << "=================================" << G4endl;
159 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> vis hash: " << sizeType << " " << std::hash<G4VisAttributes>{}(*pVA) << G4endl;
160 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> sizeType: " << sizeType << G4endl;
161 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> isVisible: " << pVA->IsVisible() << G4endl;
162 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> isDaughtersInvisible: " << pVA->IsDaughtersInvisible() << G4endl;
163 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> colour: " << colour.GetRed() << " " << colour.GetGreen() << " " << colour.GetBlue() << G4endl;
164 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> alpha: " << colour.GetAlpha() << G4endl;
165 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyline& polyline) called> lineWidth: " << lineWidth << G4endl;
166#endif
167
168 if(sizeType == world) {
169 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
170 if (polylineVisAttributesMap.find(hash) == polylineVisAttributesMap.end()) {
171 polylineVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
172
173 vtkSmartPointer <vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
174 vtkSmartPointer <vtkCellArray> lines = vtkSmartPointer<vtkCellArray>::New();
175 vtkSmartPointer <vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
176 vtkSmartPointer <vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
177 vtkSmartPointer <vtkActor> actor = vtkSmartPointer<vtkActor>::New();
178
179 polyData->SetPoints(data);
180 polyData->SetLines(lines);
181 mapper->SetInputData(polyData);
182 actor->SetMapper(mapper);
183
184 // Setup actor and mapper
185 actor->GetProperty()->SetLineWidth(lineWidth);
186 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(),
187 colour.GetBlue());
188 actor->GetProperty()->SetOpacity(opacity);
189 actor->SetVisibility(1);
190 // actor->GetProperty()->BackfaceCullingOn();
191 // actor->GetProperty()->FrontfaceCullingOn();
192
193 auto pVtkViewer = dynamic_cast<G4VtkViewer*>(fpViewer);
194 pVtkViewer->renderer->AddActor(actor);
195
196 polylineDataMap.insert(
197 std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, data));
198 polylineLineMap.insert(
199 std::pair<std::size_t, vtkSmartPointer<vtkCellArray>>(hash, lines));
200 polylinePolyDataMap.insert(
201 std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash, polyData));
203 std::pair<std::size_t, vtkSmartPointer<vtkPolyDataMapper>>(hash,
204 mapper));
206 std::pair<std::size_t, vtkSmartPointer<vtkActor>>(hash, actor));
207 }
208
209 // Data data
210 const size_t nLines = polyline.size();
211
212 for (size_t i = 0; i < nLines; ++i) {
213 auto id = polylineDataMap[hash]->InsertNextPoint(polyline[i].x(),
214 polyline[i].y(),
215 polyline[i].z());
216 if (i < nLines - 1) {
217 vtkSmartPointer <vtkLine> line = vtkSmartPointer<vtkLine>::New();
218 line->GetPointIds()->SetId(0, id);
219 line->GetPointIds()->SetId(1, id + 1);
220 polylineLineMap[hash]->InsertNextCell(line);
221 }
222 }
223 }
224 else if (sizeType == screen ) {
225
226 }
227}
228
230
232 if(fProcessing2D) {sizeType = screen;}
233 else {sizeType = world;}
234
235 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(text.GetVisAttributes());
236 G4Color colour = pVA->GetColour();
237 G4double opacity = colour.GetAlpha();
238 // G4Text::Layout layout = text.GetLayout();
239 // G4double xOffset = text.GetXOffset();
240 // G4double yOffset = text.GetYOffset();
241
242 double x = text.GetPosition().x();
243 double y = text.GetPosition().y();
244 double z = text.GetPosition().z();
245
246#ifdef G4VTKDEBUG
247 G4cout << "=================================" << G4endl;
248 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> text: " << text.GetText() << " sizeType:" << sizeType << " " << fProcessing2D << G4endl;
249 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> colour: " << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
250 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> alpha: " << colour.GetAlpha() << G4endl;
251 G4cout << "G4VtkSeneHandler::AddPrimitive(const G4Text& text) called> position: " << x << " " << y << " " << z << G4endl;
252#endif
253
254 switch (sizeType) {
255 default:
256 case (screen): {
257 vtkSmartPointer <vtkTextActor> actor = vtkSmartPointer<vtkTextActor>::New();
258 actor->SetInput(text.GetText().c_str());
259 actor->GetPositionCoordinate()->SetCoordinateSystemToNormalizedViewport();
260 // actor->SetTextScaleModeToViewport();
261 actor->SetPosition((x+1.)/2.0, (y+1.)/2.);
262 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
263 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetBlue(), colour.GetGreen());
264 actor->GetTextProperty()->SetOpacity(opacity);
265
266 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
267 pVtkViewer->renderer->AddActor(actor);
268 break;
269 }
270 case world: {
271 vtkSmartPointer <vtkBillboardTextActor3D> actor = vtkSmartPointer<vtkBillboardTextActor3D>::New();
272 actor->SetInput(text.GetText().c_str());
273 actor->SetPosition(x, y, z);
274 actor->GetTextProperty()->SetFontSize(text.GetScreenSize());
275 actor->GetTextProperty()->SetColor(colour.GetRed(), colour.GetBlue(), colour.GetGreen());
276 actor->GetTextProperty()->SetOpacity(opacity);
277
278 auto *pVtkViewer = dynamic_cast<G4VtkViewer*>(fpViewer);
279 pVtkViewer->renderer->AddActor(actor);
280 break;
281 }
282 }
283}
284
286
287 MarkerSizeType sizeType;
288 G4double size= GetMarkerSize(circle, sizeType);
289 if(fProcessing2D) {sizeType = screen;}
290 else {sizeType = world;}
291
292 // Get vis attributes - pick up defaults if none.
294 G4Color colour = pVA->GetColour();
295 G4double opacity = colour.GetAlpha();
296 // G4bool isVisible = pVA->IsVisible();
297
298#ifdef G4VTKDEBUG
299 G4cout << "=================================" << G4endl;
300 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> " << " radius:" << size << " sizeType:" << sizeType << G4endl;
301 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> colour: " << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
302 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Circle& circle) called> alpha: " << colour.GetAlpha() << G4endl;
303#endif
304
305 if (sizeType == world) {
306 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
307 if (circleVisAttributesMap.find(hash) == circleVisAttributesMap.end()) {
308 circleVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
309
310 vtkSmartPointer<vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
311 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
312 vtkSmartPointer<vtkVertexGlyphFilter> filter = vtkSmartPointer<vtkVertexGlyphFilter>::New();
313 vtkSmartPointer<vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
314 vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
315
316 polyData->SetPoints(data);
317 filter->SetInputData(polyData);
318 mapper->SetInputConnection(filter->GetOutputPort());
319 actor->SetMapper(mapper);
320
321 // Setup actor and mapper
322 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
323 actor->GetProperty()->SetOpacity(opacity);
324 actor->SetVisibility(1);
325 actor->GetProperty()->SetRenderPointsAsSpheres(true);
326 actor->GetProperty()->SetPointSize(size*5);
327
328 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
329 pVtkViewer->renderer->AddActor(actor);
330
331 circleDataMap.insert(std::pair<std::size_t,vtkSmartPointer<vtkPoints>>(hash, data));
332 circlePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyData>>(hash, polyData));
333 circleFilterMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>>(hash, filter));
334 circlePolyDataMapperMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyDataMapper>>(hash, mapper));
335 circlePolyDataActorMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkActor>>(hash, actor));
336 }
337
338 // Data data point
340 G4Point3D posPrime = rot*circle.GetPosition();
341 circleDataMap[hash]->InsertNextPoint(fObjectTransformation.dx() + posPrime.x(),
342 fObjectTransformation.dy() + posPrime.y(),
343 fObjectTransformation.dz() + posPrime.z());
344 }
345 else if (sizeType == screen) {
346
347 }
348}
349
351
352 MarkerSizeType sizeType;
353 G4double size = GetMarkerSize (square, sizeType);
354
355 // Get vis attributes - pick up defaults if none.
356 const G4VisAttributes* pVA = fpViewer -> GetApplicableVisAttributes(square.GetVisAttributes());
357 G4Color colour = pVA->GetColour();
358 G4double opacity = colour.GetAlpha();
359 // G4bool isVisible = pVA->IsVisible();
360
361 // Draw in world coordinates.
362 vtkSmartPointer<vtkRegularPolygonSource> polygonSource = vtkSmartPointer<vtkRegularPolygonSource>::New();
363 polygonSource->SetNumberOfSides(4);
364 polygonSource->SetRadius(size);
365
366
367#ifdef G4VTKDEBUG
368 G4cout << "=================================" << G4endl;
369 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Square& square) called" << G4endl;
370 G4cout << square.GetPosition().x() << " " << square.GetPosition().y() << " " << square.GetPosition().z() << G4endl;
371 //PrintThings();
372#endif
373
374 if (sizeType == world) {
375 std::size_t hash = std::hash<G4VisAttributes>{}(*pVA);
376 if (squareVisAttributesMap.find(hash) == squareVisAttributesMap.end()) {
377 squareVisAttributesMap.insert(std::pair<std::size_t, const G4VisAttributes *>(hash, pVA));
378
379 vtkSmartPointer<vtkPoints> data = vtkSmartPointer<vtkPoints>::New();
380 vtkSmartPointer<vtkPolyData> polyData = vtkSmartPointer<vtkPolyData>::New();
381 vtkSmartPointer<vtkVertexGlyphFilter> filter = vtkSmartPointer<vtkVertexGlyphFilter>::New();
382 vtkSmartPointer<vtkPolyDataMapper> mapper = vtkSmartPointer<vtkPolyDataMapper>::New();
383 vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
384
385 polyData->SetPoints(data);
386 filter->SetInputData(polyData);
387 mapper->SetInputConnection(filter->GetOutputPort());
388 actor->SetMapper(mapper);
389
390 // Setup actor and mapper
391 actor->GetProperty()->SetColor(colour.GetRed(), colour.GetGreen(), colour.GetBlue());
392 actor->GetProperty()->SetOpacity(opacity);
393 actor->SetVisibility(1);
394 actor->GetProperty()->SetPointSize(size*5);
395
396 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
397 pVtkViewer->renderer->AddActor(actor);
398
399 squareDataMap.insert(std::pair<std::size_t,vtkSmartPointer<vtkPoints>>(hash, data));
400 squarePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyData>>(hash, polyData));
401 squareFilterMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkVertexGlyphFilter>>(hash, filter));
402 squarePolyDataMapperMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkPolyDataMapper>>(hash, mapper));
403 squarePolyDataActorMap.insert(std::pair<std::size_t, vtkSmartPointer < vtkActor>>(hash, actor));
404 }
405
406 // Data data point
408 G4Point3D posPrime = rot*square.GetPosition();
409 squareDataMap[hash]->InsertNextPoint(fObjectTransformation.dx() + posPrime.x(),
410 fObjectTransformation.dy() + posPrime.y(),
411 fObjectTransformation.dz() + posPrime.z());
412 }
413 else if (sizeType == screen) {
414
415 }
416}
418 AddPrimitiveTensorGlyph(polyhedron);
419}
420
422 // only have a single polydata object for each LV but bake in the PV
423 // transformation so clippers and cutters can be implemented
424
425 //Get colour, etc..
426 if (polyhedron.GetNoFacets() == 0) {return;}
427
428 // Get vis attributes - pick up defaults if none.
430 G4Color colour = pVA->GetColour();
431 // G4bool isVisible = pVA->IsVisible();
432 // G4double lineWidth = pVA->GetLineWidth();
433 // G4VisAttributes::LineStyle lineStyle = pVA->GetLineStyle();
434
435 // G4double lineWidthScale = drawing_style.GetGlobalLineWidthScale();
436
437 // Get view parameters that the user can force through the vis attributes, thereby over-riding the current view parameter.
439 //G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
440
441#ifdef G4VTKDEBUG
442 G4cout << "=================================" << G4endl;
443 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> " << G4endl;
444 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> colour:" << colour.GetRed() << " " << colour.GetBlue() << " " << colour.GetGreen() << G4endl;
445 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> alpha:" << colour.GetAlpha() << G4endl;
446 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> lineWidth:" << lineWidth << G4endl;
447 G4cout << "G4VtkSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called> lineStyle:" << lineStyle << G4endl;
448#endif
449
450 //std::size_t vhash = std::hash<G4VisAttributes>{}(*pVA);
451
452 std::size_t vhash = 0;
453 std::hash_combine(vhash,static_cast<int>(drawing_style));
454
455 std::size_t phash = std::hash<G4Polyhedron>{}(polyhedron);
456
457 std::size_t hash = 0;
458 std::hash_combine(hash, phash);
459 std::hash_combine(hash, vhash);
460
461 if (polyhedronPolyDataMap.find(hash) == polyhedronPolyDataMap.end()) {
462
463 vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
464 vtkSmartPointer<vtkCellArray> polys = vtkSmartPointer<vtkCellArray>::New();
465 vtkSmartPointer<vtkPolyData> polydata = vtkSmartPointer<vtkPolyData>::New();
466
467 G4bool notLastFace;
468 int iVert = 0;
469 do {
470 G4Point3D vertex[4];
471 G4int edgeFlag[4];
472 G4Normal3D normals[4];
473 G4int nEdges;
474 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normals);
475
476 vtkSmartPointer<vtkIdList> poly = vtkSmartPointer<vtkIdList>::New();
477 // loop over vertices
478 for (int i = 0; i < nEdges; i++) {
479 points->InsertNextPoint(vertex[i].x(), vertex[i].y(), vertex[i].z());
480 poly->InsertNextId(iVert);
481 iVert++;
482 }
483 polys->InsertNextCell(poly);
484
485 } while (notLastFace);
486
487 polydata->SetPoints(points);
488 polydata->SetPolys(polys);
489
490 polyhedronDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, points));
491 polyhedronPolyMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkCellArray>>(hash, polys));
492 polyhedronPolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash, polydata));
493 polyhedronPolyDataCountMap.insert(std::pair<std::size_t, std::size_t>(hash, 0));
494
495 vtkSmartPointer<vtkPoints> instancePosition = vtkSmartPointer<vtkPoints>::New();
496 vtkSmartPointer<vtkDoubleArray> instanceRotation = vtkSmartPointer<vtkDoubleArray>::New();
497 vtkSmartPointer<vtkDoubleArray> instanceColors = vtkSmartPointer<vtkDoubleArray>::New();
498 vtkSmartPointer<vtkPolyDataMapper> instanceMapper = vtkSmartPointer<vtkPolyDataMapper>::New();
499 vtkSmartPointer<vtkActor> instanceActor = vtkSmartPointer<vtkActor>::New();
500 instanceColors->SetName("colors");
501
502 instanceColors->SetNumberOfComponents(4);
503 instanceRotation->SetNumberOfComponents(9);
504
505 vtkSmartPointer<vtkPolyData> instancePolyData = vtkSmartPointer<vtkPolyData>::New();
506 instancePolyData->SetPoints(instancePosition);
507 instancePolyData->GetPointData()->SetTensors(instanceRotation);
508 instancePolyData->GetPointData()->SetVectors(instanceColors);
509 instancePolyData->GetPointData()->SetScalars(instanceColors);
510
511 vtkSmartPointer<vtkCleanPolyData> filterClean = vtkSmartPointer<vtkCleanPolyData>::New();
512 filterClean->PointMergingOn();
513 filterClean->AddInputData(polydata);
514
515 vtkSmartPointer<vtkTriangleFilter> filterTriangle = vtkSmartPointer<vtkTriangleFilter>::New();
516 filterTriangle->SetInputConnection(filterClean->GetOutputPort());
517
518 vtkSmartPointer<vtkPolyDataNormals> filterNormals = vtkSmartPointer<vtkPolyDataNormals>::New();
519 filterNormals->SetFeatureAngle(45);
520 filterNormals->SetInputConnection(filterTriangle->GetOutputPort());
521
522 vtkSmartPointer<vtkFeatureEdges> filterEdge = vtkSmartPointer<vtkFeatureEdges>::New();
523 filterEdge->SetFeatureEdges(1);
524 filterEdge->SetManifoldEdges(0);
525 filterEdge->SetBoundaryEdges(0);
526 filterEdge->SetFeatureAngle(45); // TODO need to have a function and command to set this
527 filterEdge->SetInputConnection(filterTriangle->GetOutputPort());
528
529 vtkSmartPointer<vtkTensorGlyphColor> tensorGlyph = vtkSmartPointer<vtkTensorGlyphColor>::New();
530 tensorGlyph->SetInputData(instancePolyData);
531 tensorGlyph->SetSourceConnection(filterNormals->GetOutputPort());
532 tensorGlyph->ColorGlyphsOn();
533 tensorGlyph->ScalingOff();
534 tensorGlyph->ThreeGlyphsOff();
535 tensorGlyph->ExtractEigenvaluesOff();
536 tensorGlyph->SetColorModeToScalars();
537 tensorGlyph->Update();
538
539 instanceMapper->SetInputData(tensorGlyph->GetOutput());
540 instanceMapper->SetColorModeToDirectScalars();
541
542 instanceActor->SetMapper(instanceMapper);
543 // instanceActor->GetProperty()->SetLineWidth(10);
544 instanceActor->SetVisibility(1);
545
546 if(drawing_style == G4ViewParameters::hsr) {
547
548 }
549
550 if(drawing_style == G4ViewParameters::hlr) {
551 }
552
553 if(drawing_style == G4ViewParameters::wireframe) {
554 instanceActor->GetProperty()->SetRepresentationToWireframe();
555 }
556
557 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
558 pVtkViewer->renderer->AddActor(instanceActor);
559
560 instancePositionMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPoints>>(hash, instancePosition));
561 instanceRotationMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkDoubleArray>>(hash, instanceRotation));
562 instanceColoursMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkDoubleArray>>(hash, instanceColors));
563 instancePolyDataMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkPolyData>>(hash,instancePolyData));
564 instanceActorMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkActor>>(hash, instanceActor));
565 instanceTensorGlyphMap.insert(std::pair<std::size_t, vtkSmartPointer<vtkTensorGlyphColor>>(hash, tensorGlyph));
566 }
567
569
570 double red = colour.GetRed();
571 double green = colour.GetGreen();
572 double blue = colour.GetBlue();
573 double alpha = colour.GetAlpha();
574
575 instanceColoursMap[hash]->InsertNextTuple4(red, green, blue, alpha);
576
577 instancePositionMap[hash]->InsertNextPoint(fObjectTransformation.dx(),
580
582
583 instanceRotationMap[hash]->InsertNextTuple9(fInvObjTrans.xx(), fInvObjTrans.xy(),fInvObjTrans.xz(),
584 fInvObjTrans.yx(), fInvObjTrans.yy(),fInvObjTrans.yz(),
585 fInvObjTrans.zx(), fInvObjTrans.zy(),fInvObjTrans.zz());
586
587}
588
590{
591 // only have a single polydata object for each LV but bake in the PV
592 // transformation so clippers and cutters can be implemented
593 AddPrimitiveTensorGlyph(polyhedron);
594}
596 for (auto it = polylineDataMap.begin(); it != polylineDataMap.end(); it++)
597 {
598 it->second->Modified();
599 }
600
601 for (auto it = polylineLineMap.begin(); it != polylineLineMap.end(); it++)
602 {
603 it->second->Modified();
604 }
605 for (auto it = circleDataMap.begin(); it != circleDataMap.end(); it++)
606 {
607 it->second->Modified();
608 }
609
610 for (auto it = squareDataMap.begin(); it != squareDataMap.end(); it++)
611 {
612 it->second->Modified();
613 }
614
615 for(auto it = instancePositionMap.begin(); it != instancePositionMap.end(); it++)
616 {
617 it->second->Modified();
618 }
619
620 for(auto it = instanceRotationMap.begin(); it != instanceRotationMap.end(); it++)
621 {
622 it->second->Modified();
623 }
624
625 for(auto it = instanceColoursMap.begin(); it != instanceColoursMap.end(); it++)
626 {
627 it->second->Modified();
628 }
629
630 for(auto it = instancePolyDataMap.begin(); it != instancePolyDataMap.end(); it++)
631 {
632 it->second->Modified();
633 }
634
635 for(auto it = instanceActorMap.begin(); it != instanceActorMap.end(); it++)
636 {
637 it->second->Modified();
638 }
639
640 for(auto it = instanceTensorGlyphMap.begin(); it != instanceTensorGlyphMap.end(); it++)
641 {
642 it->second->Update();
643 }
644
645#ifdef G4VTKDEBUG
646 G4cout << "G4VtkSceneHandler::Modified() polyline styles: " << polylineVisAttributesMap.size() << G4endl;
647 for (auto it = polylineLineMap.begin(); it != polylineLineMap.end(); it++)
648 {
649 G4cout << "G4VtkSceneHandler::Modified() polyline segments: "
650 << it->second->GetNumberOfCells() << G4endl;
651 }
652 G4cout << "G4VtkSceneHandler::Modified() circle styles: " << circleVisAttributesMap.size() << G4endl;
653 for (auto it = circleDataMap.begin(); it != circleDataMap.end(); it++)
654 {
655 G4cout << "G4VtkSceneHandler::Modified() circles: "
656 << it->second->GetNumberOfPoints() << G4endl;
657 }
658
659 G4cout << "G4VtkSceneHandler::Modified() square styles: " << squareVisAttributesMap.size() << G4endl;
660 for (auto it = squareDataMap.begin(); it != squareDataMap.end(); it++)
661 {
662 G4cout << "G4VtkSceneHanler::Modified() squares: "
663 << it->second->GetNumberOfPoints() << G4endl;
664 }
665
666 G4cout << "G4VtkSceneHandler::Modified() unique polyhedra: " << polyhedronDataMap.size() << G4endl;
667
668 int nPlacements = 0;
669 int nCells = 0;
670 for (auto it = polyhedronPolyDataMap.begin(); it != polyhedronPolyDataMap.end(); it++) {
671 G4cout << "G4VtkSceneHandler::Modified() polyhedronPolyData: " << it->second->GetPoints()->GetNumberOfPoints() << " " << it->second->GetPolys()->GetNumberOfCells() << " " << polyhedronPolyDataCountMap[it->first] <<G4endl;
672 nCells += it->second->GetPolys()->GetNumberOfCells()*polyhedronPolyDataCountMap[it->first];
673 nPlacements += polyhedronPolyDataCountMap[it->first];
674 }
675 G4cout << "G4VtkSceneHandler::Modified() polyhedronPolyData: " << nPlacements << " " << nCells << G4endl;
676#endif
677}
678
680
682 polylineDataMap.clear();
683 polylineLineMap.clear();
684 polylinePolyDataMap.clear();
687
688 for (auto& v : polylineDataMap)
689 {v.second->Reset();}
690 for (auto& v : polylineLineMap)
691 {v.second->Reset();}
692
694 circleDataMap.clear();
695 circlePolyDataMap.clear();
696 circleFilterMap.clear();
699
701 squareDataMap.clear();
702 squarePolyDataMap.clear();
703 squareFilterMap.clear();
706
707 for (auto& v : squareDataMap)
708 {v.second->Reset();}
709
711 polyhedronDataMap.clear();
712 polyhedronPolyMap.clear();
713 polyhedronPolyDataMap.clear();
715
716 instancePositionMap.clear();
717 instanceRotationMap.clear();
718 instanceColoursMap.clear();
719 instancePolyDataMap.clear();
721 instanceActorMap.clear();
722}
723
725
727
728#if 0
729 return;
730
731 const G4VModel* pv_model = GetModel();
732 if (!pv_model) { return ; }
733
734 G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
735 if (!pPVModel) { return ; }
736
737 //-- debug information
738 if(1) {
739 G4VPhysicalVolume *pv = pPVModel->GetCurrentPV();
741 G4cout << "name=" << box.GetName() << " volumeType=" << pv->VolumeType() << " pvName=" << pv->GetName() << " lvName=" << lv->GetName() << " multiplicity=" << pv->GetMultiplicity() << " isparametrised=" << pv->IsParameterised() << " isreplicated=" << pv->IsReplicated() << " parametrisation=" << pv->GetParameterisation() << G4endl;
742 }
743
744 if(0) {
745 G4Material *mat = pPVModel->GetCurrentMaterial();
746 G4String name = mat->GetName();
747 G4double dens = mat->GetDensity()/(g/cm3);
748 G4int copyNo = pPVModel->GetCurrentPV()->GetCopyNo();
749 G4int depth = pPVModel->GetCurrentDepth();
750 G4cout << " name : " << box.GetName() << G4endl;
751 G4cout << " copy no.: " << copyNo << G4endl;
752 G4cout << " depth : " << depth << G4endl;
753 G4cout << " density : " << dens << " [g/cm3]" << G4endl;
754 G4cout << " location: " << pPVModel->GetCurrentPV()->GetObjectTranslation() << G4endl;
755 G4cout << " Multiplicity : " << pPVModel->GetCurrentPV()->GetMultiplicity() << G4endl;
756 G4cout << " Is replicated? : " << pPVModel->GetCurrentPV()->IsReplicated() << G4endl;
757 G4cout << " Is parameterised? : " << pPVModel->GetCurrentPV()->IsParameterised() << G4endl;
758 G4cout << " top phys. vol. name : " << pPVModel->GetTopPhysicalVolume()->GetName() << G4endl;
759 }
760#endif
761}
762
764{
765#ifdef G4VTKDEBUG
766 G4cout << "G4VtkSceneHandler::AddCompound" << G4endl;
767#endif
768
769 if(mesh.GetMeshType() != G4Mesh::rectangle || mesh.GetMeshDepth() != 3)
770 {
772 }
773
774 auto container = mesh.GetContainerVolume();
775 auto rep1 = container->GetLogicalVolume()->GetDaughter(0);
776 auto rep2 = rep1->GetLogicalVolume()->GetDaughter(0);
777 auto rep3 = rep2->GetLogicalVolume()->GetDaughter(0);
778
779 EAxis rep1_axis, rep2_axis, rep3_axis;
780 G4int rep1_nReplicas, rep2_nReplicas, rep3_nReplicas;
781 G4double rep1_width, rep2_width, rep3_width;
782 G4double rep1_offset, rep2_offset, rep3_offset;
783 G4bool rep1_consuming, rep2_consuming, rep3_consuming;
784
785 rep1->GetReplicationData(rep1_axis,rep1_nReplicas, rep1_width, rep1_offset, rep1_consuming);
786 rep2->GetReplicationData(rep2_axis,rep2_nReplicas, rep2_width, rep2_offset, rep2_consuming);
787 rep3->GetReplicationData(rep3_axis,rep3_nReplicas, rep3_width, rep3_offset, rep3_consuming);
788
789 // Instantiate a temporary G4PhysicalVolumeModel
791 tmpMP.SetCulling(false); // This avoids drawing transparent...
792 tmpMP.SetCullingInvisible(false); // ... or invisble volumes.
793 const G4bool useFullExtent = false; // To avoid calculating the extent
794
796 G4Transform3D(), &tmpMP, useFullExtent);
797
798 // Instantiate a pseudo scene so that we can make a "private" descent and fill vtkImageData
799 vtkSmartPointer<vtkImageData> imagedata = vtkSmartPointer<vtkImageData>::New();
800 imagedata->SetDimensions(rep1_nReplicas+1,rep2_nReplicas+1,rep3_nReplicas+1);
801 imagedata->SetSpacing(rep1_width,rep2_width,rep2_width);
802 imagedata->SetOrigin(0,0,0);
803 imagedata->AllocateScalars(VTK_DOUBLE, 1);
804
805 G4double halfX = 0., halfY = 0., halfZ = 0.;
806 struct PseudoScene : public G4PseudoScene
807 {
808 PseudoScene(
809 G4PhysicalVolumeModel* pvModel, // input...the following are output
810 vtkImageData *vtkId,
811 G4int nx, G4int ny, G4int nz,
812 G4double& halfX, G4double& halfY, G4double& halfZ) : fpPVModel(pvModel)
813 , fVtkId(vtkId)
814 , fNx(nx)
815 , fNy(ny)
816 , fNz(nz)
817 , fHalfX(halfX)
818 , fHalfY(halfY)
819 , fHalfZ(halfZ)
820 {}
821 using G4PseudoScene::AddSolid; // except for...
822 void AddSolid(const G4Box& box)
823 {
824 const G4Colour& colour = fpPVModel->GetCurrentLV()->GetVisAttributes()->GetColour();
825 // G4double density = fpPVModel->GetCurrentLV()->GetMaterial()->GetDensity();
826 const G4ThreeVector& position = fpCurrentObjectTransformation->getTranslation();
827 fHalfX = box.GetXHalfLength();
828 fHalfY = box.GetYHalfLength();
829 fHalfZ = box.GetZHalfLength();
830 fVtkId->SetScalarComponentFromDouble(int(position.x()/(2*fHalfX))+fNx,
831 int(position.y()/(2*fHalfY))+fNy,
832 int(position.z()/(2*fHalfZ))+fNz,
833 0,
834 (colour.GetRed() + colour.GetBlue() + colour.GetGreen())/3.0);
835 }
836 G4PhysicalVolumeModel* fpPVModel;
837 vtkImageData *fVtkId;
838 G4int fNx, fNy, fNz;
839 G4double &fHalfX, &fHalfY, &fHalfZ;
840 }
841 // construct the pseudoScene
842 pseudoScene(&tmpPVModel, imagedata, rep1_nReplicas, rep2_nReplicas, rep3_nReplicas, halfX, halfY, halfZ);
843
844 // Make private descent into the nested parameterisation
845 tmpPVModel.DescribeYourselfTo(pseudoScene);
846
847 vtkSmartPointer<vtkOpenGLGPUVolumeRayCastMapper> volumeMapper = vtkSmartPointer<vtkOpenGLGPUVolumeRayCastMapper>::New();
848 volumeMapper->SetInputData(imagedata);
849
850 vtkNew<vtkVolume> volume;
851 volume->SetMapper(volumeMapper);
852
853 vtkSmartPointer<vtkMatrix4x4> vtkMatrix = vtkSmartPointer<vtkMatrix4x4>::New();
854
855 vtkMatrix->SetElement(0,0,fObjectTransformation.xx());
856 vtkMatrix->SetElement(0,1,fObjectTransformation.xy());
857 vtkMatrix->SetElement(0,2,fObjectTransformation.xz());
858
859 vtkMatrix->SetElement(1,0,fObjectTransformation.yx());
860 vtkMatrix->SetElement(1,1,fObjectTransformation.yy());
861 vtkMatrix->SetElement(1,2,fObjectTransformation.yz());
862
863 vtkMatrix->SetElement(2,0,fObjectTransformation.zx());
864 vtkMatrix->SetElement(2,1,fObjectTransformation.zy());
865 vtkMatrix->SetElement(2,2,fObjectTransformation.zz());
866
867 vtkMatrix->SetElement(3,0, fObjectTransformation.dx());
868 vtkMatrix->SetElement(3,1, fObjectTransformation.dy());
869 vtkMatrix->SetElement(3,2, fObjectTransformation.dz());
870 vtkMatrix->SetElement(3,3, 1);
871
872 volume->SetUserMatrix(vtkMatrix);
873
874 auto *pVtkViewer = dynamic_cast<G4VtkViewer *>(fpViewer);
875 pVtkViewer->renderer->AddVolume(volume);
876}
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition: G4Colour.cc:155
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
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
const G4String & GetName() const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4String & GetName() const
Definition: G4Material.hh:172
Definition: G4Mesh.hh:48
MeshType GetMeshType() const
Definition: G4Mesh.hh:75
G4VPhysicalVolume * GetContainerVolume() const
Definition: G4Mesh.hh:73
@ rectangle
Definition: G4Mesh.hh:53
G4int GetMeshDepth() const
Definition: G4Mesh.hh:76
void SetCulling(G4bool)
void SetCullingInvisible(G4bool)
G4VPhysicalVolume * GetCurrentPV() const
void DescribeYourselfTo(G4VGraphicsScene &)
G4LogicalVolume * GetCurrentLV() const
G4VPhysicalVolume * GetTopPhysicalVolume() const
G4Material * GetCurrentMaterial() const
void AddSolid(const G4Box &solid)
Definition: G4Text.hh:72
G4String GetText() const
G4double GetScreenSize() const
G4Point3D GetPosition() const
virtual G4String GetCurrentDescription() const
Definition: G4VModel.cc:51
virtual G4String GetCurrentTag() const
Definition: G4VModel.cc:46
virtual G4bool IsReplicated() const =0
virtual EVolume VolumeType() const =0
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
G4VModel * GetModel() const
G4Transform3D fObjectTransformation
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
G4VViewer * fpViewer
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
virtual void AddSolid(const G4Box &)
virtual void AddCompound(const G4VTrajectory &)
G4String GetName() const
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
G4double GetLineWidth() const
G4bool IsDaughtersInvisible() const
LineStyle GetLineStyle() const
const G4Colour & GetColour() const
G4bool IsVisible() const
const G4VisAttributes * GetVisAttributes() const
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > polylinePolyDataMapperMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polylineLineMap
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > circleFilterMap
std::map< std::size_t, vtkSmartPointer< vtkVertexGlyphFilter > > squareFilterMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > instancePolyDataMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > circlePolyDataMap
std::map< std::size_t, const G4VisAttributes * > polylineVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceColoursMap
std::map< std::size_t, vtkSmartPointer< vtkDoubleArray > > instanceRotationMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > instancePositionMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > circleDataMap
void AddPrimitiveBakedTransform(const G4Polyhedron &)
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polylinePolyDataMap
static G4int fSceneIdCount
std::map< std::size_t, vtkSmartPointer< vtkPoints > > squareDataMap
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > polyhedronPolyDataMap
void AddPrimitive(const G4Polyline &)
std::map< std::size_t, std::size_t > polyhedronPolyDataCountMap
std::map< std::size_t, const G4VisAttributes * > squareVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > circlePolyDataMapperMap
G4VtkSceneHandler(G4VGraphicsSystem &system, const G4String &name)
std::map< std::size_t, vtkSmartPointer< vtkActor > > instanceActorMap
std::map< std::size_t, vtkSmartPointer< vtkCellArray > > polyhedronPolyMap
std::map< std::size_t, const G4VisAttributes * > polyhedronVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > squarePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPolyDataMapper > > squarePolyDataMapperMap
void AddSolid(const G4Box &box)
std::map< std::size_t, vtkSmartPointer< vtkPolyData > > squarePolyDataMap
void AddPrimitiveTensorGlyph(const G4Polyhedron &)
std::map< std::size_t, vtkSmartPointer< vtkTensorGlyphColor > > instanceTensorGlyphMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > polylinePolyDataActorMap
std::map< std::size_t, const G4VisAttributes * > circleVisAttributesMap
std::map< std::size_t, vtkSmartPointer< vtkActor > > circlePolyDataActorMap
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polylineDataMap
void AddCompound(const G4Mesh &mesh)
std::map< std::size_t, vtkSmartPointer< vtkPoints > > polyhedronDataMap
vtkNew< vtkRenderer > renderer
Definition: G4VtkViewer.hh:186
double dy() const
Definition: Transform3D.h:287
double zz() const
Definition: Transform3D.h:281
double yz() const
Definition: Transform3D.h:272
double dz() const
Definition: Transform3D.h:290
CLHEP::HepRotation getRotation() const
double dx() const
Definition: Transform3D.h:284
double xy() const
Definition: Transform3D.h:260
double zx() const
Definition: Transform3D.h:275
double yx() const
Definition: Transform3D.h:266
Transform3D inverse() const
Definition: Transform3D.cc:141
double zy() const
Definition: Transform3D.h:278
double xx() const
Definition: Transform3D.h:257
double yy() const
Definition: Transform3D.h:269
double xz() const
Definition: Transform3D.h:263
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const
EAxis
Definition: geomdefs.hh:54
void hash_combine(std::size_t)
std::size_t operator()(const G4Polyhedron &ph) const
std::size_t operator()(const G4VisAttributes &va) const