Geant4 11.2.2
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#include "G4VViewer.hh"
38#include "G4VtkStore.hh"
39#include "G4VtkUtility.hh"
40
41#ifndef WIN32
42#pragma GCC diagnostic push
43#pragma GCC diagnostic ignored "-Wextra-semi"
44#endif
45
46#include <vtkAutoInit.h>
47#include <vtkCamera.h>
48#include <vtkCameraOrientationWidget.h>
49#include <vtkImplicitPlaneRepresentation.h>
50#include <vtkImplicitPlaneWidget2.h>
51#include <vtkInteractorStyleTerrain.h>
52#include <vtkInteractorStyleTrackballCamera.h>
53#include <vtkLight.h>
54#include <vtkNew.h>
55#include <vtkObject.h>
56#include <vtkPlane.h>
57#include <vtkRenderWindow.h>
58#include <vtkRenderWindowInteractor.h>
59#include <vtkRenderer.h>
60#include <vtkTextActor.h>
61
62#ifndef WIN32
63#pragma GCC diagnostic pop
64#endif
65
66VTK_MODULE_INIT(vtkRenderingOpenGL2)
67VTK_MODULE_INIT(vtkInteractionStyle);
68VTK_MODULE_INIT(vtkRenderingFreeType)
69
70class vtkGeant4Callback : public vtkCommand
71{
72 public:
73 static vtkGeant4Callback* New() { return new vtkGeant4Callback; }
74
75 vtkGeant4Callback() { fVP = nullptr; }
77 void SetVtkInitialValues(G4double parallelScaleIn, G4double cameraDistanceIn)
78 {
79 parallelScale = parallelScaleIn;
80 cameraDistance = cameraDistanceIn;
81 }
82 void Execute(vtkObject* caller, unsigned long, void*) override
83 {
84 auto ren = static_cast<vtkRenderer*>(caller);
85 vtkCamera* cam = ren->GetActiveCamera();
86
87 auto cp = cam->GetPosition();
88 auto fp = cam->GetFocalPoint();
89 auto ud = cam->GetViewUp();
90 auto ps = cam->GetParallelScale();
91 auto cd = std::sqrt(std::pow(cp[0] - fp[0], 2) + std::pow(cp[1] - fp[1], 2)
92 + std::pow(cp[2] - fp[2], 2));
93
94 fVP->SetCurrentTargetPoint(G4Point3D(fp[0], fp[1], fp[2]));
95 fVP->SetViewpointDirection((G4Point3D(cp[0], cp[1], cp[2]) - G4Point3D(fp[0], fp[1], fp[2])).unit());
96 fVP->SetUpVector(G4Vector3D(ud[0], ud[1], ud[2]));
97
98 if (cam->GetParallelProjection() != 0) {
99 fVP->SetZoomFactor(parallelScale / ps);
100 }
101 else {
102 fVP->SetZoomFactor(cameraDistance / cd);
103 }
104 }
105
106 protected:
110};
111
112class vtkInfoCallback : public vtkCommand
113{
114 public:
115 static vtkInfoCallback* New() { return new vtkInfoCallback; }
116
118 {
119 t1 = std::chrono::steady_clock::now();
120 t2 = std::chrono::steady_clock::now();
121 }
122 void SetTextActor(vtkTextActor* txt) { this->TextActor = txt; }
123
124 void Execute(vtkObject* caller, unsigned long, void*) override
125 {
126 auto ren = static_cast<vtkRenderer*>(caller);
127 int nActors = ren->GetActors()->GetNumberOfItems();
128 vtkCamera* cam = ren->GetActiveCamera();
129 if (cam == nullptr) return;
130
131 double* pos = cam->GetPosition();
132 double* foc = cam->GetFocalPoint();
133 double viewAngle = cam->GetViewAngle();
134 double distance = cam->GetDistance();
135 double near;
136 double far;
137 cam->GetClippingRange(near,far);
138 double parallelScale = cam->GetParallelScale();
139
140
141 if (pos == nullptr) return;
142
143 // Get current time
144 t2 = std::chrono::steady_clock::now();
145
146 // Frame rate calculation
147 std::chrono::duration<double> tdiff = t2 - t1;
148 t1 = t2;
149 float fps = 1.0 / tdiff.count();
150
151 // String for display
152 snprintf(this->TextBuff, sizeof this->TextBuff,
153 "camera position : %.1f %.1f %.1f \n"
154 "camera focal point : %.1f %.1f %.1f \n"
155 "view angle : %.1f\n"
156 "distance : %.1f\n"
157 "clip near/far : %.1f %.1f\n"
158 "parallel scale : %.1f\n"
159 "number actors : %i\n"
160 "fps : %.1f",
161 pos[0], pos[1], pos[2], foc[0], foc[1], foc[2], viewAngle, distance, near, far, parallelScale,
162 nActors, fps);
163 if (this->TextActor != nullptr) {
164 this->TextActor->SetInput(this->TextBuff);
165 }
166 }
167
168 protected:
169 vtkTextActor* TextActor;
170 char TextBuff[256];
171 std::chrono::time_point<std::chrono::steady_clock> t1;
172 std::chrono::time_point<std::chrono::steady_clock> t2;
173};
174
175class vtkIPWCallback : public vtkCommand
176{
177 public:
178 static vtkIPWCallback* New() { return new vtkIPWCallback; }
179
180 vtkIPWCallback() = default;
181
182 void SetStore(G4VtkStore* storeIn) { store = storeIn; }
183
185 {
186 pipelineToUpdateName = nameIn;
187 pipelineToUpdateType = typeIn;
188 }
189
190 void Execute(vtkObject* caller, unsigned long, void*) override
191 {
192 if (this->plane == nullptr) {
193 this->plane = vtkPlane::New();
194 }
195
196 if (pipelineToUpdateName.empty() || pipelineToUpdateType.empty()) {
197 return;
198 }
199
200 auto planeWidget = reinterpret_cast<vtkImplicitPlaneWidget2*>(caller);
201 auto rep =
202 reinterpret_cast<vtkImplicitPlaneRepresentation*>(planeWidget->GetRepresentation());
203 rep->GetPlane(this->plane);
206 }
207
208 vtkPlane* plane{nullptr};
209 G4VtkStore* store{nullptr};
212};
213
214class G4VtkViewer : public G4VViewer
215{
216 public:
217 G4VtkViewer(G4VSceneHandler&, const G4String& name);
218 void Initialise() override;
219 ~G4VtkViewer() override;
220
221 void SetView() override;
222 void ClearView() override;
223 void DrawView() override;
224 void ShowView() override;
225 void FinishView() override;
226
233 void ExportVTPCutter(G4String fileName);
234 void ExportFormatStore(G4String fileName, G4String store);
235
236 void DrawShadows();
237
238 void EnableShadows();
239 void DisableShadows();
240
241 void AddViewHUD();
242 void EnableHUD();
243 void DisableHUD();
244
245 virtual void AddClipperPlaneWidget(const G4Plane3D& plane);
246 void EnableClipper(const G4Plane3D& plane, G4bool widget);
247 void DisableClipper();
248 virtual void EnableClipperWidget();
249 virtual void DisableClipperWidget();
250
251 virtual void AddCutterPlaneWidget(const G4Plane3D& plane);
252 void EnableCutter(const G4Plane3D& plane, G4bool bWidget);
253 void DisableCutter(G4String name);
254 virtual void EnableCutterWidget();
255 virtual void DisableCutterWidget();
256
257 virtual void AddCameraOrientationWidget();
258 virtual void EnableCameraOrientationWidget();
259 virtual void DisableCameraOrientationWidget();
260
261 void AddImageOverlay(const G4String& fileName, const G4double alpha,
262 const G4double imageBottomLeft[2], const G4double worldBottomLeft[2],
263 const G4double imageTopRight[2], const G4double worldTopRight[2],
264 const G4double rot[3], const G4double trans[3]);
265 void AddGeometryOverlay(const G4String& fileName, const G4double colour[3], const G4double alpha,
266 const G4String& representation,
267 const G4double scale[3], const G4double rotation[3],
268 const G4double translation[3]);
269
270 void Render() {_renderWindow->Render();}
272 G4cout << "StartInteractor" << G4endl;
273 _renderWindow->GetInteractor()->Start();}
274
275 void Print();
276
277 void SetPolyhedronPipeline(const G4String& t);
278
279 virtual void SetWidgetInteractor(vtkAbstractWidget* widget);
280
281 void ExportView(){};
283
284 vtkNew<vtkTextActor> infoTextActor;
285 vtkNew<vtkInfoCallback> infoCallback;
286 vtkNew<vtkGeant4Callback> geant4Callback;
288 vtkNew<vtkCamera> camera;
289 vtkNew<vtkRenderer> renderer;
290 vtkRenderWindow* _renderWindow;
291 vtkRenderWindowInteractor* renderWindowInteractor;
292
293 protected:
297
298 vtkNew<vtkImplicitPlaneRepresentation> cutterPlaneRepresentation;
299 vtkNew<vtkImplicitPlaneWidget2> cutterPlaneWidget;
300
301 vtkNew<vtkImplicitPlaneRepresentation> clipperPlaneRepresentation;
302 vtkNew<vtkImplicitPlaneWidget2> clipperPlaneWidget;
303
304 vtkNew<vtkCameraOrientationWidget> camOrientWidget;
305
306 bool bCutter = false;
307 bool bClipper = false;
308 bool bHud = false;
309 bool bOrientation = false;
310};
311
312#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
G4Plane3D VtkPlaneToG4Plane3D(vtkPlane *vtkPlane)
VTK_MODULE_INIT(vtkInteractionStyle)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void UpdatePlanePipelines(G4String name, G4String type, const G4Plane3D)
vtkNew< vtkCamera > camera
void DisableClipper()
vtkSmartPointer< vtkLight > light
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 SetGeant4View()
void EnableHUD()
vtkNew< vtkRenderer > renderer
G4double cameraDistance
void AddViewHUD()
void ExportView()
void EnableCutter(const G4Plane3D &plane, G4bool bWidget)
void StartInteractor()
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 * fVP
void SetGeant4ViewParameters(G4ViewParameters *VP)
static vtkGeant4Callback * New()
void Execute(vtkObject *caller, unsigned long, void *) override
void SetVtkInitialValues(G4double parallelScaleIn, G4double cameraDistanceIn)
G4double cameraDistance
G4String pipelineToUpdateType
vtkIPWCallback()=default
void Execute(vtkObject *caller, unsigned long, void *) override
void SetUpdatePipelineName(G4String nameIn, G4String typeIn)
vtkPlane * plane
G4String pipelineToUpdateName
G4VtkStore * store
void SetStore(G4VtkStore *storeIn)
static vtkIPWCallback * New()
char TextBuff[256]
void SetTextActor(vtkTextActor *txt)
vtkTextActor * TextActor
std::chrono::time_point< std::chrono::steady_clock > t2
static vtkInfoCallback * New()
std::chrono::time_point< std::chrono::steady_clock > t1
void Execute(vtkObject *caller, unsigned long, void *) override