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