Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkViewer.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 <cmath>
27
28#include "G4VtkViewer.hh"
29
31
32#include "G4Transform3D.hh"
33#include "G4VSceneHandler.hh"
39#include "G4VtkSceneHandler.hh"
40#include "G4VtkUtility.hh"
41#include "G4VtkVisContext.hh"
42
43#include <vtk3DSImporter.h>
44#include <vtkBMPWriter.h>
45#include <vtkIVExporter.h> // open inventor
46#include <vtkImageWriter.h>
47#include <vtkImplicitPlaneRepresentation.h>
48#include <vtkImplicitPlaneWidget2.h>
49#include <vtkJPEGWriter.h>
50#include <vtkLightCollection.h>
51#include <vtkOBJExporter.h>
52#include <vtkOBJImporter.h>
53#include <vtkGLTFExporter.h>
54#include <vtkOOGLExporter.h>
55#include <vtkX3DExporter.h>
56#include <vtkJSONRenderWindowExporter.h>
57#include <vtkVtkJSSceneGraphSerializer.h>
58//#include <vtkBufferedArchiver.h>
59#include <vtkPNGWriter.h>
60#include <vtkPNMWriter.h>
61#include <vtkPOVExporter.h>
62#include <vtkPostScriptWriter.h>
63#include <vtkRIBExporter.h> // Renderman
64#include <vtkRendererCollection.h>
65#include <vtkSingleVTPExporter.h>
66#include <vtkTIFFWriter.h>
67#include <vtkVRMLExporter.h>
68#include <vtkVRMLImporter.h>
69#include <vtkWindowToImageFilter.h>
70#include <vtkX3DExporter.h>
71
72// Readers (vtkDataReader)
73
74#include <vtkCameraPass.h>
75#include <vtkOpenGLRenderer.h>
76#include <vtkRenderPass.h>
77#include <vtkRenderPassCollection.h>
78#include <vtkSequencePass.h>
79#include <vtkShadowMapBakerPass.h>
80#include <vtkShadowMapPass.h>
81
83 : G4VViewer(sceneHandler, sceneHandler.IncrementViewCount(), name)
84{
85 vtkObject::GlobalWarningDisplayOff();
86
87 // Set default and current view parameters
88 fVP.SetAutoRefresh(true);
89 fDefaultVP.SetAutoRefresh(true);
90}
91
93{
94 _renderWindow = vtkRenderWindow::New();
95 renderWindowInteractor = vtkRenderWindowInteractor::New();
96
97#ifdef G4VTKDEBUG
98 G4cout << "G4VtkViewer::G4VtkViewer" << G4endl;
99 G4cout << "G4VtkViewer::G4VtkViewer> " << fVP.GetWindowSizeHintX() << " "
100 << fVP.GetWindowSizeHintY() << G4endl;
101 G4cout << "G4VtkViewer::G4VtkViewer> " << fVP.GetWindowLocationHintX() << " "
102 << fVP.GetWindowLocationHintY() << G4endl;
103#endif
104
105 // Need windowSizeX/Y - obtain from _renderWindow?
106 G4int screenSizeX = _renderWindow->GetScreenSize()[0];
107 G4int screenSizeY = _renderWindow->GetScreenSize()[1];
108 G4int positionX = fVP.GetWindowLocationHintX();
109 if (fVP.IsWindowLocationHintXNegative()) {
110 positionX = screenSizeX + positionX - fVP.GetWindowSizeHintX();
111 }
112 G4int positionY = fVP.GetWindowLocationHintY();
113 if (!fVP.IsWindowLocationHintYNegative()) {
114 positionY = screenSizeY + positionY - fVP.GetWindowSizeHintY();
115 }
116 _renderWindow->SetPosition(positionX, positionY);
117#ifdef __APPLE__
118 // Adjust window size for Apple to make it correspond to OpenGL.
119 // Maybe it's OpenGL that shoud be adjusted.
120 const G4double pixelFactor = 2.;
121#else
122 const G4double pixelFactor = 1.;
123#endif
124 _renderWindow->SetSize(pixelFactor * fVP.GetWindowSizeHintX(),
125 pixelFactor * fVP.GetWindowSizeHintY());
126 _renderWindow->SetWindowName("Vtk viewer");
127
128 _renderWindow->AddRenderer(renderer);
129 renderWindowInteractor->SetRenderWindow(_renderWindow);
130
131 // TODO proper camera parameter settings
132 camera->SetPosition(0, 0, 1000);
133 camera->SetFocalPoint(0, 0, 0);
134 renderer->SetActiveCamera(camera);
135
136 // Hidden line removal
137 renderer->SetUseHiddenLineRemoval(0);
138
139 // Shadows
140 renderer->SetUseShadows(0);
141
142 // Set callback to match VTK parameters to Geant4
143 geant4Callback->SetGeant4ViewParameters(&fVP);
144 renderer->AddObserver(vtkCommand::EndEvent, geant4Callback);
145
147 renderWindowInteractor->SetInteractorStyle(style);
148}
149
151{
152#ifdef G4VTKDEBUG
153 G4cout << "G4VtkViewer::~G4VtkViewer()" << G4endl;
154#endif
155}
156
158{
159#ifdef G4VTKDEBUG
160 G4cout << "G4VtkViewer::SetView()" << G4endl;
161#endif
162 // background colour
163 const G4Colour backgroundColour = fVP.GetBackgroundColour();
164 renderer->SetBackground(backgroundColour.GetRed(), backgroundColour.GetGreen(),
165 backgroundColour.GetBlue());
166
167 // target and camera positions
168 G4double radius = fSceneHandler.GetExtent().GetExtentRadius();
169 if (radius <= 0.) {
170 radius = 1.;
171 }
172
173 if (firstSetView) cameraDistance = fVP.GetCameraDistance(radius);
174 G4Point3D viewpointDirection = fVP.GetViewpointDirection();
175 G4Point3D targetPosition = fVP.GetCurrentTargetPoint();
176 G4double zoomFactor = fVP.GetZoomFactor();
177 G4double fieldHalfAngle = fVP.GetFieldHalfAngle();
178
179 G4Point3D cameraPosition = targetPosition + viewpointDirection.unit() /zoomFactor * cameraDistance;
180
181 vtkCamera* activeCamera = renderer->GetActiveCamera();
182 activeCamera->SetFocalPoint(targetPosition.x(), targetPosition.y(), targetPosition.z());
183 activeCamera->SetViewAngle(2*fieldHalfAngle / CLHEP::pi * 180);
184 activeCamera->SetPosition(cameraPosition.x(), cameraPosition.y(), cameraPosition.z());
185 activeCamera->SetParallelScale(cameraDistance / zoomFactor);
186
187 if (fieldHalfAngle == 0) {
188 activeCamera->SetParallelProjection(1);
189 }
190 else {
191 activeCamera->SetParallelProjection(0);
192 }
193
194 // need to set camera distance and parallel scale on first set view
195 if (firstSetView) {
196 geant4Callback->SetVtkInitialValues(cameraDistance, cameraDistance);
197 activeCamera->SetParallelScale(cameraDistance);
198 firstSetView = false;
199 }
200
201 // camera up direction
202 const G4Vector3D upVector = fVP.GetUpVector();
203 renderer->GetActiveCamera()->SetViewUp(upVector.x(), upVector.y(), upVector.z());
204
205 // Light
206 const G4Vector3D lightDirection = fVP.GetLightpointDirection();
207 G4bool lightsMoveWithCamera = fVP.GetLightsMoveWithCamera();
208 G4Vector3D lightPosition = targetPosition + lightDirection.unit() * cameraDistance;
209
210 vtkLightCollection* currentLights = renderer->GetLights();
211 if (currentLights->GetNumberOfItems() != 0) {
212 auto currentLight = dynamic_cast<vtkLight*>(currentLights->GetItemAsObject(0));
213 if (currentLight != nullptr) {
214 currentLight->SetPosition(lightPosition.x(), lightPosition.y(), lightPosition.z());
215 if (lightsMoveWithCamera) {
216 currentLight->SetLightTypeToCameraLight();
217 }
218 else {
219 currentLight->SetLightTypeToSceneLight();
220 }
221 }
222 }
223
224 // cut away
225 if (fVP.IsCutaway()) {
226 G4cout << "Add cutaway planes" << G4endl;
227 }
228
229 // section
230 if (fVP.IsSection()) {
231 G4cout << "Add section" << G4endl;
232 }
233}
234
236{
237#ifdef G4VTKDEBUG
238 G4cout << "G4VtkViewer::ClearView()" << G4endl;
239#endif
240
241 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
242 G4VtkStore& ts = fVtkSceneHandler.GetTransientStore();
243 ts.Clear();
244
245 G4VtkStore& s = fVtkSceneHandler.GetStore();
246 s.Clear();
247}
248
250{
251#ifdef G4VTKDEBUG
252 G4cout << "G4VtkViewer::DrawView()" << G4endl;
253#endif
254
255 // First, a view should decide when to re-visit the G4 kernel.
256 // Sometimes it might not be necessary, e.g., if the scene is stored
257 // in a graphical database (e.g., OpenGL's display lists) and only
258 // the viewing angle has changed. But graphics systems without a
259 // graphical database will always need to visit the G4 kernel.
260
261 NeedKernelVisit(); // Default is - always visit G4 kernel.
262
263 // Note: this routine sets the fNeedKernelVisit flag of *all* the
264 // views of the scene.
265
266 ProcessView(); // The basic logic is here.
267
268 // Add HUD
269 AddViewHUD();
270
271 // Add clipper and cutter widgets
272 auto g4p = G4Plane3D();
275
276 // Add camera orientation widget
278
279 // ...before finally...
280 FinishView(); // Flush streams and/or swap buffers.
281}
282
284{
285 _renderWindow->SetMultiSamples(0);
286
287 vtkNew<vtkShadowMapPass> shadows;
288 vtkNew<vtkSequencePass> seq;
289
290 vtkNew<vtkRenderPassCollection> passes;
291 passes->AddItem(shadows->GetShadowMapBakerPass());
292 passes->AddItem(shadows);
293 seq->SetPasses(passes);
294
295 vtkNew<vtkCameraPass> cameraP;
296 cameraP->SetDelegatePass(seq);
297
298 // tell the renderer to use our render pass pipeline
299 auto glrenderer = dynamic_cast<vtkOpenGLRenderer*>(renderer.GetPointer());
300 glrenderer->SetPass(cameraP);
301}
302
304{
305#ifdef G4VTKDEBUG
306 G4cout << "G4VtkViewer::ShowView()" << G4endl;
307#endif
308
309 infoTextActor->GetTextProperty()->SetFontSize(28);
310 G4Colour colour = fVP.GetBackgroundColour();
311
312 // make sure text is always visible
313 infoTextActor->GetTextProperty()->SetColor(std::fmod(colour.GetRed() + 0.5, 1.0),
314 std::fmod(colour.GetGreen() + 0.5, 1.0),
315 std::fmod(colour.GetBlue() + 0.5, 1.0));
316 infoTextActor->GetTextProperty()->SetFontSize(20);
317 infoCallback->SetTextActor(infoTextActor);
318 renderer->AddObserver(vtkCommand::EndEvent, infoCallback);
319 geant4Callback->SetGeant4ViewParameters(&fVP);
320 renderer->AddObserver(vtkCommand::EndEvent, geant4Callback);
321}
322
324{
325#ifdef G4VTKDEBUG
326 G4cout << "G4VtkViewer::FinishView()" << G4endl;
327#endif
328
329 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
330 fVtkSceneHandler.Modified();
331
332 _renderWindow->GetInteractor()->Initialize();
333 _renderWindow->Render();
334
335 if (firstFinishView) {
336 firstFinishView = false;
337 }
338 else {
339 _renderWindow->GetInteractor()->Start();
340 }
341}
342
344{
345 vtkImageWriter* imWriter = nullptr;
346
347 if (format == "bmp") {
348 imWriter = vtkBMPWriter::New();
349 }
350 else if (format == "jpg") {
351 imWriter = vtkJPEGWriter::New();
352 }
353 else if (format == "pnm") {
354 imWriter = vtkPNMWriter::New();
355 }
356 else if (format == "png") {
357 imWriter = vtkPNGWriter::New();
358 }
359 else if (format == "tiff") {
360 imWriter = vtkTIFFWriter::New();
361 }
362 else if (format == "ps") {
363 imWriter = vtkPostScriptWriter::New();
364 }
365 else {
366 return;
367 }
368
369 _renderWindow->Render();
370
373 winToImage->SetInput(_renderWindow);
374 winToImage->SetScale(1);
375 if (format == "ps") {
376 winToImage->SetInputBufferTypeToRGB();
377 winToImage->ReadFrontBufferOff();
378 winToImage->Update();
379 }
380 else {
381 winToImage->SetInputBufferTypeToRGBA();
382 }
383
384 imWriter->SetFileName((path + "." + format).c_str());
385 imWriter->SetInputConnection(winToImage->GetOutputPort());
386 imWriter->Write();
387}
388
390{
392 exporter->SetRenderWindow(_renderWindow);
393 exporter->SetFilePrefix(path.c_str());
394 exporter->Write();
395}
396
398{
400 exporter->SetRenderWindow(_renderWindow);
401 exporter->SetFileName((path + ".vrml").c_str());
402 exporter->Write();
403}
404
406{
408 exporter->SetRenderWindow(_renderWindow);
409 exporter->SetFileName((path + ".vtp").c_str());
410 exporter->Write();
411}
412
415 exporter->SetRenderWindow(_renderWindow);
416 exporter->SetFileName((fileName+".gltf").c_str());
417 exporter->InlineDataOn();
418 exporter->Write();
419}
420
423 exporter->SetRenderWindow(_renderWindow);
424 exporter->SetFileName((fileName+".x3d").c_str());
425 exporter->Write();
426}
427
428
436
438{
439 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
440 G4VtkStore& s = fVtkSceneHandler.GetStore();
441
442 // create new renderer
443 vtkNew<vtkRenderer> tempRenderer;
444
445 // loop over pipelines
446 auto separate = s.GetSeparatePipeMap();
447 for (const auto& i : separate) {
448 i.second->GetActor();
449 auto children = i.second->GetChildPipelines();
450 for (auto child : children) {
451 if (child->GetTypeName() == "G4VtkCutterPipeline") {
452 auto childCutter = dynamic_cast<G4VtkCutterPipeline*>(child);
453 tempRenderer->AddActor(childCutter->GetActor());
454 }
455 }
456 }
457
458 auto tensor = s.GetTensorPipeMap();
459 for (const auto& i : tensor) {
460 i.second->GetActor();
461 auto children = i.second->GetChildPipelines();
462 for (auto child : children) {
463 if (child->GetTypeName() == "G4VtkCutterPipeline") {
464 auto childCutter = dynamic_cast<G4VtkCutterPipeline*>(child);
465 tempRenderer->AddActor(childCutter->GetActor());
466 }
467 }
468 }
469
470 auto append = s.GetAppendPipeMap();
471 for (const auto& i : append) {
472 i.second->GetActor();
473 auto children = i.second->GetChildPipelines();
474 for (auto child : children) {
475 if (child->GetTypeName() == "G4VtkCutterPipeline") {
476 auto childCutter = dynamic_cast<G4VtkCutterPipeline*>(child);
477 tempRenderer->AddActor(childCutter->GetActor());
478 }
479 }
480 }
481
482 auto baked = s.GetBakePipeMap();
483 for (const auto& i : baked) {
484 i.second->GetActor();
485 auto children = i.second->GetChildPipelines();
486 for (auto child : children) {
487 if (child->GetTypeName() == "G4VtkCutterPipeline") {
488 auto childCutter = dynamic_cast<G4VtkCutterPipeline*>(child);
489 tempRenderer->AddActor(childCutter->GetActor());
490 }
491 }
492 }
493
494 vtkNew<vtkRenderWindow> tempRenderWindow;
495 tempRenderWindow->AddRenderer(tempRenderer);
496 vtkNew<vtkSingleVTPExporter> exporter;
497 exporter->SetRenderWindow(tempRenderWindow);
498 exporter->SetFileName(fileName.c_str());
499 exporter->Write();
500}
501
503{
504 vtkSmartPointer<vtkRenderWindow> tempRenderWindow;
505 vtkNew<vtkRenderer> tempRenderer;
506 tempRenderWindow = vtkSmartPointer<vtkRenderWindow>::New();
507 tempRenderWindow->AddRenderer(tempRenderer);
508
509 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
510
511 if (storeName == "transient") {
512 G4VtkStore& store = fVtkSceneHandler.GetTransientStore();
513 store.AddToRenderer(tempRenderer);
514 }
515 else {
516 G4VtkStore& store = fVtkSceneHandler.GetStore();
517 store.AddToRenderer(tempRenderer);
518 }
519
520 if (fileName.find("obj") != std::string::npos) {
521 vtkNew<vtkOBJExporter> exporter;
522 exporter->SetRenderWindow(tempRenderWindow);
523 exporter->SetFilePrefix(fileName.c_str());
524 exporter->Write();
525 }
526 else if (fileName.find("vrml") != std::string::npos) {
527 vtkNew<vtkVRMLExporter> exporter;
528 exporter->SetRenderWindow(tempRenderWindow);
529 exporter->SetFileName(fileName.c_str());
530 exporter->Write();
531 }
532 else if (fileName.find("vtp") != std::string::npos) {
533 vtkNew<vtkSingleVTPExporter> exporter;
534 exporter->SetRenderWindow(tempRenderWindow);
535 exporter->SetFileName(fileName.c_str());
536 exporter->Write();
537 }
538 else if (fileName.find("gltf") != std::string::npos) {
539 vtkNew<vtkGLTFExporter> exporter;
540 exporter->SetRenderWindow(tempRenderWindow);
541 exporter->SetFileName(fileName.c_str());
542 exporter->InlineDataOn();
543 exporter->Write();
544 }
545}
546
548{
549 // make sure text is always visible
550 G4Colour colour = fVP.GetBackgroundColour();
551 infoTextActor->GetTextProperty()->SetColor(std::fmod(colour.GetRed() + 0.5, 1.0),
552 std::fmod(colour.GetGreen() + 0.5, 1.0),
553 std::fmod(colour.GetBlue() + 0.5, 1.0));
554 infoTextActor->GetTextProperty()->SetFontSize(20);
555 infoCallback->SetTextActor(infoTextActor);
556 renderer->AddObserver(vtkCommand::EndEvent, infoCallback);
557 renderer->AddActor(infoTextActor);
558 infoTextActor->SetVisibility(0);
559}
560
562{
563 vtkNew<vtkIPWCallback> clipperCallback;
564 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
565 G4VtkStore& store = fVtkSceneHandler.GetStore();
566 clipperCallback->SetStore(&store);
567 clipperCallback->SetUpdatePipelineName("clipper", "clipper");
568
569 G4double bounds[6];
570 store.GetBounds(bounds);
571 auto vplane = G4Plane3DToVtkPlane(plane);
572 clipperPlaneRepresentation->SetPlaceFactor(
573 1.25); // This must be set prior to placing the widget.
574 clipperPlaneRepresentation->PlaceWidget(bounds);
575 clipperPlaneRepresentation->SetNormal(vplane->GetNormal());
576
577 vtkNew<vtkPropCollection> planeRepActors;
578 clipperPlaneRepresentation->GetActors(planeRepActors);
579 planeRepActors->InitTraversal();
580
583 clipperPlaneWidget->AddObserver(vtkCommand::InteractionEvent, clipperCallback);
584
585 clipperPlaneWidget->SetEnabled(0);
586}
587
589{
590 vtkNew<vtkIPWCallback> cutterCallback;
591 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
592 G4VtkStore& store = fVtkSceneHandler.GetStore();
593 cutterCallback->SetStore(&store);
594 cutterCallback->SetUpdatePipelineName("cutter", "cutter");
595
596 G4double bounds[6];
597 store.GetBounds(bounds);
598 auto vplane = G4Plane3DToVtkPlane(plane);
599 cutterPlaneRepresentation->SetPlaceFactor(1.25); // This must be set prior to placing the widget.
600 cutterPlaneRepresentation->PlaceWidget(bounds);
601 cutterPlaneRepresentation->SetNormal(vplane->GetNormal());
602
605 cutterPlaneWidget->AddObserver(vtkCommand::InteractionEvent, cutterCallback);
606
607 cutterPlaneWidget->SetEnabled(0);
608}
609
611{
612 renderer->SetUseShadows(1);
613}
614
616{
617 renderer->SetUseShadows(0);
618}
619
621{
622 infoTextActor->SetVisibility(1);
623}
624
626{
627 infoTextActor->SetVisibility(0);
628}
629
630void G4VtkViewer::EnableClipper(const G4Plane3D& plane, G4bool bWidget)
631{
632 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
633 G4VtkStore& s = fVtkSceneHandler.GetStore();
634 G4String name = G4String("clipper");
635 s.AddClipper(name, plane);
636 if (bWidget) {
638 }
639}
640
642{
643 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
644 G4VtkStore& s = fVtkSceneHandler.GetStore();
645 s.RemoveClipper("clipper");
646}
647
649{
650 clipperPlaneWidget->SetEnabled(1);
651}
652
654{
655 clipperPlaneWidget->SetEnabled(0);
656}
657
658void G4VtkViewer::EnableCutter(const G4Plane3D& plane, G4bool bWidget)
659{
660 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
661 G4VtkStore& s = fVtkSceneHandler.GetStore();
662 G4String name = G4String("cutter");
663 s.AddCutter(name, plane);
664 if (bWidget) {
666 }
667}
668
670{
671 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
672 G4VtkStore& s = fVtkSceneHandler.GetStore();
673 s.RemoveCutter("cutter");
674}
675
677{
678 G4cout << "enable cutter widget" << G4endl;
679 cutterPlaneWidget->SetEnabled(1);
680}
681
683{
684 cutterPlaneWidget->SetEnabled(0);
685}
686
688{
689 camOrientWidget->SetParentRenderer(renderer);
690 // Enable the widget.
691 camOrientWidget->Off();
692}
693
698
703
704void G4VtkViewer::AddImageOverlay(const G4String& fileName, const G4double alpha,
705 const G4double imageBottomLeft[2],
706 const G4double worldBottomLeft[2],
707 const G4double imageTopRight[2], const G4double worldTopRight[2],
708 const G4double rotation[3], const G4double translation[3])
709{
710 auto xScale = (worldTopRight[0] - worldBottomLeft[0]) / (imageTopRight[0] - imageBottomLeft[0]);
711 auto yScale = -(worldTopRight[1] - worldBottomLeft[1]) / (imageTopRight[1] - imageBottomLeft[1]);
712
713 G4cout << xScale << " " << yScale << G4endl;
714 auto transformation = G4Transform3D::Identity;
715 auto scal = G4Scale3D(xScale, yScale, 1);
716 auto rotx = G4RotateX3D(rotation[0]/180*CLHEP::pi);
717 auto roty = G4RotateY3D(rotation[1]/180*CLHEP::pi);
718 auto rotz = G4RotateZ3D(rotation[2]/180*CLHEP::pi);
719 auto tranImg = G4Translate3D( -std::fabs(imageBottomLeft[0] + imageTopRight[0]) / 2.0,
720 -std::fabs(imageBottomLeft[1] + imageTopRight[1]) / 2.0,
721 0);
722 auto tran = G4Translate3D(translation[0],
723 translation[1],
724 translation[2]);
725
726 G4cout << translation[0] << " " << translation[1] << " " << translation[2] << G4endl;
727 transformation = tran * rotz * roty * rotx * scal * tranImg * transformation;
728
729 G4cout << transformation.dx() << " " << transformation.dy() << " " << transformation.dz() << G4endl;
730 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
731 G4VtkStore& st = fVtkSceneHandler.GetTransientStore();
732
733 G4VtkVisContext vc = G4VtkVisContext(this, nullptr, false, transformation);
734 vc.alpha = alpha;
735 st.AddNonG4ObjectImage(fileName, vc);
736}
737
738void G4VtkViewer::AddGeometryOverlay(const G4String& fileName, const G4double colour[3],
739 const G4double alpha, const G4String& representation,
740 const G4double scale[3], const G4double rotation[3],
741 const G4double translation[3])
742{
743 auto transformation = G4Transform3D::Identity;
744 auto scal = G4Scale3D(scale[0], scale[1], scale[2]);
745 auto rotx = G4RotateX3D(rotation[0]);
746 auto roty = G4RotateY3D(rotation[1]);
747 auto rotz = G4RotateZ3D(rotation[2]);
748 auto tran = G4Translate3D(translation[0], translation[1], translation[2]);
749
750 transformation = tran * rotz * roty * rotx * scal * transformation;
751
752 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
753 G4VtkStore& st = fVtkSceneHandler.GetTransientStore();
754
755
756 G4VtkVisContext vc = G4VtkVisContext(this, nullptr, false, transformation);
757 if (representation == "w")
759 else if (representation == "s")
761 vc.alpha = alpha;
762 vc.red = colour[0];
763 vc.green = colour[1];
764 vc.blue = colour[2];
765 st.AddNonG4ObjectPolydata(fileName, vc);
766}
767
769{
770 cutterPlaneRepresentation->VisibilityOff();
771
772 G4cout << "Number of VTK props> " << renderer->GetNumberOfPropsRendered() << G4endl;
773 G4cout << "Number of VTK actors> " << renderer->GetActors()->GetNumberOfItems() << G4endl;
774 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
775 G4VtkStore& s = fVtkSceneHandler.GetStore();
776 G4VtkStore& st = fVtkSceneHandler.GetTransientStore();
777 s.Print();
778 st.Print();
779}
780
782{
783 // Get the scene handler
784 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
785 fVtkSceneHandler.SetPolyhedronPipeline(type);
786}
787
788void G4VtkViewer::SetWidgetInteractor(vtkAbstractWidget* widget)
789{
790 widget->SetInteractor(_renderWindow->GetInteractor());
791}
HepGeom::Plane3D< G4double > G4Plane3D
Definition G4Plane3D.hh:34
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::RotateX3D G4RotateX3D
HepGeom::RotateY3D G4RotateY3D
HepGeom::Translate3D G4Translate3D
HepGeom::Scale3D G4Scale3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
vtkSmartPointer< vtkPlane > G4Plane3DToVtkPlane(const G4Plane3D &g4plane)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetBlue() const
Definition G4Colour.hh:172
G4double GetRed() const
Definition G4Colour.hh:170
G4double GetGreen() const
Definition G4Colour.hh:171
void ProcessView()
Definition G4VViewer.cc:109
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:253
void NeedKernelVisit()
Definition G4VViewer.cc:82
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:258
G4ViewParameters fVP
Definition G4VViewer.hh:257
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition G4VViewer.cc:49
void SetPolyhedronPipeline(const G4String &str)
std::map< std::size_t, std::shared_ptr< G4VtkPolydataPipeline > > & GetSeparatePipeMap()
void AddNonG4ObjectImage(const G4String &fileName, const G4VtkVisContext &vc)
void AddToRenderer(vtkRenderer *renderer)
void RemoveCutter(G4String name)
void GetBounds(G4double maxBound[6])
void Clear()
void Print()
void RemoveClipper(G4String name)
void AddNonG4ObjectPolydata(const G4String fileName, const G4VtkVisContext &vc)
vtkNew< vtkCamera > camera
void DisableClipper()
vtkNew< vtkCameraOrientationWidget > camOrientWidget
void DisableCutter(G4String name)
void ExportVRMLScene(G4String)
void SetView() override
void AddImageOverlay(const G4String &fileName, const G4double alpha, const G4double imageBottomLeft[2], const G4double worldBottomLeft[2], const G4double imageTopRight[2], const G4double worldTopRight[2], const G4double rot[3], const G4double trans[3])
vtkNew< vtkImplicitPlaneWidget2 > cutterPlaneWidget
vtkNew< vtkGeant4Callback > geant4Callback
void Initialise() override
G4bool firstFinishView
void EnableShadows()
void ShowView() override
void SetPolyhedronPipeline(const G4String &t)
void ExportScreenShot(G4String, G4String)
void DrawShadows()
void FinishView() override
void ExportGLTFScene(G4String)
~G4VtkViewer() override
vtkNew< vtkTextActor > infoTextActor
virtual void SetWidgetInteractor(vtkAbstractWidget *widget)
virtual void AddClipperPlaneWidget(const G4Plane3D &plane)
vtkRenderWindowInteractor * renderWindowInteractor
vtkNew< vtkImplicitPlaneRepresentation > cutterPlaneRepresentation
void DisableShadows()
void ExportJSONRenderWindowScene(G4String)
virtual void EnableClipperWidget()
void DisableHUD()
G4VtkViewer(G4VSceneHandler &, const G4String &name)
virtual void AddCameraOrientationWidget()
void ExportFormatStore(G4String fileName, G4String store)
void EnableClipper(const G4Plane3D &plane, G4bool widget)
virtual void AddCutterPlaneWidget(const G4Plane3D &plane)
virtual void EnableCameraOrientationWidget()
vtkNew< vtkInfoCallback > infoCallback
virtual void DisableCutterWidget()
virtual void DisableCameraOrientationWidget()
vtkNew< vtkImplicitPlaneWidget2 > clipperPlaneWidget
void EnableHUD()
vtkNew< vtkRenderer > renderer
G4double cameraDistance
void AddViewHUD()
void EnableCutter(const G4Plane3D &plane, G4bool bWidget)
G4bool firstSetView
void ExportVTPScene(G4String)
vtkNew< vtkImplicitPlaneRepresentation > clipperPlaneRepresentation
virtual void DisableClipperWidget()
void ExportOBJScene(G4String)
void ExportX3DScene(G4String)
void ClearView() override
virtual void EnableCutterWidget()
void AddGeometryOverlay(const G4String &fileName, const G4double colour[3], const G4double alpha, const G4String &representation, const G4double scale[3], const G4double rotation[3], const G4double translation[3])
vtkRenderWindow * _renderWindow
void ExportVTPCutter(G4String fileName)
void DrawView() override
G4ViewParameters::DrawingStyle fDrawingStyle
BasicVector3D< T > unit() const
static DLL_API const Transform3D Identity