Geant4-11
G4Qt3DViewer.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// John Allison 17th June 2019
27
28#include "G4Qt3DViewer.hh"
29
30#include "G4Qt3DSceneHandler.hh"
31#include "G4Qt3DUtils.hh"
32
33#include "G4Scene.hh"
34#include "G4UImanager.hh"
35#include "G4UIQt.hh"
36#include "G4SystemOfUnits.hh"
37
39(G4Qt3DSceneHandler& sceneHandler, const G4String& name)
40: G4VViewer(sceneHandler, sceneHandler.IncrementViewCount(), name)
41, fQt3DSceneHandler(sceneHandler)
42, fKeyPressed(false)
43, fMousePressed(false)
44, fMousePressedX(0.)
45, fMousePressedY(0.)
46{}
47
49{
50 setObjectName(fName.c_str());
51
52 fVP.SetAutoRefresh(true);
54
55 auto UI = G4UImanager::GetUIpointer();
56 auto uiQt = dynamic_cast<G4UIQt*>(UI->GetG4UIWindow());
57 if (!uiQt) {
58 fViewId = -1; // This flags an error.
59 G4cerr << "G4Qt3DViewer::G4Qt3DViewer requires G4UIQt"
60 << G4endl;
61 return;
62 }
63 fUIWidget = QWidget::createWindowContainer(this);
64 fUIWidget->setMinimumSize(QSize(200, 100));
65 fUIWidget->setMaximumSize(screen()->size());
66// fUIWidget->setFocusPolicy(Qt::NoFocus); //??
67 uiQt->AddTabWidget(fUIWidget,QString(fName));
68
69 setRootEntity(fQt3DSceneHandler.fpQt3DScene);
70}
71
73{}
74
76{
77 // Background colour
78 defaultFrameGraph()->setClearColor(G4Qt3DUtils::ConvertToQColor(fVP.GetBackgroundColour()));
79
80 // Get radius of scene, etc.
81 // Note that this procedure properly takes into account zoom, dolly and pan.
82 const G4Point3D targetPoint
86 if(radius<=0.) radius = 1.;
87 const G4double cameraDistance = fVP.GetCameraDistance (radius);
88 const G4Point3D cameraPosition =
89 targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
90 const GLdouble pnear = fVP.GetNearDistance (cameraDistance, radius);
91 const GLdouble pfar = fVP.GetFarDistance (cameraDistance, pnear, radius);
92 const GLdouble right = fVP.GetFrontHalfHeight (pnear, radius);
93 const GLdouble left = -right;
94 const GLdouble top = fVP.GetFrontHalfHeight (pnear, radius);
95 const GLdouble bottom = -top;
96
97 camera()->setObjectName((fName + " camera").c_str());
98 camera()->setViewCenter(G4Qt3DUtils::ConvertToQVector3D(targetPoint));
99 camera()->setPosition(G4Qt3DUtils::ConvertToQVector3D(cameraPosition));
100 camera()->setUpVector(G4Qt3DUtils::ConvertToQVector3D(fVP.GetUpVector()));
101
102// auto lightEntity = new Qt3DCore::QEntity(fQt3DSceneHandler.fpQt3DScene);
103// auto directionalLight = new Qt3DRender::QDirectionalLight(lightEntity);
106// directionalLight->setWorldDirection(G4Qt3DUtils::ConvertToQVector3D(fVP.GetActualLightpointDirection()));
107// lightEntity->addComponent(directionalLight);
108
109 const auto& size = fUIWidget->size();
110 G4double w = size.width();
111 G4double h = size.height();
112#ifdef G4QT3DDEBUG
113 // Curiously w,h are wrong first time - 640,480 instead of (my Mac) 991,452.
114 G4cout << "W,H: " << w << ',' << h << G4endl;
115#endif
116 const G4double aspectRatio = w/h;
117 if (fVP.GetFieldHalfAngle() == 0.) {
118 camera()->lens()->setOrthographicProjection
119 (left*aspectRatio,right*aspectRatio,bottom,top,pnear,pfar);
120 } else {
121 camera()->lens()->setPerspectiveProjection
122 (2.*fVP.GetFieldHalfAngle()/deg,aspectRatio,pnear,pfar);
123 }
124}
125
127{}
128
130{
131 // First, a view should decide when to re-visit the G4 kernel.
132 // Sometimes it might not be necessary, e.g., if the scene is stored
133 // in a graphical database (e.g., OpenGL's display lists) and only
134 // the viewing angle has changed. But graphics systems without a
135 // graphical database will always need to visit the G4 kernel.
136
137 // The fNeedKernelVisit flag might have been set by the user in
138 // /vis/viewer/rebuild, but if not, make decision and set flag only
139 // if necessary...
141 G4bool kernelVisitWasNeeded = fNeedKernelVisit; // Keep (ProcessView resets).
142 fLastVP = fVP;
143
144 ProcessView (); // Clears store and processes scene only if necessary.
145
146 if (kernelVisitWasNeeded) {
147 // We might need to do something if the kernel was visited.
148 } else {
149 }
150
151 // ...before finally...
152 FinishView (); // Flush streams and/or swap buffers.
153}
154
156{
157 // show() may only be called from master thread
159 show();
160 }
161 // The way Qt seems to work, we don't seem to need a show() anyway, but
162 // we'll leave it in - it seems not to have any effect, good or bad.
163}
164
166{
168 show();
169 }
170}
171
172namespace {
173 QThread* masterQThread = nullptr;
174 QThread* visSubThreadQThread = nullptr;
175}
176
177# include <chrono>
178
180{
181#ifdef G4QT3DDEBUG
182 G4cout << "G4Qt3DViewer::MovingToVisSubThread" << G4endl;
183#endif
184 // Still on master thread but vis thread has been launched
185 // Make note of master QThread
186 masterQThread = QThread::currentThread();
187 // Wait until SwitchToVisSubThread has found vis sub-thread QThread
188 std::this_thread::sleep_for(std::chrono::milliseconds(100));
189 // Move relevant stuff to vis sub-thread QThread
190 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
191 auto p2 = p1->parent();
192 p2->moveToThread(visSubThreadQThread);
193}
194
196{
197#ifdef G4QT3DDEBUG
198 G4cout << "G4Qt3DViewer::SwitchToVisSubThread" << G4endl;
199#endif
200 // On vis sub-thread before any drawing
201 // Make note of vis-subthread QThread for MovingToVisSubThread
202 visSubThreadQThread = QThread::currentThread();
203 // Wait until SwitchToVisSubThread has moved stuff
204 std::this_thread::sleep_for(std::chrono::milliseconds(1000));
205}
206
208{
209#ifdef G4QT3DDEBUG
210 G4cout << "G4Qt3DViewer::MovingToMasterThread" << G4endl;
211#endif
212 // On vis sub-thread just before exit
213 // Move relevant stuff to master QThread.
214 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
215 auto p2 = p1->parent();
216 p2->moveToThread(masterQThread);
217 // Zero - will be different next run
218 visSubThreadQThread = nullptr;
219}
220
222{
223#ifdef G4QT3DDEBUG
224 G4cout << "G4Qt3DViewer::SwitchToMasterThread" << G4endl;
225#endif
226 // On master thread after vis sub-thread has terminated
227 // Nothing to do
228}
229
231
232 // If there's a significant difference with the last view parameters
233 // of either the scene handler or this viewer, trigger a rebuild.
234
236 NeedKernelVisit (); // Sets fNeedKernelVisit.
237 }
238}
239
241{
242 // Typical comparison. Taken from OpenGL.
243 if (
244 (lastVP.GetDrawingStyle () != fVP.GetDrawingStyle ()) ||
246 (lastVP.IsAuxEdgeVisible () != fVP.IsAuxEdgeVisible ()) ||
247 (lastVP.IsCulling () != fVP.IsCulling ()) ||
248 (lastVP.IsCullingInvisible () != fVP.IsCullingInvisible ()) ||
249 (lastVP.IsDensityCulling () != fVP.IsDensityCulling ()) ||
250 (lastVP.IsCullingCovered () != fVP.IsCullingCovered ()) ||
251 (lastVP.GetCBDAlgorithmNumber() !=
253 (lastVP.IsSection () != fVP.IsSection ()) ||
254 (lastVP.IsCutaway () != fVP.IsCutaway ()) ||
255 (lastVP.IsExplode () != fVP.IsExplode ()) ||
256 (lastVP.GetNoOfSides () != fVP.GetNoOfSides ()) ||
259 (lastVP.IsMarkerNotHidden () != fVP.IsMarkerNotHidden ()) ||
260 (lastVP.GetDefaultVisAttributes()->GetColour() !=
265 (lastVP.IsPicking () != fVP.IsPicking ()) ||
266 (lastVP.GetVisAttributesModifiers() !=
268 (lastVP.IsSpecialMeshRendering() !=
270 ) {
271 return true;
272 }
273
274 if (lastVP.IsDensityCulling () &&
275 (lastVP.GetVisibleDensity () != fVP.GetVisibleDensity ()))
276 return true;
277
278 if (lastVP.GetCBDAlgorithmNumber() > 0) {
279 if (lastVP.GetCBDParameters().size() != fVP.GetCBDParameters().size()) return true;
280 else if (lastVP.GetCBDParameters() != fVP.GetCBDParameters()) return true;
281 }
282
283 if (lastVP.IsExplode () &&
284 (lastVP.GetExplodeFactor () != fVP.GetExplodeFactor ()))
285 return true;
286
287 if (lastVP.IsSpecialMeshRendering() &&
289 return true;
290
291 return false;
292}
293
295{
296 fKeyPressed = true;
297 fKey = ev->key();
298}
299
300void G4Qt3DViewer::keyReleaseEvent(QKeyEvent* /*ev*/)
301{
302 fKeyPressed = false;
303}
304
305void G4Qt3DViewer::mouseDoubleClickEvent(QMouseEvent* /*ev*/) {}
306
307void G4Qt3DViewer::mouseMoveEvent(QMouseEvent* ev)
308{
309 // I think we only want these if a mouse button is pressed.
310 // But they come even when not pressed (on my MacBook Pro trackpad).
311 // Documentation says:
312 /* Mouse move events will occur only when a mouse button is pressed down,
313 unless mouse tracking has been enabled with QWidget::setMouseTracking().*/
314 // But this is a window not a widget.
315 // As a workaround we maintain a flag changed by mousePress/ReleaseEvent.
316
317 G4double x = ev->x();
318 G4double y = ev->y();
321 fMousePressedX = x;
322 fMousePressedY = y;
323
324 if (fMousePressed) {
325
326 if (fKeyPressed && fKey == Qt::Key_Shift) { // Translation (pan)
327
329 const G4double scale = 300; // Roughly pixels per window, empirically chosen
330 const G4double dxScene = dx*sceneRadius/scale;
331 const G4double dyScene = dy*sceneRadius/scale;
332 fVP.IncrementPan(-dxScene,dyScene);
333
334 } else { // Rotation
335
336 // Simple ad-hoc algorithms
338 const G4Vector3D& y_prime = x_prime.cross(fVP.GetViewpointDirection());
339 const G4double scale = 200; // Roughly pixels per window, empirically chosen
340 G4Vector3D newViewpointDirection = fVP.GetViewpointDirection();
341 newViewpointDirection += dx*x_prime/scale;
342 newViewpointDirection += dy*y_prime/scale;
343 fVP.SetViewpointDirection(newViewpointDirection.unit());
344
346 G4Vector3D newUpVector = fVP.GetUpVector();
347 newUpVector += dx*x_prime/scale;
348 newUpVector += dy*y_prime/scale;
349 fVP.SetUpVector(newUpVector.unit());
350 }
351 }
352 }
353
354 SetView();
355 DrawView();
356}
357
358void G4Qt3DViewer::mousePressEvent(QMouseEvent* ev)
359{
360 fMousePressed = true;
361 fMousePressedX = ev->x();
362 fMousePressedY = ev->y();
363}
364
365void G4Qt3DViewer::mouseReleaseEvent(QMouseEvent* /*ev*/)
366{
367 fMousePressed = false;
368}
369
370void G4Qt3DViewer::wheelEvent(QWheelEvent* ev)
371{
372 // Take note of up-down motion only
373 const G4double angleY = ev->angleDelta().y();
374
375 if (fVP.GetFieldHalfAngle() == 0.) { // Orthographic projection
376 const G4double scale = 500; // Empirically chosen
377 fVP.MultiplyZoomFactor(1.+angleY/scale);
378 } else { // Perspective projection
379 const G4double scale = fVP.GetFieldHalfAngle()/(10.*deg); // Empirical
380 fVP.SetDolly(fVP.GetDolly()+angleY/scale);
381 }
382
383 SetView();
384 DrawView();
385}
static constexpr double deg
Definition: G4SIunits.hh:132
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Qt3DCore::QEntity * fpQt3DScene
void mousePressEvent(QMouseEvent *)
void KernelVisitDecision()
G4Qt3DSceneHandler & fQt3DSceneHandler
Definition: G4Qt3DViewer.hh:69
G4double fMousePressedY
Definition: G4Qt3DViewer.hh:76
void Initialise()
Definition: G4Qt3DViewer.cc:48
G4bool fMousePressed
Definition: G4Qt3DViewer.hh:75
void ClearView()
void SwitchToMasterThread()
void SwitchToVisSubThread()
void MovingToMasterThread()
G4double fMousePressedX
Definition: G4Qt3DViewer.hh:76
void mouseReleaseEvent(QMouseEvent *)
virtual ~G4Qt3DViewer()
Definition: G4Qt3DViewer.cc:72
G4bool fKeyPressed
Definition: G4Qt3DViewer.hh:73
void wheelEvent(QWheelEvent *)
void keyPressEvent(QKeyEvent *)
QWidget * fUIWidget
Definition: G4Qt3DViewer.hh:71
G4Qt3DViewer(G4Qt3DSceneHandler &, const G4String &name)
Definition: G4Qt3DViewer.cc:39
void mouseDoubleClickEvent(QMouseEvent *)
void MovingToVisSubThread()
void mouseMoveEvent(QMouseEvent *)
void FinishView()
G4bool CompareForKernelVisit(G4ViewParameters &)
G4ViewParameters fLastVP
Definition: G4Qt3DViewer.hh:68
void SetView()
Definition: G4Qt3DViewer.cc:75
void keyReleaseEvent(QKeyEvent *)
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4Scene * GetScene() const
G4bool fNeedKernelVisit
Definition: G4VViewer.hh:224
void ProcessView()
Definition: G4VViewer.cc:105
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:215
G4String fName
Definition: G4VViewer.hh:217
void NeedKernelVisit()
Definition: G4VViewer.cc:78
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:220
G4int fViewId
Definition: G4VViewer.hh:216
G4ViewParameters fVP
Definition: G4VViewer.hh:219
void SetViewpointDirection(const G4Vector3D &viewpointDirection)
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
void SetAutoRefresh(G4bool)
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
G4double GetCameraDistance(G4double radius) const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsMarkerNotHidden() const
G4double GetGlobalLineWidthScale() const
G4bool IsCutaway() const
const G4Colour & GetBackgroundColour() const
const G4Vector3D & GetViewpointDirection() const
G4bool IsSection() const
const G4Point3D & GetCurrentTargetPoint() const
G4bool IsPicking() const
G4double GetFarDistance(G4double cameraDistance, G4double nearDistance, G4double radius) const
G4double GetFieldHalfAngle() const
G4bool IsCulling() const
G4double GetFrontHalfHeight(G4double nearDistance, G4double radius) const
const G4VisAttributes * GetDefaultTextVisAttributes() const
void SetDolly(G4double dolly)
G4bool IsExplode() const
void IncrementPan(G4double right, G4double up)
const G4Vector3D & GetUpVector() const
const std::vector< G4double > & GetCBDParameters() const
G4int GetCBDAlgorithmNumber() const
const std::vector< G4ModelingParameters::PVNameCopyNo > & GetSpecialMeshVolumes() const
G4double GetGlobalMarkerScale() const
G4bool IsCullingInvisible() const
const G4VisAttributes * GetDefaultVisAttributes() const
void SetUpVector(const G4Vector3D &upVector)
RotationStyle GetRotationStyle() const
G4bool IsDensityCulling() const
void MultiplyZoomFactor(G4double zoomFactorMultiplier)
G4double GetVisibleDensity() const
G4bool IsCullingCovered() const
G4double GetNearDistance(G4double cameraDistance, G4double radius) const
DrawingStyle GetDrawingStyle() const
G4bool IsAuxEdgeVisible() const
G4double GetDolly() const
const G4Colour & GetColour() const
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > unit() const
const char * name(G4int ptype)
QColor ConvertToQColor(const G4Colour &c)
Definition: G4Qt3DUtils.cc:46
QVector3D ConvertToQVector3D(const G4ThreeVector &v)
Definition: G4Qt3DUtils.cc:52
G4bool IsMasterThread()
Definition: G4Threading.cc:124