Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkViewer.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//
29// John Allison 5th April 2001
30// A template for a simplest possible graphics driver.
31//?? Lines or sections marked like this require specialisation for your driver.
32
33#ifndef G4VTKVIEWER_HH
34#define G4VTKVIEWER_HH
35
36// #define G4VTKDEBUG // Comment this out to suppress debug code.
37
38#include "G4VViewer.hh"
39
40#pragma GCC diagnostic push
41#pragma GCC diagnostic ignored "-Wextra-semi"
42#include "vtkObject.h"
43#include "vtkAutoInit.h"
44#include "vtkCamera.h"
45#include "vtkLight.h"
46#include "vtkRenderer.h"
47#include "vtkRenderWindow.h"
48#include "vtkRenderWindowInteractor.h"
49#include "vtkInteractorStyleTrackballCamera.h"
50#include "vtkInteractorStyleTerrain.h"
51#include "vtkTextActor.h"
52#pragma GCC diagnostic pop
53
54VTK_MODULE_INIT(vtkRenderingOpenGL2)
55VTK_MODULE_INIT(vtkInteractionStyle);
56VTK_MODULE_INIT(vtkRenderingFreeType)
57
58class vtkGeant4Callback : public vtkCommand
59{
60public:
61 static vtkGeant4Callback* New() {return new vtkGeant4Callback;}
62
63 vtkGeant4Callback() { fVP = nullptr; }
65 void SetVtkInitialValues(G4double parallelScaleIn, G4double cameraDistanceIn)
66 {
67 parallelScale = parallelScaleIn;
68 cameraDistance = cameraDistanceIn;
69 }
70 virtual void Execute(vtkObject *caller, unsigned long, void*)
71 {
72 vtkRenderer *ren = static_cast<vtkRenderer *>(caller);
73 vtkCamera *cam = ren->GetActiveCamera();
74 //G4cout << cam->GetFocalPoint()[0] << " " << cam->GetFocalPoint()[1] << " " << cam->GetFocalPoint()[2] << G4endl;
75 //
76 // G4cout << cam->GetPosition()[0] << " " << cam->GetPosition()[1] << " " << cam->GetPosition()[2] << G4endl;
77
78 auto cp = cam->GetPosition();
79 auto fp = cam->GetFocalPoint();
80 auto ud = cam->GetViewUp();
81
82 fVP->SetCurrentTargetPoint(G4Point3D(fp[0],fp[1],fp[2]));
83 fVP->SetViewpointDirection((G4Point3D(cp[0],cp[1],cp[2]) -
84 G4Point3D(fp[0],fp[1],fp[2])).unit());
85 fVP->SetUpVector(G4Vector3D(ud[0],ud[1], ud[2]));
86
87 if(cam->GetParallelProjection()) {
88 fVP->SetZoomFactor(parallelScale/cam->GetParallelScale());
89 }
90 else {
91 auto cd = std::sqrt(std::pow(cp[0]-fp[0],2) +
92 std::pow(cp[1]-cp[1],2) +
93 std::pow(cp[2]-cp[2],2));
94 fVP->SetZoomFactor(cameraDistance/cd);
95 }
96 }
97
98protected:
102};
103
104class vtkInfoCallback : public vtkCommand
105{
106public:
107 static vtkInfoCallback *New() { return new vtkInfoCallback; }
108
110 t1 = std::chrono::steady_clock::now();
111 t2 = std::chrono::steady_clock::now();
112 }
113 void SetTextActor(vtkTextActor *txt) { this->TextActor = txt; }
114
115 virtual void Execute(vtkObject *caller, unsigned long, void*)
116 {
117 vtkRenderer *ren = static_cast<vtkRenderer *>(caller);
118 int nActors = ren->GetActors()->GetNumberOfItems();
119 vtkCamera *cam = ren->GetActiveCamera();
120 if(!cam) return;
121
122 double *pos = cam->GetPosition();
123 double *foc = cam->GetFocalPoint();
124 double viewAngle = cam->GetViewAngle();
125 double distance = cam->GetDistance();
126 double parallelScale = cam->GetParallelScale();
127
128 if(!pos) return;
129
130 // Get current time
131 t2 = std::chrono::steady_clock::now();
132
133 // Frame rate calculation
134 std::chrono::duration<double> tdiff = t2-t1;
135 t1 = t2;
136 float fps = 1.0/tdiff.count();
137
138 // String for display
139 snprintf(this->TextBuff,sizeof this->TextBuff,
140 "camera position : %.1f %.1f %.1f \n"
141 "camera focal point : %.1f %.1f %.1f \n"
142 "view angle : %.1f\n"
143 "distance : %.1f\n"
144 "parallel scale : %.1f\n"
145 "number actors : %i\n"
146 "fps : %.1f",pos[0], pos[1], pos[2], foc[0], foc[1], foc[2], viewAngle, distance, parallelScale, nActors, fps);
147 if(this->TextActor) {
148 this->TextActor->SetInput(this->TextBuff);
149 }
150 }
151protected:
152 vtkTextActor *TextActor;
153 char TextBuff[256];
154 std::chrono::time_point<std::chrono::steady_clock> t1;
155 std::chrono::time_point<std::chrono::steady_clock> t2;
156};
157
158class G4VtkViewer: public G4VViewer {
159 public:
160 G4VtkViewer(G4VSceneHandler&, const G4String& name);
161 void Initialise();
162 virtual ~G4VtkViewer();
163
164 void SetView();
165 void ClearView();
166 void DrawView();
167 void ShowView();
168 void FinishView();
169
170 void DrawViewHUD();
171 void DrawShadows();
172
177
178 void ExportView(){};
180
181 vtkNew<vtkTextActor> infoTextActor;
182 vtkNew<vtkInfoCallback> infoCallback;
183 vtkNew<vtkGeant4Callback> geant4Callback;
184 vtkSmartPointer<vtkLight> light;
185 vtkNew<vtkCamera> camera;
186 vtkNew<vtkRenderer> renderer;
187 vtkRenderWindow* _renderWindow;
188 vtkRenderWindowInteractor* renderWindowInteractor;
189
190 private:
191 G4bool firstSetView = true;
192};
193
194#endif
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:34
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
VTK_MODULE_INIT(vtkInteractionStyle)
vtkNew< vtkCamera > camera
Definition: G4VtkViewer.hh:185
vtkSmartPointer< vtkLight > light
Definition: G4VtkViewer.hh:184
void ExportVRMLScene(G4String)
Definition: G4VtkViewer.cc:387
void DrawView()
Definition: G4VtkViewer.cc:243
void SetView()
Definition: G4VtkViewer.cc:125
vtkNew< vtkGeant4Callback > geant4Callback
Definition: G4VtkViewer.hh:183
void ExportScreenShot(G4String, G4String)
Definition: G4VtkViewer.cc:331
void DrawShadows()
Definition: G4VtkViewer.cc:276
vtkNew< vtkTextActor > infoTextActor
Definition: G4VtkViewer.hh:181
vtkRenderWindowInteractor * renderWindowInteractor
Definition: G4VtkViewer.hh:188
void ShowView()
Definition: G4VtkViewer.cc:296
vtkNew< vtkInfoCallback > infoCallback
Definition: G4VtkViewer.hh:182
void SetGeant4View()
Definition: G4VtkViewer.hh:179
vtkNew< vtkRenderer > renderer
Definition: G4VtkViewer.hh:186
void Initialise()
Definition: G4VtkViewer.cc:65
void ExportView()
Definition: G4VtkViewer.hh:178
void ExportVTPScene(G4String)
Definition: G4VtkViewer.cc:397
void ExportOBJScene(G4String)
Definition: G4VtkViewer.cc:377
virtual ~G4VtkViewer()
Definition: G4VtkViewer.cc:123
void FinishView()
Definition: G4VtkViewer.cc:321
vtkRenderWindow * _renderWindow
Definition: G4VtkViewer.hh:187
void ClearView()
Definition: G4VtkViewer.cc:219
void DrawViewHUD()
Definition: G4VtkViewer.cc:263
G4ViewParameters * fVP
Definition: G4VtkViewer.hh:99
void SetGeant4ViewParameters(G4ViewParameters *VP)
Definition: G4VtkViewer.hh:64
G4double parallelScale
Definition: G4VtkViewer.hh:100
static vtkGeant4Callback * New()
Definition: G4VtkViewer.hh:61
virtual void Execute(vtkObject *caller, unsigned long, void *)
Definition: G4VtkViewer.hh:70
void SetVtkInitialValues(G4double parallelScaleIn, G4double cameraDistanceIn)
Definition: G4VtkViewer.hh:65
G4double cameraDistance
Definition: G4VtkViewer.hh:101
char TextBuff[256]
Definition: G4VtkViewer.hh:153
void SetTextActor(vtkTextActor *txt)
Definition: G4VtkViewer.hh:113
vtkTextActor * TextActor
Definition: G4VtkViewer.hh:152
virtual void Execute(vtkObject *caller, unsigned long, void *)
Definition: G4VtkViewer.hh:115
std::chrono::time_point< std::chrono::steady_clock > t2
Definition: G4VtkViewer.hh:155
static vtkInfoCallback * New()
Definition: G4VtkViewer.hh:107
std::chrono::time_point< std::chrono::steady_clock > t1
Definition: G4VtkViewer.hh:154