Geant4 11.1.1
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)
46, fMousePressedX(0.)
47, fMousePressedY(0.)
48{}
49
51{
52 setObjectName(fName.c_str());
53
54 fVP.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 fUIWidget->setMinimumSize(QSize(200, 100));
67 fUIWidget->setMaximumSize(screen()->size());
68// fUIWidget->setFocusPolicy(Qt::NoFocus); //??
69 uiQt->AddTabWidget(fUIWidget,QString(fName));
70
71 setRootEntity(fQt3DSceneHandler.fpQt3DScene);
72}
73
75{}
76
78{
79 // Background colour
80 defaultFrameGraph()->setClearColor(G4Qt3DUtils::ConvertToQColor(fVP.GetBackgroundColour()));
81
82 // Get radius of scene, etc.
83 // Note that this procedure properly takes into account zoom, dolly and pan.
84 const G4Point3D targetPoint
88 if(radius<=0.) radius = 1.;
89 const G4double cameraDistance = fVP.GetCameraDistance (radius);
90 const G4Point3D cameraPosition =
91 targetPoint + cameraDistance * fVP.GetViewpointDirection().unit();
92 const GLdouble pnear = fVP.GetNearDistance (cameraDistance, radius);
93 const GLdouble pfar = fVP.GetFarDistance (cameraDistance, pnear, radius);
94 const GLdouble right = fVP.GetFrontHalfHeight (pnear, radius);
95 const GLdouble left = -right;
96 const GLdouble top = fVP.GetFrontHalfHeight (pnear, radius);
97 const GLdouble bottom = -top;
98
99 camera()->setObjectName((fName + " camera").c_str());
100 camera()->setViewCenter(G4Qt3DUtils::ConvertToQVector3D(targetPoint));
101 camera()->setPosition(G4Qt3DUtils::ConvertToQVector3D(cameraPosition));
102 camera()->setUpVector(G4Qt3DUtils::ConvertToQVector3D(fVP.GetUpVector()));
103
104// auto lightEntity = new Qt3DCore::QEntity(fQt3DSceneHandler.fpQt3DScene);
105// auto directionalLight = new Qt3DRender::QDirectionalLight(lightEntity);
106//// directionalLight->setColor("white");
107//// directionalLight->setIntensity(1.);
108// directionalLight->setWorldDirection(G4Qt3DUtils::ConvertToQVector3D(fVP.GetActualLightpointDirection()));
109// lightEntity->addComponent(directionalLight);
110
111 const auto& size = fUIWidget->size();
112 G4double w = size.width();
113 G4double h = size.height();
114#ifdef G4QT3DDEBUG
115 // Curiously w,h are wrong first time - 640,480 instead of (my Mac) 991,452.
116 G4cout << "W,H: " << w << ',' << h << G4endl;
117#endif
118 const G4double aspectRatio = w/h;
119 if (fVP.GetFieldHalfAngle() == 0.) {
120 camera()->lens()->setOrthographicProjection
121 (left*aspectRatio,right*aspectRatio,bottom,top,pnear,pfar);
122 } else {
123 camera()->lens()->setPerspectiveProjection
124 (2.*fVP.GetFieldHalfAngle()/deg,aspectRatio,pnear,pfar);
125 }
126}
127
129{}
130
132{
133 // First, a view should decide when to re-visit the G4 kernel.
134 // Sometimes it might not be necessary, e.g., if the scene is stored
135 // in a graphical database (e.g., OpenGL's display lists) and only
136 // the viewing angle has changed. But graphics systems without a
137 // graphical database will always need to visit the G4 kernel.
138
139 // The fNeedKernelVisit flag might have been set by the user in
140 // /vis/viewer/rebuild, but if not, make decision and set flag only
141 // if necessary...
143 G4bool kernelVisitWasNeeded = fNeedKernelVisit; // Keep (ProcessView resets).
144 fLastVP = fVP;
145
146 ProcessView (); // Clears store and processes scene only if necessary.
147
148 if (kernelVisitWasNeeded) {
149 // We might need to do something if the kernel was visited.
150 } else {
151 }
152
153 // ...before finally...
154 FinishView (); // Flush streams and/or swap buffers.
155}
156
158{
159 // show() may only be called from master thread
161 show();
162 }
163 // The way Qt seems to work, we don't seem to need a show() anyway, but
164 // we'll leave it in - it seems not to have any effect, good or bad.
165}
166
168{
170 show();
171 }
172}
173
174namespace {
175 QThread* masterQThread = nullptr;
176 QThread* visSubThreadQThread = nullptr;
177}
178
179# include <chrono>
180
182{
183#ifdef G4QT3DDEBUG
184 G4cout << "G4Qt3DViewer::MovingToVisSubThread" << G4endl;
185#endif
186 // Still on master thread but vis thread has been launched
187 // Make note of master QThread
188 masterQThread = QThread::currentThread();
189 // Wait until SwitchToVisSubThread has found vis sub-thread QThread
190 std::this_thread::sleep_for(std::chrono::milliseconds(100));
191 // Move relevant stuff to vis sub-thread QThread
192 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
193 auto p2 = p1->parent();
194 p2->moveToThread(visSubThreadQThread);
195}
196
198{
199#ifdef G4QT3DDEBUG
200 G4cout << "G4Qt3DViewer::SwitchToVisSubThread" << G4endl;
201#endif
202 // On vis sub-thread before any drawing
203 // Make note of vis-subthread QThread for MovingToVisSubThread
204 visSubThreadQThread = QThread::currentThread();
205 // Wait until SwitchToVisSubThread has moved stuff
206 std::this_thread::sleep_for(std::chrono::milliseconds(1000));
207}
208
210{
211#ifdef G4QT3DDEBUG
212 G4cout << "G4Qt3DViewer::MovingToMasterThread" << G4endl;
213#endif
214 // On vis sub-thread just before exit
215 // Move relevant stuff to master QThread.
216 auto p1 = fQt3DSceneHandler.fpQt3DScene->parent();
217 auto p2 = p1->parent();
218 p2->moveToThread(masterQThread);
219 // Zero - will be different next run
220 visSubThreadQThread = nullptr;
221}
222
224{
225#ifdef G4QT3DDEBUG
226 G4cout << "G4Qt3DViewer::SwitchToMasterThread" << G4endl;
227#endif
228 // On master thread after vis sub-thread has terminated
229 // Nothing to do
230}
231
233
234 // If there's a significant difference with the last view parameters
235 // of either the scene handler or this viewer, trigger a rebuild.
236
238 NeedKernelVisit (); // Sets fNeedKernelVisit.
239 }
240}
241
243{
244 // Typical comparison. Taken from OpenInventor.
245 if (
246 (vp.GetDrawingStyle () != fVP.GetDrawingStyle ()) ||
248 (vp.IsAuxEdgeVisible () != fVP.IsAuxEdgeVisible ()) ||
249 (vp.IsCulling () != fVP.IsCulling ()) ||
251 (vp.IsDensityCulling () != fVP.IsDensityCulling ()) ||
252 (vp.IsCullingCovered () != fVP.IsCullingCovered ()) ||
253 (vp.GetCBDAlgorithmNumber() !=
255 (vp.IsSection () != fVP.IsSection ()) ||
256 (vp.IsCutaway () != fVP.IsCutaway ()) ||
257 // This assumes use of generic clipping (sectioning, slicing,
258 // DCUT, cutaway). If a decision is made to implement locally,
259 // this will need changing. See G4OpenGLViewer::SetView,
260 // G4OpenGLStoredViewer.cc::CompareForKernelVisit and
261 // G4OpenGLStoredSceneHander::CreateSection/CutawayPolyhedron.
262 (vp.IsExplode () != fVP.IsExplode ()) ||
263 (vp.GetNoOfSides () != fVP.GetNoOfSides ()) ||
272 (vp.IsPicking () != fVP.IsPicking ()) ||
273 // Scaling for Open Inventor is done by the scene handler so it
274 // needs a kernel visit. (In this respect, it differs from the
275 // OpenGL drivers, where it's done in SetView.)
276 (vp.GetScaleFactor () != fVP.GetScaleFactor ()) ||
283 )
284 return true;
285
286 if (vp.IsDensityCulling () &&
288 return true;
289
290 if (vp.GetCBDAlgorithmNumber() > 0) {
291 if (vp.GetCBDParameters().size() != fVP.GetCBDParameters().size()) return true;
292 else if (vp.GetCBDParameters() != fVP.GetCBDParameters()) return true;
293 }
294
295 if (vp.IsSection () &&
296 (vp.GetSectionPlane () != fVP.GetSectionPlane ()))
297 return true;
298
299 if (vp.IsCutaway ()) {
300 if (vp.GetCutawayMode() != fVP.GetCutawayMode()) return true;
301 if (vp.GetCutawayPlanes ().size () !=
302 fVP.GetCutawayPlanes ().size ()) return true;
303 for (size_t i = 0; i < vp.GetCutawayPlanes().size(); ++i)
304 if (vp.GetCutawayPlanes()[i] != fVP.GetCutawayPlanes()[i])
305 return true;
306 }
307
308 if (vp.IsExplode () &&
310 return true;
311
312 if (vp.IsSpecialMeshRendering() &&
314 return true;
315
316 return false;
317}
318
320{
321 fKeyPressed = true;
322 fKey = ev->key();
323}
324
325void G4Qt3DViewer::keyReleaseEvent(QKeyEvent* /*ev*/)
326{
327 fKeyPressed = false;
328}
329
330void G4Qt3DViewer::mouseDoubleClickEvent(QMouseEvent* /*ev*/) {}
331
332void G4Qt3DViewer::mouseMoveEvent(QMouseEvent* ev)
333{
334 // I think we only want these if a mouse button is pressed.
335 // But they come even when not pressed (on my MacBook Pro trackpad).
336 // Documentation says:
337 /* Mouse move events will occur only when a mouse button is pressed down,
338 unless mouse tracking has been enabled with QWidget::setMouseTracking().*/
339 // But this is a window not a widget.
340 // As a workaround we maintain a flag changed by mousePress/ReleaseEvent.
341
342 G4double x = ev->x();
343 G4double y = ev->y();
346 fMousePressedX = x;
347 fMousePressedY = y;
348
349 if (fMousePressed) {
350
351 if (fKeyPressed && fKey == Qt::Key_Shift) { // Translation (pan)
352
354 const G4double scale = 300; // Roughly pixels per window, empirically chosen
355 const G4double dxScene = dx*sceneRadius/scale;
356 const G4double dyScene = dy*sceneRadius/scale;
357 fVP.IncrementPan(-dxScene,dyScene);
358
359 } else { // Rotation
360
361 // Simple ad-hoc algorithms
363 const G4Vector3D& y_prime = x_prime.cross(fVP.GetViewpointDirection());
364 const G4double scale = 200; // Roughly pixels per window, empirically chosen
365 G4Vector3D newViewpointDirection = fVP.GetViewpointDirection();
366 newViewpointDirection += dx*x_prime/scale;
367 newViewpointDirection += dy*y_prime/scale;
368 fVP.SetViewpointDirection(newViewpointDirection.unit());
369
371 G4Vector3D newUpVector = fVP.GetUpVector();
372 newUpVector += dx*x_prime/scale;
373 newUpVector += dy*y_prime/scale;
374 fVP.SetUpVector(newUpVector.unit());
375 }
376 }
377 }
378
379 SetView();
380 DrawView();
381}
382
383void G4Qt3DViewer::mousePressEvent(QMouseEvent* ev)
384{
385 fMousePressed = true;
386 fMousePressedX = ev->x();
387 fMousePressedY = ev->y();
388}
389
390void G4Qt3DViewer::mouseReleaseEvent(QMouseEvent* /*ev*/)
391{
392 fMousePressed = false;
393}
394
395void G4Qt3DViewer::wheelEvent(QWheelEvent* ev)
396{
397 // Take note of up-down motion only
398 const G4double angleY = ev->angleDelta().y();
399
400 if (fVP.GetFieldHalfAngle() == 0.) { // Orthographic projection
401 const G4double scale = 500; // Empirically chosen
402 fVP.MultiplyZoomFactor(1.+angleY/scale);
403 } else { // Perspective projection
404 const G4double delta = fSceneHandler.GetExtent().GetExtentRadius()/200.; // Empirical
405 fVP.SetDolly(fVP.GetDolly()+angleY*delta);
406 }
407
408 SetView();
409 DrawView();
410}
#define G4warn
Definition: G4Scene.cc:41
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#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:50
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:74
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:41
void mouseDoubleClickEvent(QMouseEvent *)
void MovingToVisSubThread()
void mouseMoveEvent(QMouseEvent *)
void FinishView()
G4bool CompareForKernelVisit(G4ViewParameters &)
G4ViewParameters fLastVP
Definition: G4Qt3DViewer.hh:68
void SetView()
Definition: G4Qt3DViewer.cc:77
void keyReleaseEvent(QKeyEvent *)
const G4VisExtent & GetExtent() const
const G4Point3D & GetStandardTargetPoint() const
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4Scene * GetScene() const
virtual const G4VisExtent & GetExtent() const
G4bool fNeedKernelVisit
Definition: G4VViewer.hh:227
void ProcessView()
Definition: G4VViewer.cc:107
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:216
G4String fName
Definition: G4VViewer.hh:218
void NeedKernelVisit()
Definition: G4VViewer.cc:80
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:221
G4int fViewId
Definition: G4VViewer.hh:217
G4ViewParameters fVP
Definition: G4VViewer.hh:220
void SetViewpointDirection(const G4Vector3D &viewpointDirection)
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
const G4Vector3D & GetScaleFactor() const
void SetAutoRefresh(G4bool)
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
CutawayMode GetCutawayMode() 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)
const G4Planes & GetCutawayPlanes() const
RotationStyle GetRotationStyle() const
G4bool IsDensityCulling() const
void MultiplyZoomFactor(G4double zoomFactorMultiplier)
G4double GetVisibleDensity() const
SMROption GetSpecialMeshRenderingOption() const
G4bool IsCullingCovered() const
const G4Plane3D & GetSectionPlane() 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
QColor ConvertToQColor(const G4Colour &c)
Definition: G4Qt3DUtils.cc:46
QVector3D ConvertToQVector3D(const G4ThreeVector &v)
Definition: G4Qt3DUtils.cc:52
G4bool IsMasterThread()
Definition: G4Threading.cc:124