Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Qt3DViewer Class Reference

#include <G4Qt3DViewer.hh>

+ Inheritance diagram for G4Qt3DViewer:

Public Member Functions

 G4Qt3DViewer (G4Qt3DSceneHandler &, const G4String &name)
 
virtual ~G4Qt3DViewer ()
 
void Initialise ()
 
void SetView ()
 
void ClearView ()
 
void DrawView ()
 
void ShowView ()
 
void FinishView ()
 
void MovingToVisSubThread ()
 
void SwitchToVisSubThread ()
 
void MovingToMasterThread ()
 
void SwitchToMasterThread ()
 
- Public Member Functions inherited from G4VViewer
 G4VViewer (G4VSceneHandler &, G4int id, const G4String &name="")
 
virtual ~G4VViewer ()
 
virtual void ResetView ()
 
void RefreshView ()
 
virtual G4bool ReadyToDraw ()
 
std::vector< G4ThreeVectorComputeFlyThrough (G4Vector3D *)
 
virtual void DoneWithMasterThread ()
 
virtual void DoneWithVisSubThread ()
 
void InsertModelInSceneTree (G4VModel *)
 
const G4SceneTreeItemGetSceneTree ()
 
G4SceneTreeItemAccessSceneTree ()
 
void UpdateGUISceneTree ()
 
const G4StringGetName () const
 
const G4StringGetShortName () const
 
void SetName (const G4String &)
 
G4int GetViewId () const
 
G4VSceneHandlerGetSceneHandler () const
 
const G4ViewParametersGetViewParameters () const
 
const G4ViewParametersGetDefaultViewParameters () const
 
G4double GetKernelVisitElapsedTimeSeconds () const
 
virtual const std::vector< G4ModelingParameters::VisAttributesModifier > * GetPrivateVisAttributesModifiers () const
 
void SetViewParameters (const G4ViewParameters &vp)
 
void SetDefaultViewParameters (const G4ViewParameters &vp)
 
const G4VisAttributesGetApplicableVisAttributes (const G4VisAttributes *) const
 
void SetNeedKernelVisit (G4bool need)
 
void NeedKernelVisit ()
 
void ProcessView ()
 

Protected Member Functions

virtual void resizeEvent (QResizeEvent *)
 
void KernelVisitDecision ()
 
G4bool CompareForKernelVisit (G4ViewParameters &)
 
void keyPressEvent (QKeyEvent *)
 
void keyReleaseEvent (QKeyEvent *)
 
void mouseDoubleClickEvent (QMouseEvent *)
 
void mouseMoveEvent (QMouseEvent *)
 
void mousePressEvent (QMouseEvent *)
 
void mouseReleaseEvent (QMouseEvent *)
 
void wheelEvent (QWheelEvent *)
 
- Protected Member Functions inherited from G4VViewer
void SetTouchable (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath)
 
void TouchableSetVisibility (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, G4bool visibility)
 
void TouchableSetColour (const std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > &fullPath, const G4Colour &)
 

Protected Attributes

G4ViewParameters fLastVP
 
G4Qt3DSceneHandlerfQt3DSceneHandler
 
QWidget * fUIWidget
 
G4bool fKeyPressed
 
int fKey
 
G4bool fMousePressed
 
G4double fMousePressedX
 
G4double fMousePressedY
 
- Protected Attributes inherited from G4VViewer
G4VSceneHandlerfSceneHandler
 
G4int fViewId
 
G4String fName
 
G4String fShortName
 
G4ViewParameters fVP
 
G4ViewParameters fDefaultVP
 
G4double fKernelVisitElapsedTimeSeconds = 999.
 
G4SceneTreeItem fSceneTree
 
G4bool fNeedKernelVisit
 

Additional Inherited Members

- Public Attributes inherited from G4VViewer
const G4int fMaxNTouchables = 10000
 
G4bool fCurtailDescent = false
 

Detailed Description

Definition at line 37 of file G4Qt3DViewer.hh.

Constructor & Destructor Documentation

◆ G4Qt3DViewer()

G4Qt3DViewer::G4Qt3DViewer ( G4Qt3DSceneHandler & sceneHandler,
const G4String & name )

Definition at line 40 of file G4Qt3DViewer.cc.

42: G4VViewer(sceneHandler, sceneHandler.IncrementViewCount(), name)
43, fQt3DSceneHandler(sceneHandler)
44, fKeyPressed(false)
45, fMousePressed(false)
48{}
G4Qt3DSceneHandler & fQt3DSceneHandler
G4double fMousePressedY
G4bool fMousePressed
G4double fMousePressedX
G4bool fKeyPressed
G4int IncrementViewCount()
G4VViewer(G4VSceneHandler &, G4int id, const G4String &name="")
Definition G4VViewer.cc:49

◆ ~G4Qt3DViewer()

G4Qt3DViewer::~G4Qt3DViewer ( )
virtual

Definition at line 71 of file G4Qt3DViewer.cc.

72{
73 setRootEntity(nullptr);
74}

Member Function Documentation

◆ ClearView()

void G4Qt3DViewer::ClearView ( void )
virtual

Implements G4VViewer.

Definition at line 131 of file G4Qt3DViewer.cc.

132{}

◆ CompareForKernelVisit()

G4bool G4Qt3DViewer::CompareForKernelVisit ( G4ViewParameters & vp)
protected

Definition at line 291 of file G4Qt3DViewer.cc.

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}
G4ViewParameters fVP
Definition G4VViewer.hh:257
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

Referenced by KernelVisitDecision().

◆ DrawView()

void G4Qt3DViewer::DrawView ( )
virtual

Implements G4VViewer.

Definition at line 134 of file G4Qt3DViewer.cc.

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}
bool G4bool
Definition G4Types.hh:86
void KernelVisitDecision()
G4ViewParameters fLastVP
G4bool fNeedKernelVisit
Definition G4VViewer.hh:265
void ProcessView()
Definition G4VViewer.cc:109

Referenced by mouseMoveEvent(), and wheelEvent().

◆ FinishView()

void G4Qt3DViewer::FinishView ( void )
virtual

Reimplemented from G4VViewer.

Definition at line 172 of file G4Qt3DViewer.cc.

173{
174#if QT_VERSION < 0x060000
176 show();
177 }
178#endif
179}
G4bool IsMasterThread()

Referenced by DrawView().

◆ Initialise()

void G4Qt3DViewer::Initialise ( )
virtual

Reimplemented from G4VViewer.

Definition at line 50 of file G4Qt3DViewer.cc.

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}
#define G4warn
Definition G4Scene.cc:41
#define G4endl
Definition G4ios.hh:67
QWidget * fUIWidget
static G4UImanager * GetUIpointer()
G4String fName
Definition G4VViewer.hh:255
G4ViewParameters fDefaultVP
Definition G4VViewer.hh:258
G4int fViewId
Definition G4VViewer.hh:254

◆ KernelVisitDecision()

void G4Qt3DViewer::KernelVisitDecision ( )
protected

Definition at line 281 of file G4Qt3DViewer.cc.

281 {
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}
G4bool CompareForKernelVisit(G4ViewParameters &)
void NeedKernelVisit()
Definition G4VViewer.cc:82

Referenced by DrawView().

◆ keyPressEvent()

void G4Qt3DViewer::keyPressEvent ( QKeyEvent * ev)
protected

Definition at line 368 of file G4Qt3DViewer.cc.

369{
370 fKeyPressed = true;
371 fKey = ev->key();
372}

◆ keyReleaseEvent()

void G4Qt3DViewer::keyReleaseEvent ( QKeyEvent * )
protected

Definition at line 374 of file G4Qt3DViewer.cc.

375{
376 fKeyPressed = false;
377}

◆ mouseDoubleClickEvent()

void G4Qt3DViewer::mouseDoubleClickEvent ( QMouseEvent * )
protected

Definition at line 379 of file G4Qt3DViewer.cc.

379{}

◆ mouseMoveEvent()

void G4Qt3DViewer::mouseMoveEvent ( QMouseEvent * ev)
protected

Definition at line 381 of file G4Qt3DViewer.cc.

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}
double G4double
Definition G4Types.hh:83
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > unit() const

◆ mousePressEvent()

void G4Qt3DViewer::mousePressEvent ( QMouseEvent * ev)
protected

Definition at line 436 of file G4Qt3DViewer.cc.

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}

◆ mouseReleaseEvent()

void G4Qt3DViewer::mouseReleaseEvent ( QMouseEvent * )
protected

Definition at line 448 of file G4Qt3DViewer.cc.

449{
450 fMousePressed = false;
451}

◆ MovingToMasterThread()

void G4Qt3DViewer::MovingToMasterThread ( )
virtual

Reimplemented from G4VViewer.

Definition at line 256 of file G4Qt3DViewer.cc.

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}

◆ MovingToVisSubThread()

void G4Qt3DViewer::MovingToVisSubThread ( )
virtual

Reimplemented from G4VViewer.

Definition at line 206 of file G4Qt3DViewer.cc.

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}
G4TemplateAutoLock< G4Mutex > G4AutoLock
#define G4CONDITIONWAITLAMBDA(cond, mutex, lambda)
#define G4CONDITIONBROADCAST(cond)

◆ resizeEvent()

void G4Qt3DViewer::resizeEvent ( QResizeEvent * )
protectedvirtual

Definition at line 76 of file G4Qt3DViewer.cc.

76 {
77 SetView();
78}

◆ SetView()

void G4Qt3DViewer::SetView ( )
virtual

Implements G4VViewer.

Definition at line 80 of file G4Qt3DViewer.cc.

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}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
G4GLOB_DLL std::ostream G4cout
G4VSceneHandler & fSceneHandler
Definition G4VViewer.hh:253
QColor ConvertToQColor(const G4Colour &c)
QVector3D ConvertToQVector3D(const G4ThreeVector &v)

Referenced by mouseMoveEvent(), resizeEvent(), and wheelEvent().

◆ ShowView()

void G4Qt3DViewer::ShowView ( void )
virtual

Reimplemented from G4VViewer.

Definition at line 160 of file G4Qt3DViewer.cc.

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}

◆ SwitchToMasterThread()

void G4Qt3DViewer::SwitchToMasterThread ( )
virtual

Reimplemented from G4VViewer.

Definition at line 275 of file G4Qt3DViewer.cc.

277{
278 visSubThreadEstablished = false;
279}

◆ SwitchToVisSubThread()

void G4Qt3DViewer::SwitchToVisSubThread ( )
virtual

Reimplemented from G4VViewer.

Definition at line 236 of file G4Qt3DViewer.cc.

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}

◆ wheelEvent()

void G4Qt3DViewer::wheelEvent ( QWheelEvent * ev)
protected

Definition at line 453 of file G4Qt3DViewer.cc.

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}

Member Data Documentation

◆ fKey

int G4Qt3DViewer::fKey
protected

Definition at line 75 of file G4Qt3DViewer.hh.

Referenced by keyPressEvent(), and mouseMoveEvent().

◆ fKeyPressed

G4bool G4Qt3DViewer::fKeyPressed
protected

Definition at line 74 of file G4Qt3DViewer.hh.

Referenced by G4Qt3DViewer(), keyPressEvent(), keyReleaseEvent(), and mouseMoveEvent().

◆ fLastVP

G4ViewParameters G4Qt3DViewer::fLastVP
protected

Definition at line 69 of file G4Qt3DViewer.hh.

Referenced by DrawView(), and KernelVisitDecision().

◆ fMousePressed

G4bool G4Qt3DViewer::fMousePressed
protected

Definition at line 76 of file G4Qt3DViewer.hh.

Referenced by G4Qt3DViewer(), mouseMoveEvent(), mousePressEvent(), and mouseReleaseEvent().

◆ fMousePressedX

G4double G4Qt3DViewer::fMousePressedX
protected

Definition at line 77 of file G4Qt3DViewer.hh.

Referenced by G4Qt3DViewer(), mouseMoveEvent(), and mousePressEvent().

◆ fMousePressedY

G4double G4Qt3DViewer::fMousePressedY
protected

Definition at line 77 of file G4Qt3DViewer.hh.

Referenced by G4Qt3DViewer(), mouseMoveEvent(), and mousePressEvent().

◆ fQt3DSceneHandler

G4Qt3DSceneHandler& G4Qt3DViewer::fQt3DSceneHandler
protected

◆ fUIWidget

QWidget* G4Qt3DViewer::fUIWidget
protected

Definition at line 72 of file G4Qt3DViewer.hh.

Referenced by Initialise(), and SetView().


The documentation for this class was generated from the following files: