Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkViewer Class Reference

#include <G4VtkViewer.hh>

+ Inheritance diagram for G4VtkViewer:

Public Member Functions

 G4VtkViewer (G4VSceneHandler &, const G4String &name)
 
void Initialise () override
 
 ~G4VtkViewer () override
 
void SetView () override
 
void ClearView () override
 
void DrawView () override
 
void ShowView () override
 
void FinishView () override
 
void ExportScreenShot (G4String, G4String)
 
void ExportOBJScene (G4String)
 
void ExportVRMLScene (G4String)
 
void ExportVTPScene (G4String)
 
void ExportGLTFScene (G4String)
 
void ExportX3DScene (G4String)
 
void ExportJSONRenderWindowScene (G4String)
 
void ExportVTPCutter (G4String fileName)
 
void ExportFormatStore (G4String fileName, G4String store)
 
void DrawShadows ()
 
void EnableShadows ()
 
void DisableShadows ()
 
void AddViewHUD ()
 
void EnableHUD ()
 
void DisableHUD ()
 
virtual void AddClipperPlaneWidget (const G4Plane3D &plane)
 
void EnableClipper (const G4Plane3D &plane, G4bool widget)
 
void DisableClipper ()
 
virtual void EnableClipperWidget ()
 
virtual void DisableClipperWidget ()
 
virtual void AddCutterPlaneWidget (const G4Plane3D &plane)
 
void EnableCutter (const G4Plane3D &plane, G4bool bWidget)
 
void DisableCutter (G4String name)
 
virtual void EnableCutterWidget ()
 
virtual void DisableCutterWidget ()
 
virtual void AddCameraOrientationWidget ()
 
virtual void EnableCameraOrientationWidget ()
 
virtual void DisableCameraOrientationWidget ()
 
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])
 
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])
 
void Render ()
 
void StartInteractor ()
 
void Print ()
 
void SetPolyhedronPipeline (const G4String &t)
 
virtual void SetWidgetInteractor (vtkAbstractWidget *widget)
 
void ExportView ()
 
void SetGeant4View ()
 
- Public Member Functions inherited from G4VViewer
 G4VViewer (G4VSceneHandler &, G4int id, const G4String &name="")
 
virtual ~G4VViewer ()
 
virtual void ResetView ()
 
void RefreshView ()
 
virtual G4bool ReadyToDraw ()
 
std::vector< G4ThreeVectorComputeFlyThrough (G4Vector3D *)
 
virtual void DoneWithMasterThread ()
 
virtual void MovingToVisSubThread ()
 
virtual void SwitchToVisSubThread ()
 
virtual void DoneWithVisSubThread ()
 
virtual void MovingToMasterThread ()
 
virtual void SwitchToMasterThread ()
 
void InsertModelInSceneTree (G4VModel *)
 
const G4SceneTreeItemGetSceneTree ()
 
G4SceneTreeItemAccessSceneTree ()
 
void UpdateGUISceneTree ()
 
const G4StringGetName () const
 
const G4StringGetShortName () const
 
void SetName (const G4String &)
 
G4int GetViewId () const
 
G4VSceneHandlerGetSceneHandler () const
 
const G4ViewParametersGetViewParameters () const
 
const G4ViewParametersGetDefaultViewParameters () const
 
G4double GetKernelVisitElapsedTimeSeconds () const
 
virtual const std::vector< G4ModelingParameters::VisAttributesModifier > * GetPrivateVisAttributesModifiers () const
 
void SetViewParameters (const G4ViewParameters &vp)
 
void SetDefaultViewParameters (const G4ViewParameters &vp)
 
const G4VisAttributesGetApplicableVisAttributes (const G4VisAttributes *) const
 
void SetNeedKernelVisit (G4bool need)
 
void NeedKernelVisit ()
 
void ProcessView ()
 

Public Attributes

vtkNew< vtkTextActor > infoTextActor
 
vtkNew< vtkInfoCallbackinfoCallback
 
vtkNew< vtkGeant4Callbackgeant4Callback
 
vtkSmartPointer< vtkLight > light
 
vtkNew< vtkCamera > camera
 
vtkNew< vtkRenderer > renderer
 
vtkRenderWindow * _renderWindow
 
vtkRenderWindowInteractor * renderWindowInteractor
 
- Public Attributes inherited from G4VViewer
const G4int fMaxNTouchables = 10000
 
G4bool fCurtailDescent = false
 

Protected Attributes

G4bool firstSetView = true
 
G4bool firstFinishView = true
 
G4double cameraDistance
 
vtkNew< vtkImplicitPlaneRepresentation > cutterPlaneRepresentation
 
vtkNew< vtkImplicitPlaneWidget2 > cutterPlaneWidget
 
vtkNew< vtkImplicitPlaneRepresentation > clipperPlaneRepresentation
 
vtkNew< vtkImplicitPlaneWidget2 > clipperPlaneWidget
 
vtkNew< vtkCameraOrientationWidget > camOrientWidget
 
bool bCutter = false
 
bool bClipper = false
 
bool bHud = false
 
bool bOrientation = false
 
- Protected Attributes inherited from G4VViewer
G4VSceneHandlerfSceneHandler
 
G4int fViewId
 
G4String fName
 
G4String fShortName
 
G4ViewParameters fVP
 
G4ViewParameters fDefaultVP
 
G4double fKernelVisitElapsedTimeSeconds = 999.
 
G4SceneTreeItem fSceneTree
 
G4bool fNeedKernelVisit
 

Additional Inherited Members

- Protected Member Functions inherited from G4VViewer
void SetTouchable (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath)
 
void TouchableSetVisibility (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, G4bool visibility)
 
void TouchableSetColour (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, const G4Colour &)
 

Detailed Description

Definition at line 214 of file G4VtkViewer.hh.

Constructor & Destructor Documentation

◆ G4VtkViewer()

G4VtkViewer::G4VtkViewer ( G4VSceneHandler & sceneHandler,
const G4String & name )

Definition at line 82 of file G4VtkViewer.cc.

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}
G4int IncrementViewCount()
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:258
G4ViewParameters fVP
Definition G4VViewer.hh:257
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition G4VViewer.cc:49

Referenced by G4VtkOffscreenViewer::G4VtkOffscreenViewer(), and G4VtkQtViewer::G4VtkQtViewer().

◆ ~G4VtkViewer()

G4VtkViewer::~G4VtkViewer ( )
override

Definition at line 150 of file G4VtkViewer.cc.

151{
152#ifdef G4VTKDEBUG
153 G4cout << "G4VtkViewer::~G4VtkViewer()" << G4endl;
154#endif
155}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

Member Function Documentation

◆ AddCameraOrientationWidget()

void G4VtkViewer::AddCameraOrientationWidget ( )
virtual

Definition at line 687 of file G4VtkViewer.cc.

688{
689 camOrientWidget->SetParentRenderer(renderer);
690 // Enable the widget.
691 camOrientWidget->Off();
692}
vtkNew< vtkCameraOrientationWidget > camOrientWidget
vtkNew< vtkRenderer > renderer

Referenced by DrawView().

◆ AddClipperPlaneWidget()

void G4VtkViewer::AddClipperPlaneWidget ( const G4Plane3D & plane)
virtual

Definition at line 561 of file G4VtkViewer.cc.

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}
double G4double
Definition G4Types.hh:83
vtkSmartPointer< vtkPlane > G4Plane3DToVtkPlane(const G4Plane3D &g4plane)
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:253
void GetBounds(G4double maxBound[6])
virtual void SetWidgetInteractor(vtkAbstractWidget *widget)
vtkNew< vtkImplicitPlaneWidget2 > clipperPlaneWidget
vtkNew< vtkImplicitPlaneRepresentation > clipperPlaneRepresentation

Referenced by DrawView().

◆ AddCutterPlaneWidget()

void G4VtkViewer::AddCutterPlaneWidget ( const G4Plane3D & plane)
virtual

Definition at line 588 of file G4VtkViewer.cc.

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}
vtkNew< vtkImplicitPlaneWidget2 > cutterPlaneWidget
vtkNew< vtkImplicitPlaneRepresentation > cutterPlaneRepresentation

Referenced by DrawView().

◆ AddGeometryOverlay()

void G4VtkViewer::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] )

Definition at line 738 of file G4VtkViewer.cc.

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}
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::RotateX3D G4RotateX3D
HepGeom::RotateY3D G4RotateY3D
HepGeom::Translate3D G4Translate3D
HepGeom::Scale3D G4Scale3D
void AddNonG4ObjectPolydata(const G4String fileName, const G4VtkVisContext &vc)
G4ViewParameters::DrawingStyle fDrawingStyle
static DLL_API const Transform3D Identity

◆ AddImageOverlay()

void G4VtkViewer::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] )

Definition at line 704 of file G4VtkViewer.cc.

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}
void AddNonG4ObjectImage(const G4String &fileName, const G4VtkVisContext &vc)

◆ AddViewHUD()

void G4VtkViewer::AddViewHUD ( )

Definition at line 547 of file G4VtkViewer.cc.

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}
G4double GetBlue() const
Definition G4Colour.hh:172
G4double GetRed() const
Definition G4Colour.hh:170
G4double GetGreen() const
Definition G4Colour.hh:171
vtkNew< vtkTextActor > infoTextActor
vtkNew< vtkInfoCallback > infoCallback

Referenced by DrawView().

◆ ClearView()

void G4VtkViewer::ClearView ( void )
overridevirtual

Implements G4VViewer.

Definition at line 235 of file G4VtkViewer.cc.

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}
void Clear()

◆ DisableCameraOrientationWidget()

void G4VtkViewer::DisableCameraOrientationWidget ( )
virtual

Definition at line 699 of file G4VtkViewer.cc.

700{
701 camOrientWidget->Off();
702}

◆ DisableClipper()

void G4VtkViewer::DisableClipper ( )

Definition at line 641 of file G4VtkViewer.cc.

642{
643 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
644 G4VtkStore& s = fVtkSceneHandler.GetStore();
645 s.RemoveClipper("clipper");
646}
void RemoveClipper(G4String name)

◆ DisableClipperWidget()

void G4VtkViewer::DisableClipperWidget ( )
virtual

Definition at line 653 of file G4VtkViewer.cc.

654{
655 clipperPlaneWidget->SetEnabled(0);
656}

◆ DisableCutter()

void G4VtkViewer::DisableCutter ( G4String name)

Definition at line 669 of file G4VtkViewer.cc.

670{
671 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
672 G4VtkStore& s = fVtkSceneHandler.GetStore();
673 s.RemoveCutter("cutter");
674}
void RemoveCutter(G4String name)

◆ DisableCutterWidget()

void G4VtkViewer::DisableCutterWidget ( )
virtual

Definition at line 682 of file G4VtkViewer.cc.

683{
684 cutterPlaneWidget->SetEnabled(0);
685}

◆ DisableHUD()

void G4VtkViewer::DisableHUD ( )

Definition at line 625 of file G4VtkViewer.cc.

626{
627 infoTextActor->SetVisibility(0);
628}

◆ DisableShadows()

void G4VtkViewer::DisableShadows ( )

Definition at line 615 of file G4VtkViewer.cc.

616{
617 renderer->SetUseShadows(0);
618}

◆ DrawShadows()

void G4VtkViewer::DrawShadows ( )

Definition at line 283 of file G4VtkViewer.cc.

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}
vtkRenderWindow * _renderWindow

◆ DrawView()

void G4VtkViewer::DrawView ( )
overridevirtual

Implements G4VViewer.

Definition at line 249 of file G4VtkViewer.cc.

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}
HepGeom::Plane3D< G4double > G4Plane3D
Definition G4Plane3D.hh:34
void ProcessView()
Definition G4VViewer.cc:109
void NeedKernelVisit()
Definition G4VViewer.cc:82
void FinishView() override
virtual void AddClipperPlaneWidget(const G4Plane3D &plane)
virtual void AddCameraOrientationWidget()
virtual void AddCutterPlaneWidget(const G4Plane3D &plane)
void AddViewHUD()

◆ EnableCameraOrientationWidget()

void G4VtkViewer::EnableCameraOrientationWidget ( )
virtual

Definition at line 694 of file G4VtkViewer.cc.

695{
696 camOrientWidget->On();
697}

◆ EnableClipper()

void G4VtkViewer::EnableClipper ( const G4Plane3D & plane,
G4bool widget )

Definition at line 630 of file G4VtkViewer.cc.

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}
void AddClipper(G4String name, const G4Plane3D &plane)
virtual void EnableClipperWidget()
const char * name(G4int ptype)

◆ EnableClipperWidget()

void G4VtkViewer::EnableClipperWidget ( )
virtual

Reimplemented in G4VtkQtViewer.

Definition at line 648 of file G4VtkViewer.cc.

649{
650 clipperPlaneWidget->SetEnabled(1);
651}

Referenced by EnableClipper(), and G4VtkQtViewer::EnableClipperWidget().

◆ EnableCutter()

void G4VtkViewer::EnableCutter ( const G4Plane3D & plane,
G4bool bWidget )

Definition at line 658 of file G4VtkViewer.cc.

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}
void AddCutter(G4String name, const G4Plane3D &plane)
virtual void EnableCutterWidget()

◆ EnableCutterWidget()

void G4VtkViewer::EnableCutterWidget ( )
virtual

Definition at line 676 of file G4VtkViewer.cc.

677{
678 G4cout << "enable cutter widget" << G4endl;
679 cutterPlaneWidget->SetEnabled(1);
680}

Referenced by EnableCutter().

◆ EnableHUD()

void G4VtkViewer::EnableHUD ( )

Definition at line 620 of file G4VtkViewer.cc.

621{
622 infoTextActor->SetVisibility(1);
623}

◆ EnableShadows()

void G4VtkViewer::EnableShadows ( )

Definition at line 610 of file G4VtkViewer.cc.

611{
612 renderer->SetUseShadows(1);
613}

◆ ExportFormatStore()

void G4VtkViewer::ExportFormatStore ( G4String fileName,
G4String store )

Definition at line 502 of file G4VtkViewer.cc.

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}
void AddToRenderer(vtkRenderer *renderer)

◆ ExportGLTFScene()

void G4VtkViewer::ExportGLTFScene ( G4String fileName)

Definition at line 413 of file G4VtkViewer.cc.

413 {
414 vtkSmartPointer<vtkGLTFExporter> exporter = vtkSmartPointer<vtkGLTFExporter>::New();
415 exporter->SetRenderWindow(_renderWindow);
416 exporter->SetFileName((fileName+".gltf").c_str());
417 exporter->InlineDataOn();
418 exporter->Write();
419}

◆ ExportJSONRenderWindowScene()

void G4VtkViewer::ExportJSONRenderWindowScene ( G4String )

Definition at line 429 of file G4VtkViewer.cc.

429 {
430 vtkSmartPointer<vtkJSONRenderWindowExporter> exporter = vtkSmartPointer<vtkJSONRenderWindowExporter>::New();
431 vtkSmartPointer<vtkVtkJSSceneGraphSerializer> serializer = vtkSmartPointer<vtkVtkJSSceneGraphSerializer>::New();
432 exporter->SetRenderWindow(_renderWindow);
433 exporter->SetSerializer(serializer);
434 exporter->Write();
435}

◆ ExportOBJScene()

void G4VtkViewer::ExportOBJScene ( G4String path)

Definition at line 389 of file G4VtkViewer.cc.

390{
391 vtkSmartPointer<vtkOBJExporter> exporter = vtkSmartPointer<vtkOBJExporter>::New();
392 exporter->SetRenderWindow(_renderWindow);
393 exporter->SetFilePrefix(path.c_str());
394 exporter->Write();
395}

◆ ExportScreenShot()

void G4VtkViewer::ExportScreenShot ( G4String path,
G4String format )

Definition at line 343 of file G4VtkViewer.cc.

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
371 vtkSmartPointer<vtkWindowToImageFilter> winToImage =
372 vtkSmartPointer<vtkWindowToImageFilter>::New();
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}

◆ ExportView()

void G4VtkViewer::ExportView ( )
inline

Definition at line 282 of file G4VtkViewer.hh.

282{};

◆ ExportVRMLScene()

void G4VtkViewer::ExportVRMLScene ( G4String path)

Definition at line 397 of file G4VtkViewer.cc.

398{
399 vtkSmartPointer<vtkVRMLExporter> exporter = vtkSmartPointer<vtkVRMLExporter>::New();
400 exporter->SetRenderWindow(_renderWindow);
401 exporter->SetFileName((path + ".vrml").c_str());
402 exporter->Write();
403}

◆ ExportVTPCutter()

void G4VtkViewer::ExportVTPCutter ( G4String fileName)

Definition at line 437 of file G4VtkViewer.cc.

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}
std::map< std::size_t, std::shared_ptr< G4VtkPolydataPipeline > > & GetSeparatePipeMap()
std::map< std::size_t, std::shared_ptr< G4VtkPolydataInstanceTensorPipeline > > & GetTensorPipeMap()
std::map< std::size_t, std::shared_ptr< G4VtkPolydataInstanceBakePipeline > > & GetBakePipeMap()
std::map< std::size_t, std::shared_ptr< G4VtkPolydataInstanceAppendPipeline > > & GetAppendPipeMap()

◆ ExportVTPScene()

void G4VtkViewer::ExportVTPScene ( G4String path)

Definition at line 405 of file G4VtkViewer.cc.

406{
407 vtkSmartPointer<vtkSingleVTPExporter> exporter = vtkSmartPointer<vtkSingleVTPExporter>::New();
408 exporter->SetRenderWindow(_renderWindow);
409 exporter->SetFileName((path + ".vtp").c_str());
410 exporter->Write();
411}

◆ ExportX3DScene()

void G4VtkViewer::ExportX3DScene ( G4String fileName)

Definition at line 421 of file G4VtkViewer.cc.

421 {
422 vtkSmartPointer<vtkX3DExporter> exporter = vtkSmartPointer<vtkX3DExporter>::New();
423 exporter->SetRenderWindow(_renderWindow);
424 exporter->SetFileName((fileName+".x3d").c_str());
425 exporter->Write();
426}

◆ FinishView()

void G4VtkViewer::FinishView ( void )
overridevirtual

Reimplemented from G4VViewer.

Definition at line 323 of file G4VtkViewer.cc.

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}
G4bool firstFinishView

Referenced by DrawView().

◆ Initialise()

void G4VtkViewer::Initialise ( )
overridevirtual

Reimplemented from G4VViewer.

Definition at line 92 of file G4VtkViewer.cc.

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
146 vtkSmartPointer<G4VtkInteractorStyle> style = vtkSmartPointer<G4VtkInteractorStyle>::New();
147 renderWindowInteractor->SetInteractorStyle(style);
148}
int G4int
Definition G4Types.hh:85
vtkNew< vtkCamera > camera
vtkNew< vtkGeant4Callback > geant4Callback
vtkRenderWindowInteractor * renderWindowInteractor

Referenced by G4VtkOffscreenViewer::Initialise().

◆ Print()

void G4VtkViewer::Print ( )

Definition at line 768 of file G4VtkViewer.cc.

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}
void Print()

◆ Render()

void G4VtkViewer::Render ( )
inline

Definition at line 271 of file G4VtkViewer.hh.

271{_renderWindow->Render();}

◆ SetGeant4View()

void G4VtkViewer::SetGeant4View ( )
inline

Definition at line 283 of file G4VtkViewer.hh.

283{};

◆ SetPolyhedronPipeline()

void G4VtkViewer::SetPolyhedronPipeline ( const G4String & t)

Definition at line 781 of file G4VtkViewer.cc.

782{
783 // Get the scene handler
784 auto& fVtkSceneHandler = dynamic_cast<G4VtkSceneHandler&>(fSceneHandler);
785 fVtkSceneHandler.SetPolyhedronPipeline(type);
786}

◆ SetView()

void G4VtkViewer::SetView ( )
overridevirtual

Implements G4VViewer.

Definition at line 157 of file G4VtkViewer.cc.

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}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
bool G4bool
Definition G4Types.hh:86
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
G4double cameraDistance
G4bool firstSetView
BasicVector3D< T > unit() const

◆ SetWidgetInteractor()

void G4VtkViewer::SetWidgetInteractor ( vtkAbstractWidget * widget)
virtual

Reimplemented in G4VtkQtViewer.

Definition at line 788 of file G4VtkViewer.cc.

789{
790 widget->SetInteractor(_renderWindow->GetInteractor());
791}

Referenced by AddClipperPlaneWidget(), and AddCutterPlaneWidget().

◆ ShowView()

void G4VtkViewer::ShowView ( void )
overridevirtual

Reimplemented from G4VViewer.

Definition at line 303 of file G4VtkViewer.cc.

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}

◆ StartInteractor()

void G4VtkViewer::StartInteractor ( )
inline

Definition at line 272 of file G4VtkViewer.hh.

272 {
273 G4cout << "StartInteractor" << G4endl;
274 _renderWindow->GetInteractor()->Start();}

Member Data Documentation

◆ _renderWindow

◆ bClipper

bool G4VtkViewer::bClipper = false
protected

Definition at line 308 of file G4VtkViewer.hh.

◆ bCutter

bool G4VtkViewer::bCutter = false
protected

Definition at line 307 of file G4VtkViewer.hh.

◆ bHud

bool G4VtkViewer::bHud = false
protected

Definition at line 309 of file G4VtkViewer.hh.

◆ bOrientation

bool G4VtkViewer::bOrientation = false
protected

Definition at line 310 of file G4VtkViewer.hh.

◆ camera

vtkNew<vtkCamera> G4VtkViewer::camera

Definition at line 289 of file G4VtkViewer.hh.

Referenced by Initialise().

◆ cameraDistance

G4double G4VtkViewer::cameraDistance
protected

Definition at line 297 of file G4VtkViewer.hh.

Referenced by SetView().

◆ camOrientWidget

vtkNew<vtkCameraOrientationWidget> G4VtkViewer::camOrientWidget
protected

◆ clipperPlaneRepresentation

vtkNew<vtkImplicitPlaneRepresentation> G4VtkViewer::clipperPlaneRepresentation
protected

Definition at line 302 of file G4VtkViewer.hh.

Referenced by AddClipperPlaneWidget().

◆ clipperPlaneWidget

vtkNew<vtkImplicitPlaneWidget2> G4VtkViewer::clipperPlaneWidget
protected

◆ cutterPlaneRepresentation

vtkNew<vtkImplicitPlaneRepresentation> G4VtkViewer::cutterPlaneRepresentation
protected

Definition at line 299 of file G4VtkViewer.hh.

Referenced by AddCutterPlaneWidget(), and Print().

◆ cutterPlaneWidget

vtkNew<vtkImplicitPlaneWidget2> G4VtkViewer::cutterPlaneWidget
protected

Definition at line 300 of file G4VtkViewer.hh.

Referenced by AddCutterPlaneWidget(), DisableCutterWidget(), and EnableCutterWidget().

◆ firstFinishView

G4bool G4VtkViewer::firstFinishView = true
protected

Definition at line 296 of file G4VtkViewer.hh.

Referenced by FinishView().

◆ firstSetView

G4bool G4VtkViewer::firstSetView = true
protected

Definition at line 295 of file G4VtkViewer.hh.

Referenced by SetView().

◆ geant4Callback

vtkNew<vtkGeant4Callback> G4VtkViewer::geant4Callback

Definition at line 287 of file G4VtkViewer.hh.

Referenced by G4VtkQtViewer::Initialise(), Initialise(), SetView(), and ShowView().

◆ infoCallback

vtkNew<vtkInfoCallback> G4VtkViewer::infoCallback

Definition at line 286 of file G4VtkViewer.hh.

Referenced by AddViewHUD(), and ShowView().

◆ infoTextActor

vtkNew<vtkTextActor> G4VtkViewer::infoTextActor

Definition at line 285 of file G4VtkViewer.hh.

Referenced by AddViewHUD(), DisableHUD(), EnableHUD(), and ShowView().

◆ light

vtkSmartPointer<vtkLight> G4VtkViewer::light

Definition at line 288 of file G4VtkViewer.hh.

◆ renderer

vtkNew<vtkRenderer> G4VtkViewer::renderer

◆ renderWindowInteractor

vtkRenderWindowInteractor* G4VtkViewer::renderWindowInteractor

Definition at line 292 of file G4VtkViewer.hh.

Referenced by Initialise().


The documentation for this class was generated from the following files: