Geant4 11.1.1
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 Initialise ()
 
virtual void ResetView ()
 
virtual void SetView ()=0
 
virtual void ClearView ()=0
 
virtual void DrawView ()=0
 
void RefreshView ()
 
virtual void ShowView ()
 
virtual void FinishView ()
 
std::vector< G4ThreeVectorComputeFlyThrough (G4Vector3D *)
 
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

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.
 
G4bool fNeedKernelVisit
 

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
Definition: G4Qt3DViewer.hh:69
G4double fMousePressedY
Definition: G4Qt3DViewer.hh:76
G4bool fMousePressed
Definition: G4Qt3DViewer.hh:75
G4double fMousePressedX
Definition: G4Qt3DViewer.hh:76
G4bool fKeyPressed
Definition: G4Qt3DViewer.hh:73
G4int IncrementViewCount()

◆ ~G4Qt3DViewer()

G4Qt3DViewer::~G4Qt3DViewer ( )
virtual

Definition at line 74 of file G4Qt3DViewer.cc.

75{}

Member Function Documentation

◆ ClearView()

void G4Qt3DViewer::ClearView ( void  )
virtual

Implements G4VViewer.

Definition at line 128 of file G4Qt3DViewer.cc.

129{}

◆ CompareForKernelVisit()

G4bool G4Qt3DViewer::CompareForKernelVisit ( G4ViewParameters vp)
protected

Definition at line 242 of file G4Qt3DViewer.cc.

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}
G4ViewParameters fVP
Definition: G4VViewer.hh:220
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 131 of file G4Qt3DViewer.cc.

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}
bool G4bool
Definition: G4Types.hh:86
void KernelVisitDecision()
void FinishView()
G4ViewParameters fLastVP
Definition: G4Qt3DViewer.hh:68
G4bool fNeedKernelVisit
Definition: G4VViewer.hh:227
void ProcessView()
Definition: G4VViewer.cc:107

Referenced by mouseMoveEvent(), and wheelEvent().

◆ FinishView()

void G4Qt3DViewer::FinishView ( void  )
virtual

Reimplemented from G4VViewer.

Definition at line 167 of file G4Qt3DViewer.cc.

168{
170 show();
171 }
172}
G4bool IsMasterThread()
Definition: G4Threading.cc:124

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);
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}
#define G4warn
Definition: G4Scene.cc:41
#define G4endl
Definition: G4ios.hh:57
Qt3DCore::QEntity * fpQt3DScene
QWidget * fUIWidget
Definition: G4Qt3DViewer.hh:71
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4String fName
Definition: G4VViewer.hh:218
G4ViewParameters fDefaultVP
Definition: G4VViewer.hh:221
G4int fViewId
Definition: G4VViewer.hh:217
void SetAutoRefresh(G4bool)

◆ KernelVisitDecision()

void G4Qt3DViewer::KernelVisitDecision ( )
protected

Definition at line 232 of file G4Qt3DViewer.cc.

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

Referenced by DrawView().

◆ keyPressEvent()

void G4Qt3DViewer::keyPressEvent ( QKeyEvent *  ev)
protected

Definition at line 319 of file G4Qt3DViewer.cc.

320{
321 fKeyPressed = true;
322 fKey = ev->key();
323}

◆ keyReleaseEvent()

void G4Qt3DViewer::keyReleaseEvent ( QKeyEvent *  )
protected

Definition at line 325 of file G4Qt3DViewer.cc.

326{
327 fKeyPressed = false;
328}

◆ mouseDoubleClickEvent()

void G4Qt3DViewer::mouseDoubleClickEvent ( QMouseEvent *  )
protected

Definition at line 330 of file G4Qt3DViewer.cc.

330{}

◆ mouseMoveEvent()

void G4Qt3DViewer::mouseMoveEvent ( QMouseEvent *  ev)
protected

Definition at line 332 of file G4Qt3DViewer.cc.

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}
double G4double
Definition: G4Types.hh:83
void SetView()
Definition: G4Qt3DViewer.cc:77
const G4VisExtent & GetExtent() const
void SetViewpointDirection(const G4Vector3D &viewpointDirection)
const G4Vector3D & GetViewpointDirection() const
void IncrementPan(G4double right, G4double up)
const G4Vector3D & GetUpVector() const
void SetUpVector(const G4Vector3D &upVector)
RotationStyle GetRotationStyle() const
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > unit() const

◆ mousePressEvent()

void G4Qt3DViewer::mousePressEvent ( QMouseEvent *  ev)
protected

Definition at line 383 of file G4Qt3DViewer.cc.

384{
385 fMousePressed = true;
386 fMousePressedX = ev->x();
387 fMousePressedY = ev->y();
388}

◆ mouseReleaseEvent()

void G4Qt3DViewer::mouseReleaseEvent ( QMouseEvent *  )
protected

Definition at line 390 of file G4Qt3DViewer.cc.

391{
392 fMousePressed = false;
393}

◆ MovingToMasterThread()

void G4Qt3DViewer::MovingToMasterThread ( )

Definition at line 209 of file G4Qt3DViewer.cc.

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}
G4GLOB_DLL std::ostream G4cout

◆ MovingToVisSubThread()

void G4Qt3DViewer::MovingToVisSubThread ( )

Definition at line 181 of file G4Qt3DViewer.cc.

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}

◆ SetView()

void G4Qt3DViewer::SetView ( )
virtual

Implements G4VViewer.

Definition at line 77 of file G4Qt3DViewer.cc.

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}
const G4Point3D & GetStandardTargetPoint() const
G4Scene * GetScene() const
G4VSceneHandler & fSceneHandler
Definition: G4VViewer.hh:216
G4double GetCameraDistance(G4double radius) const
const G4Point3D & GetCurrentTargetPoint() const
G4double GetFarDistance(G4double cameraDistance, G4double nearDistance, G4double radius) const
G4double GetFieldHalfAngle() const
G4double GetFrontHalfHeight(G4double nearDistance, G4double radius) const
G4double GetNearDistance(G4double cameraDistance, G4double radius) const
QColor ConvertToQColor(const G4Colour &c)
Definition: G4Qt3DUtils.cc:46
QVector3D ConvertToQVector3D(const G4ThreeVector &v)
Definition: G4Qt3DUtils.cc:52

Referenced by mouseMoveEvent(), and wheelEvent().

◆ ShowView()

void G4Qt3DViewer::ShowView ( void  )
virtual

Reimplemented from G4VViewer.

Definition at line 157 of file G4Qt3DViewer.cc.

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}

◆ SwitchToMasterThread()

void G4Qt3DViewer::SwitchToMasterThread ( )

Definition at line 223 of file G4Qt3DViewer.cc.

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}

◆ SwitchToVisSubThread()

void G4Qt3DViewer::SwitchToVisSubThread ( )

Definition at line 197 of file G4Qt3DViewer.cc.

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}

◆ wheelEvent()

void G4Qt3DViewer::wheelEvent ( QWheelEvent *  ev)
protected

Definition at line 395 of file G4Qt3DViewer.cc.

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}
virtual const G4VisExtent & GetExtent() const
void SetDolly(G4double dolly)
void MultiplyZoomFactor(G4double zoomFactorMultiplier)
G4double GetDolly() const

Member Data Documentation

◆ fKey

int G4Qt3DViewer::fKey
protected

Definition at line 74 of file G4Qt3DViewer.hh.

Referenced by keyPressEvent(), and mouseMoveEvent().

◆ fKeyPressed

G4bool G4Qt3DViewer::fKeyPressed
protected

Definition at line 73 of file G4Qt3DViewer.hh.

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

◆ fLastVP

G4ViewParameters G4Qt3DViewer::fLastVP
protected

Definition at line 68 of file G4Qt3DViewer.hh.

Referenced by DrawView(), and KernelVisitDecision().

◆ fMousePressed

G4bool G4Qt3DViewer::fMousePressed
protected

Definition at line 75 of file G4Qt3DViewer.hh.

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

◆ fMousePressedX

G4double G4Qt3DViewer::fMousePressedX
protected

Definition at line 76 of file G4Qt3DViewer.hh.

Referenced by mouseMoveEvent(), and mousePressEvent().

◆ fMousePressedY

G4double G4Qt3DViewer::fMousePressedY
protected

Definition at line 76 of file G4Qt3DViewer.hh.

Referenced by mouseMoveEvent(), and mousePressEvent().

◆ fQt3DSceneHandler

G4Qt3DSceneHandler& G4Qt3DViewer::fQt3DSceneHandler
protected

◆ fUIWidget

QWidget* G4Qt3DViewer::fUIWidget
protected

Definition at line 71 of file G4Qt3DViewer.hh.

Referenced by Initialise(), and SetView().


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