Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
38#define G4warn G4cout
39
41(G4Qt3DSceneHandler& sceneHandler, const G4String& name)
42: G4VViewer(sceneHandler, sceneHandler.IncrementViewCount(), name)
43, fQt3DSceneHandler(sceneHandler)
44, fKeyPressed(false)
45, fMousePressed(false)
48{}
49
51{
52 setObjectName(fName.c_str());
53
54 fVP.SetAutoRefresh(true);
55 fDefaultVP.SetAutoRefresh(true);
56
57 auto UI = G4UImanager::GetUIpointer();
58 auto uiQt = dynamic_cast<G4UIQt*>(UI->GetG4UIWindow());
59 if (!uiQt) {
60 fViewId = -1; // This flags an error.
61 G4warn << "G4Qt3DViewer::G4Qt3DViewer requires G4UIQt"
62 << G4endl;
63 return;
64 }
65 fUIWidget = QWidget::createWindowContainer(this);
66 uiQt->AddTabWidget(fUIWidget,QString(fName));
67
68 setRootEntity(fQt3DSceneHandler.fpQt3DScene);
69}
70
72{
73 setRootEntity(nullptr);
74}
75
76void G4Qt3DViewer::resizeEvent(QResizeEvent*) {
77 SetView();
78}
79
81{
82 // Background colour
83 defaultFrameGraph()->setClearColor(G4Qt3DUtils::ConvertToQColor(fVP.GetBackgroundColour()));
84
85 // Get radius of scene, etc.
86 // Note that this procedure properly takes into account zoom, dolly and pan.
87 const G4Point3D targetPoint
88 = fSceneHandler.GetScene()->GetStandardTargetPoint()
89 + fVP.GetCurrentTargetPoint ();
90 G4double radius = fSceneHandler.GetScene()->GetExtent().GetExtentRadius();
91 if(radius<=0.) radius = 1.;
92 const G4double cameraDistance = fVP.GetCameraDistance (radius);
93 const G4Point3D cameraPosition =
94 targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
95 const GLdouble pnear = fVP.GetNearDistance (cameraDistance, radius);
96 const GLdouble pfar = fVP.GetFarDistance (cameraDistance, pnear, radius);
97 const GLdouble right = fVP.GetFrontHalfHeight (pnear, radius);
98 const GLdouble left = -right;
99 const GLdouble top = fVP.GetFrontHalfHeight (pnear, radius);
100 const GLdouble bottom = -top;
101
102 camera()->setObjectName((fName + " camera").c_str());
103 camera()->setViewCenter(G4Qt3DUtils::ConvertToQVector3D(targetPoint));
104 camera()->setPosition(G4Qt3DUtils::ConvertToQVector3D(cameraPosition));
105 camera()->setUpVector(G4Qt3DUtils::ConvertToQVector3D(fVP.GetUpVector()));
106
107// auto lightEntity = new Qt3DCore::QEntity(fQt3DSceneHandler.fpQt3DScene);
108// auto directionalLight = new Qt3DRender::QDirectionalLight(lightEntity);
109//// directionalLight->setColor("white");
110//// directionalLight->setIntensity(1.);
111// directionalLight->setWorldDirection(G4Qt3DUtils::ConvertToQVector3D(fVP.GetActualLightpointDirection()));
112// lightEntity->addComponent(directionalLight);
113
114 const auto& size = fUIWidget->size();
115 G4double w = size.width();
116 G4double h = size.height();
117#ifdef G4QT3DDEBUG
118 // Curiously w,h are wrong first time - 640,480 instead of (my Mac) 991,452.
119 G4cout << "W,H: " << w << ',' << h << G4endl;
120#endif
121 const G4double aspectRatio = w/h;
122 if (fVP.GetFieldHalfAngle() == 0.) {
123 camera()->lens()->setOrthographicProjection
124 (left*aspectRatio,right*aspectRatio,bottom,top,pnear,pfar);
125 } else {
126 camera()->lens()->setPerspectiveProjection
127 (2.*fVP.GetFieldHalfAngle()/deg,aspectRatio,pnear,pfar);
128 }
129}
130
133
135{
136 // First, a view should decide when to re-visit the G4 kernel.
137 // Sometimes it might not be necessary, e.g., if the scene is stored
138 // in a graphical database (e.g., OpenGL's display lists) and only
139 // the viewing angle has changed. But graphics systems without a
140 // graphical database will always need to visit the G4 kernel.
141
142 // The fNeedKernelVisit flag might have been set by the user in
143 // /vis/viewer/rebuild, but if not, make decision and set flag only
144 // if necessary...
146 G4bool kernelVisitWasNeeded = fNeedKernelVisit; // Keep (ProcessView resets).
147 fLastVP = fVP;
148
149 ProcessView (); // Clears store and processes scene only if necessary.
150
151 if (kernelVisitWasNeeded) {
152 // We might need to do something if the kernel was visited.
153 } else {
154 }
155
156 // ...before finally...
157 FinishView (); // Flush streams and/or swap buffers.
158}
159
161{
162#if QT_VERSION < 0x060000
163 // show() may only be called from master thread
165 show();
166 }
167 // The way Qt seems to work, we don't seem to need a show() anyway, but
168 // we'll leave it in - it seems not to have any effect, good or bad.
169#endif
170}
171
173{
174#if QT_VERSION < 0x060000
176 show();
177 }
178#endif
179}
180
181// Note: the order of calling of MovingToVisSubThread and SwitchToVisSubThread
182// is undefined. The order of calling is
183// DoneWithMasterThread
184// MovingToVisSubThread ) or ( SwitchToVisSubThread
185// SwitchToVisSubThread ) ( MovingToVisSubThread
186// DoneWithVisSubThread
187// MovingToMasterThread
188// SwitchToMasterThread
189// So regarding the move/switch to the vis sub-thread, we have to employ mutex locks and conditions.
190// If the viewer wishes to accept drawing from the vis sub-thread, Qt3D has to moveToThread.
191// But at this point we are still on the master thread and the value of the sub-thread's QThread is
192// not known. So it has to wait - a conditional wait, conditional on establishment of the sub-thread
193// and the provision of a pointer to the QThread version of the vis sub-thread. In turn, the
194// sub-thread has to wait until the QObjects have been moved to the sub-thread.
195
196namespace {
197 QThread* masterQThread = nullptr;
198 QThread* visSubThreadQThread = nullptr;
199
200 G4Mutex visSubThreadMutex = G4MUTEX_INITIALIZER;
201 G4Condition waitForVisSubThreadInitialized = G4CONDITION_INITIALIZER;
202 G4bool visSubThreadEstablished = false;
203 G4bool qObjectsSwitched = false;
204}
205
207// Still on master thread but vis thread has been launched
208{
209 // Make note of master QThread
210 masterQThread = QThread::currentThread();
211
212 // Wait until SwitchToVisSubThread has found vis sub-thread QThread
213 {
214 G4AutoLock lock(&visSubThreadMutex);
215 G4CONDITIONWAITLAMBDA(&waitForVisSubThreadInitialized, &lock, []{return visSubThreadEstablished;})
216 }
217
218 // Move relevant stuff to vis sub-thread QThread
219 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
220 if(p1) {
221 auto p2 = p1->parent();
222 if(p2) {
223 p2->moveToThread(visSubThreadQThread);
224 } else {
225 p1->moveToThread(visSubThreadQThread);
226 }
227 }
228
229 // Inform sub-thread
230 G4AutoLock lock(&visSubThreadMutex);
231 qObjectsSwitched = true;
232 lock.unlock();
233 G4CONDITIONBROADCAST(&waitForVisSubThreadInitialized);
234}
235
237// On vis sub-thread before any drawing
238{
239 // Make note of vis-subthread QThread for MovingToVisSubThread
240 visSubThreadQThread = QThread::currentThread();
241
242 // Let MovingToVisSubThread know we have the QThread
243 {
244 G4AutoLock lock(&visSubThreadMutex);
245 visSubThreadEstablished = true;
246 G4CONDITIONBROADCAST(&waitForVisSubThreadInitialized);
247 }
248
249 // Wait until MovingToVisSubThread has moved stuff
250 {
251 G4AutoLock lock(&visSubThreadMutex);
252 G4CONDITIONWAITLAMBDA(&waitForVisSubThreadInitialized, &lock, []{return qObjectsSwitched;})
253 }
254}
255
257// On vis sub-thread just before exit
258{
259 // Move relevant stuff to master QThread.
260 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
261 if(p1) {
262 auto p2 = p1->parent();
263 if(p2) {
264 p2->moveToThread(masterQThread);
265 } else {
266 p1->moveToThread(masterQThread);
267 }
268 }
269
270 // Reset
271 visSubThreadQThread = nullptr;
272 qObjectsSwitched = false;
273}
274
276// On master thread after vis sub-thread has terminated
277{
278 visSubThreadEstablished = false;
279}
280
282
283 // If there's a significant difference with the last view parameters
284 // of either the scene handler or this viewer, trigger a rebuild.
285
287 NeedKernelVisit (); // Sets fNeedKernelVisit.
288 }
289}
290
292{
293 // Typical comparison. Taken from OpenInventor.
294 if (
295 (vp.GetDrawingStyle () != fVP.GetDrawingStyle ()) ||
296 (vp.GetNumberOfCloudPoints() != fVP.GetNumberOfCloudPoints()) ||
297 (vp.IsAuxEdgeVisible () != fVP.IsAuxEdgeVisible ()) ||
298 (vp.IsCulling () != fVP.IsCulling ()) ||
299 (vp.IsCullingInvisible () != fVP.IsCullingInvisible ()) ||
300 (vp.IsDensityCulling () != fVP.IsDensityCulling ()) ||
301 (vp.IsCullingCovered () != fVP.IsCullingCovered ()) ||
302 (vp.GetCBDAlgorithmNumber() !=
303 fVP.GetCBDAlgorithmNumber()) ||
304 (vp.IsSection () != fVP.IsSection ()) ||
305 (vp.IsCutaway () != fVP.IsCutaway ()) ||
306 // This assumes use of generic clipping (sectioning, slicing,
307 // DCUT, cutaway). If a decision is made to implement locally,
308 // this will need changing. See G4OpenGLViewer::SetView,
309 // G4OpenGLStoredViewer.cc::CompareForKernelVisit and
310 // G4OpenGLStoredSceneHander::CreateSection/CutawayPolyhedron.
311 (vp.IsExplode () != fVP.IsExplode ()) ||
312 (vp.GetNoOfSides () != fVP.GetNoOfSides ()) ||
313 (vp.GetGlobalMarkerScale() != fVP.GetGlobalMarkerScale()) ||
314 (vp.GetGlobalLineWidthScale() != fVP.GetGlobalLineWidthScale()) ||
315 (vp.IsMarkerNotHidden () != fVP.IsMarkerNotHidden ()) ||
317 fVP.GetDefaultVisAttributes()->GetColour()) ||
319 fVP.GetDefaultTextVisAttributes()->GetColour()) ||
320 (vp.GetBackgroundColour ()!= fVP.GetBackgroundColour ())||
321 (vp.IsPicking () != fVP.IsPicking ()) ||
322 // Scaling for Open Inventor is done by the scene handler so it
323 // needs a kernel visit. (In this respect, it differs from the
324 // OpenGL drivers, where it's done in SetView.)
325 (vp.GetScaleFactor () != fVP.GetScaleFactor ()) ||
327 fVP.GetVisAttributesModifiers()) ||
329 fVP.IsSpecialMeshRendering()) ||
331 fVP.GetSpecialMeshRenderingOption())
332 )
333 return true;
334
335 if (vp.IsDensityCulling () &&
336 (vp.GetVisibleDensity () != fVP.GetVisibleDensity ()))
337 return true;
338
339 if (vp.GetCBDAlgorithmNumber() > 0) {
340 if (vp.GetCBDParameters().size() != fVP.GetCBDParameters().size()) return true;
341 else if (vp.GetCBDParameters() != fVP.GetCBDParameters()) return true;
342 }
343
344 if (vp.IsSection () &&
345 (vp.GetSectionPlane () != fVP.GetSectionPlane ()))
346 return true;
347
348 if (vp.IsCutaway ()) {
349 if (vp.GetCutawayMode() != fVP.GetCutawayMode()) return true;
350 if (vp.GetCutawayPlanes ().size () !=
351 fVP.GetCutawayPlanes ().size ()) return true;
352 for (size_t i = 0; i < vp.GetCutawayPlanes().size(); ++i)
353 if (vp.GetCutawayPlanes()[i] != fVP.GetCutawayPlanes()[i])
354 return true;
355 }
356
357 if (vp.IsExplode () &&
358 (vp.GetExplodeFactor () != fVP.GetExplodeFactor ()))
359 return true;
360
361 if (vp.IsSpecialMeshRendering() &&
362 (vp.GetSpecialMeshVolumes() != fVP.GetSpecialMeshVolumes()))
363 return true;
364
365 return false;
366}
367
369{
370 fKeyPressed = true;
371 fKey = ev->key();
372}
373
374void G4Qt3DViewer::keyReleaseEvent(QKeyEvent* /*ev*/)
375{
376 fKeyPressed = false;
377}
378
379void G4Qt3DViewer::mouseDoubleClickEvent(QMouseEvent* /*ev*/) {}
380
381void G4Qt3DViewer::mouseMoveEvent(QMouseEvent* ev)
382{
383 // I think we only want these if a mouse button is pressed.
384 // But they come even when not pressed (on my MacBook Pro trackpad).
385 // Documentation says:
386 /* Mouse move events will occur only when a mouse button is pressed down,
387 unless mouse tracking has been enabled with QWidget::setMouseTracking().*/
388 // But this is a window not a widget.
389 // As a workaround we maintain a flag changed by mousePress/ReleaseEvent.
390#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
391 G4double x = ev->x();
392 G4double y = ev->y();
393#else
394 G4double x = ev->position().x();
395 G4double y = ev->position().y();
396#endif
399 fMousePressedX = x;
400 fMousePressedY = y;
401
402 if (fMousePressed) {
403
404 if (fKeyPressed && fKey == Qt::Key_Shift) { // Translation (pan)
405
406 const G4double sceneRadius = fQt3DSceneHandler.fpScene->GetExtent().GetExtentRadius();
407 const G4double scale = 300; // Roughly pixels per window, empirically chosen
408 const G4double dxScene = dx*sceneRadius/scale;
409 const G4double dyScene = dy*sceneRadius/scale;
410 fVP.IncrementPan(-dxScene,dyScene);
411
412 } else { // Rotation
413
414 // Simple ad-hoc algorithms
415 const G4Vector3D& x_prime = fVP.GetViewpointDirection().cross(fVP.GetUpVector());
416 const G4Vector3D& y_prime = x_prime.cross(fVP.GetViewpointDirection());
417 const G4double scale = 200; // Roughly pixels per window, empirically chosen
418 G4Vector3D newViewpointDirection = fVP.GetViewpointDirection();
419 newViewpointDirection += dx*x_prime/scale;
420 newViewpointDirection += dy*y_prime/scale;
421 fVP.SetViewpointDirection(newViewpointDirection.unit());
422
423 if (fVP.GetRotationStyle() == G4ViewParameters::freeRotation) {
424 G4Vector3D newUpVector = fVP.GetUpVector();
425 newUpVector += dx*x_prime/scale;
426 newUpVector += dy*y_prime/scale;
427 fVP.SetUpVector(newUpVector.unit());
428 }
429 }
430 }
431
432 SetView();
433 DrawView();
434}
435
436void G4Qt3DViewer::mousePressEvent(QMouseEvent* ev)
437{
438 fMousePressed = true;
439#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
440 fMousePressedX = ev->x();
441 fMousePressedY = ev->y();
442#else
443 fMousePressedX = ev->position().x();
444 fMousePressedY = ev->position().y();
445#endif
446}
447
448void G4Qt3DViewer::mouseReleaseEvent(QMouseEvent* /*ev*/)
449{
450 fMousePressed = false;
451}
452
453void G4Qt3DViewer::wheelEvent(QWheelEvent* ev)
454{
455 // Take note of up-down motion only
456 const G4double angleY = ev->angleDelta().y();
457
458 if (fVP.GetFieldHalfAngle() == 0.) { // Orthographic projection
459 const G4double scale = 500; // Empirically chosen
460 fVP.MultiplyZoomFactor(1.+angleY/scale);
461 } else { // Perspective projection
462 const G4double delta = fSceneHandler.GetExtent().GetExtentRadius()/200.; // Empirical
463 fVP.SetDolly(fVP.GetDolly()+angleY*delta);
464 }
465
466 SetView();
467 DrawView();
468}
G4TemplateAutoLock< G4Mutex > G4AutoLock
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
#define G4warn
Definition G4Scene.cc:41
#define G4CONDITION_INITIALIZER
#define G4MUTEX_INITIALIZER
#define G4CONDITIONWAITLAMBDA(cond, mutex, lambda)
G4int G4Condition
#define G4CONDITIONBROADCAST(cond)
std::mutex G4Mutex
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void mousePressEvent(QMouseEvent *)
void KernelVisitDecision()
G4Qt3DSceneHandler & fQt3DSceneHandler
G4double fMousePressedY
void Initialise()
G4bool fMousePressed
void SwitchToMasterThread()
void SwitchToVisSubThread()
void MovingToMasterThread()
G4double fMousePressedX
void mouseReleaseEvent(QMouseEvent *)
virtual ~G4Qt3DViewer()
G4bool fKeyPressed
virtual void resizeEvent(QResizeEvent *)
void wheelEvent(QWheelEvent *)
void keyPressEvent(QKeyEvent *)
QWidget * fUIWidget
G4Qt3DViewer(G4Qt3DSceneHandler &, const G4String &name)
void mouseDoubleClickEvent(QMouseEvent *)
void MovingToVisSubThread()
void mouseMoveEvent(QMouseEvent *)
G4bool CompareForKernelVisit(G4ViewParameters &)
G4ViewParameters fLastVP
void keyReleaseEvent(QKeyEvent *)
static G4UImanager * GetUIpointer()
G4bool fNeedKernelVisit
Definition G4VViewer.hh:265
void ProcessView()
Definition G4VViewer.cc:109
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:253
G4String fName
Definition G4VViewer.hh:255
void NeedKernelVisit()
Definition G4VViewer.cc:82
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:258
G4int fViewId
Definition G4VViewer.hh:254
G4ViewParameters fVP
Definition G4VViewer.hh:257
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition G4VViewer.cc:49
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
const G4Vector3D & GetScaleFactor() const
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
CutawayMode GetCutawayMode() const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsMarkerNotHidden() const
G4double GetGlobalLineWidthScale() const
G4bool IsCutaway() const
const G4Colour & GetBackgroundColour() const
G4bool IsSection() const
G4bool IsPicking() const
G4bool IsCulling() const
const G4VisAttributes * GetDefaultTextVisAttributes() const
G4bool IsExplode() 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
const G4Planes & GetCutawayPlanes() const
G4bool IsDensityCulling() const
G4double GetVisibleDensity() const
SMROption GetSpecialMeshRenderingOption() const
G4bool IsCullingCovered() const
const G4Plane3D & GetSectionPlane() const
DrawingStyle GetDrawingStyle() const
G4bool IsAuxEdgeVisible() const
const G4Colour & GetColour() const
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > unit() const
QColor ConvertToQColor(const G4Colour &c)
QVector3D ConvertToQVector3D(const G4ThreeVector &v)
G4bool IsMasterThread()