Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/visualization/Qt3D/src/G4Qt3DViewer.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 40 G4Qt3DViewer::G4Qt3DViewer
 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 
 50 void G4Qt3DViewer::Initialise()
 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 
 71 G4Qt3DViewer::~G4Qt3DViewer()
 72 {
 73   setRootEntity(nullptr);
 74 }
 75 
 76 void G4Qt3DViewer::resizeEvent(QResizeEvent*) {
 77   SetView();
 78 }
 79 
 80 void G4Qt3DViewer::SetView()
 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 
131 void G4Qt3DViewer::ClearView()
132 {}
133 
134 void G4Qt3DViewer::DrawView()
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...
145   if (!fNeedKernelVisit) KernelVisitDecision();
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 
160 void G4Qt3DViewer::ShowView()
161 {
162 #if QT_VERSION < 0x060000
163   // show() may only be called from master thread
164   if (G4Threading::IsMasterThread()) {
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 
172 void G4Qt3DViewer::FinishView()
173 {
174 #if QT_VERSION < 0x060000
175   if (G4Threading::IsMasterThread()) {
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 
196 namespace {
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 
206 void G4Qt3DViewer::MovingToVisSubThread()
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 
236 void G4Qt3DViewer::SwitchToVisSubThread()
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 
256 void G4Qt3DViewer::MovingToMasterThread()
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 
275 void G4Qt3DViewer::SwitchToMasterThread()
276 // On master thread after vis sub-thread has terminated
277 {
278   visSubThreadEstablished = false;
279 }
280 
281 void G4Qt3DViewer::KernelVisitDecision () {
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 
286   if (CompareForKernelVisit(fLastVP)) {
287     NeedKernelVisit ();  // Sets fNeedKernelVisit.
288   }
289 }
290 
291 G4bool G4Qt3DViewer::CompareForKernelVisit(G4ViewParameters& vp)
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 ())  ||
316      (vp.GetDefaultVisAttributes()->GetColour() !=
317       fVP.GetDefaultVisAttributes()->GetColour())            ||
318      (vp.GetDefaultTextVisAttributes()->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 ())     ||
326      (vp.GetVisAttributesModifiers() !=
327       fVP.GetVisAttributesModifiers())                       ||
328      (vp.IsSpecialMeshRendering() !=
329       fVP.IsSpecialMeshRendering())                          ||
330      (vp.GetSpecialMeshRenderingOption() !=
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 
368 void G4Qt3DViewer::keyPressEvent(QKeyEvent* ev)
369 {
370   fKeyPressed = true;
371   fKey = ev->key();
372 }
373 
374 void G4Qt3DViewer::keyReleaseEvent(QKeyEvent* /*ev*/)
375 {
376   fKeyPressed = false;
377 }
378 
379 void G4Qt3DViewer::mouseDoubleClickEvent(QMouseEvent* /*ev*/) {}
380 
381 void 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
397   G4double dx = x-fMousePressedX;
398   G4double dy = y-fMousePressedY;
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 
436 void 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 
448 void G4Qt3DViewer::mouseReleaseEvent(QMouseEvent* /*ev*/)
449 {
450   fMousePressed = false;
451 }
452 
453 void 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 }
469