Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpenInventorQtExaminerViewer.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
27// Frederick Jones TRIUMF 07 January 2018
28
30
31#include "ui_OIQtListsDialog.h"
32
33#include "saveViewPt.h"
34#include "pickext.h"
35#include "pickref.h"
36#include "wireframe.h"
37
38#include <algorithm> // For using sort on a vector
39
40#include "G4ios.hh"
41#include "G4UImanager.hh"
42#include "G4UIQt.hh"
43
44#include <Inventor/Qt/SoQt.h>
45#include <Inventor/Qt/SoQtCursor.h>
46#include <Inventor/events/SoKeyboardEvent.h>
47#include <Inventor/events/SoMouseButtonEvent.h>
48#include <Inventor/events/SoLocation2Event.h>
49#include <Inventor/nodes/SoSeparator.h>
50#include <Inventor/nodes/SoOrthographicCamera.h>
51#include <Inventor/nodes/SoPerspectiveCamera.h>
52
53// FWJ moved to header file
54//#include <Inventor/nodes/SoEventCallback.h>
55#include <Inventor/nodes/SoLineSet.h>
56#include <Inventor/nodes/SoMaterial.h>
57#include <Inventor/errors/SoDebugError.h>
58#include <Inventor/SoPickedPoint.h>
59#include <Inventor/actions/SoWriteAction.h>
60#include <Inventor/projectors/SbPlaneProjector.h>
61
62#include <Inventor/sensors/SoTimerSensor.h> // Animation
63#include <Inventor/sensors/SoNodeSensor.h> // Detect start of run
64
65#include "Geant4_SoPolyhedron.h"
66#include "G4TrajectoryPoint.hh"
67#include "G4AttHolder.hh"
68#include "G4AttCheck.hh"
69
70#include <Inventor/nodes/SoCallback.h>
71#include <Inventor/nodes/SoSwitch.h>
72#include <Inventor/nodes/SoScale.h>
73#include <Inventor/nodes/SoTranslation.h>
74#include <Inventor/actions/SoSearchAction.h>
75#include <Inventor/actions/SoGetBoundingBoxAction.h>
76
77#include <Inventor/nodes/SoCoordinate3.h>
78// For rendering distance during animation:
79#include <Inventor/nodes/SoText2.h>
80#include <Inventor/nodes/SoFont.h>
81#include <Inventor/nodes/SoPointSet.h>
82#include <Inventor/nodes/SoDrawStyle.h>
83#include <Inventor/nodes/SoBaseColor.h>
84
85// For searching for nodes within kits:
86#include <Inventor/nodekits/SoBaseKit.h>
87
88#include <QMenuBar>
89#include <QPushButton>
90#include <QRadioButton>
91#include <QToolButton>
92#include <QListWidget>
93#include <QListWidgetItem>
94#include <QInputDialog>
95#include <QMessageBox>
96#include <QFileDialog>
97#include <QStyle>
98#include <QCommonStyle>
99//#include <QMainWindow>
100
101#ifndef G4GMAKE
102#include "moc_G4OpenInventorQtExaminerViewer.cpp"
103#endif
104
105#define G4warn G4cout
106
107G4OpenInventorQtExaminerViewer* G4OpenInventorQtExaminerViewer::viewer = 0;
108
109#define MIN_SPEED 2.1 // Lower number means faster
110#define START_STEP 0.3
111#define SPEED_INDICATOR_STEP 0.045
112#define MAX_SPEED_INDICATOR 0.81
113// Number of steps 90 degree rotation around an element is split into
114#define ROT_CNT 6
115
116
117// Constructor
119G4OpenInventorQtExaminerViewer(QWidget* parent, const char* name, SbBool embed,
120 SoQtFullViewer::BuildFlag flag,
121 SoQtViewer::Type type)
122 : SoQtExaminerViewer(parent, name, embed, flag, type),
123 externalQtApp(0), processSoEventCount(0)
124{
125 // FWJ DEBUG
126 // G4cout << "G4OpenInventorQtExaminerViewer CONSTRUCTOR CALLED" << G4endl;
127 // G4cout << "G4OpenInventorQtExaminerViewer parent=" << parent << G4endl;
128
129 // FWJ THIS DOESN'T WORK APPARENTLY NO MAINWINDOW
130 // QMenuBar* menubar = ((QMainWindow*)parent)->menuBar();
131
132 fName = new QString(name);
133 viewer = this;
135}
136
137// Destructor
139{
140 // if (superimposition != NULL) {
141 // removeSuperimposition(superimposition);
142 // superimposition->unref();
143 // superimposition = NULL;
144 // }
145 // if (animateSensor->isScheduled())
146 // animateSensor->unschedule();
147 // delete animateSensor;
148 // delete sceneChangeSensor;
149 // delete[] curViewPtName;
150 // delete searcher;
151
152 viewer = 0;
153}
154
155
157{
158 setFeedbackSize(40);
159
160 hookBeamOn = new HookEventProcState(this);
161 newEvents = false;
162
163 buildWidget(getParentWidget());
164
165 fileName = "bookmarkFile"; // Default viewpoint file name
166 viewPtIdx = -1; // index of the most recent viewpoint in viewPtList vector
167
168 animateSensor = new SoTimerSensor(animateSensorCB, this);
169 animateSensorRotation = new SoTimerSensor(animateSensorRotationCB, this);
171
173 myCam = new SoPerspectiveCamera;
174 MAX_VP_IDX = 3;
175 MAX_VP_NAME = 35; // Max length of a viewpoint name, padded with spaces
176 curViewPtName = new char[MAX_VP_NAME + 1];
177 left_right = up_down = 0; // For movements around the beam during animation
178 speedStep = START_STEP; // For smoother animation speed increase/decrease
179 rotUpVec = false; // Used during scene element rotations
180 step = 1; //By default
181 // Used for moving along the beam with the
182 // mouse instead of rotating the view
183 lshiftdown = rshiftdown = false;
184 // Used for rotating the view with the camera
185 // staying in place
186 lctrldown = rctrldown = false;
187 // Used to send abbreviated output to the console when
188 abbrOutputFlag = false;
189 pickRefPathFlag = false;
190 prevColorField = NULL;
191 // warningFlag = false; // We come from the warning dialog
192 // myElementList = NULL;
193 // FWJ default path look-ahead
194 pathLookahead = 5;
195
196 newSceneGraph = NULL;
197 zcoordSetFlag = false;
198
199 //////////////////////////SUPERIMPOSED SCENE//////////////////////////
200 searcher = NULL;
201 // Used in animation; progressively scaled for gradual speed change
202 maxSpeed = 0.0f;
203
204 static const char * superimposed[] = {
205 "#Inventor V2.1 ascii", "",
206 "Separator ",
207 "{",
208 " MaterialBinding ",
209 " {",
210 " value OVERALL",
211 " }",
212 " OrthographicCamera ",
213 " {",
214 " height 1",
215 " nearDistance 0",
216 " farDistance 1",
217 " }",
218 " DEF soxt->callback Callback { }",
219 " Separator ",
220 " {",
221 " DEF soxt->translation Translation ",
222 " {",
223 " translation 0 0 0",
224 " }",
225 " DEF soxt->scale Scale ",
226 " {",
227 " scaleFactor 1 1 1",
228 " }",
229 " DEF soxt->geometry Coordinate3 ",
230 " {",
231 " point ",
232 " [",
233 " -0.81 -0.04 0, -0.81 0 0,",
234 " -0.81 0.04 0, 0 -0.04 0,",
235 " 0 0 0, 0 0.04 0,",
236 " 0.81 -0.04 0, 0.81 0 0,",
237 " 0.81 0.04 0,",
238 " 0 0.02 0,", // idx 9
239 " 0.81 0.02 0, 0.81 -0.02 0,",
240 " 0 -0.02 0,",
241 " 0 0.01 0,", // idx 13
242 " 0.4 0.01 0, 0.4 -0.01 0,",
243 " 0 -0.01 0",
244 " ]",
245 " }",
246 // current speed indicator (outline)
247 " DEF soxt->animSpeedOutlineSwitch Switch ",
248 " {",
249 " whichChild -3",
250 " Material ",
251 " {",
252 " emissiveColor 0 0 0",
253 " }",
254 " IndexedFaceSet ",
255 " {",
256 " coordIndex ",
257 " [",
258 " 12, 11, 10, 9, -1",
259 " ]",
260 " }",
261 " }",
262 // the coordinate system
263 " DEF soxt->axisSwitch Switch ",
264 " {",
265 " whichChild -3",
266 " BaseColor ",
267 " {",
268 " rgb 1 1 1",
269 " }",
270 " IndexedLineSet ",
271 " {",
272 " coordIndex ",
273 " [",
274 " 0, 2, -1,",
275 " 3, 5, -1,",
276 " 6, 8, -1,",
277 " 1, 7, -1",
278 " ]",
279 " }",
280 " }",
281 // current speed indicator
282 " DEF soxt->animSpeedSwitch Switch ",
283 " {",
284 " whichChild -3",
285 " Material ",
286 " {",
287 " emissiveColor 0 1 0",
288 " }",
289 " IndexedFaceSet ",
290 " {",
291 " coordIndex ",
292 " [",
293 " 16, 15, 14, 13, -1",
294 " ]",
295 " }",
296 " }",
297 " }",
298 // For displaying either z position (during animation) or current viewpoint name
299 " DEF soxt->curInfoSwitch Switch ",
300 " {",
301 " whichChild -3",
302 " DEF soxt->curInfoTrans Translation ",
303 " {",
304 " translation 0 0 0 ",
305 // " translation 10 20 30 ",
306 " }",
307 " DEF soxt->curInfoFont Font ",
308 " {",
309 " name defaultFont:Bold",
310 " size 16",
311 " }",
312 " DEF soxt->curInfoText Text2 ",
313 " {",
314 " string Hello",
315 " }",
316 " }",
317 // Need to use different fields for mouseover
318 // because newlines are ignored when the scene is rendered
319 " Separator ",
320 " {",
321 " DEF soxt->mouseOverTransLogName Translation ",
322 " {",
323 " translation 0 0 0 ",
324 " }",
325 " DEF soxt->mouseOverFontLogName Font ",
326 " {",
327 " name defaultFont:Bold",
328 " size 16",
329 " }",
330 " DEF soxt->mouseOverTextLogName Text2 { } ",
331 " }",
332 " Separator ",
333 " {",
334 " DEF soxt->mouseOverTransSolid Translation ",
335 " {",
336 " translation 0 0 0 ",
337 " }",
338 " DEF soxt->mouseOverFontSolid Font ",
339 " {",
340 " name defaultFont:Bold",
341 " size 16",
342 " }",
343 " DEF soxt->mouseOverTextSolid Text2 { } ",
344 " }",
345 " Separator ",
346 " {",
347 " DEF soxt->mouseOverTransMaterial Translation ",
348 " {",
349 " translation 0 0 0 ",
350 " }",
351 " DEF soxt->mouseOverFontMaterial Font ",
352 " {",
353 " name defaultFont:Bold",
354 " size 16",
355 " }",
356 " DEF soxt->mouseOverTextMaterial Text2 { } ",
357 " }",
358 " Separator ",
359 " {",
360 " DEF soxt->mouseOverTransZPos Translation ",
361 " {",
362 " translation 0 0 0 ",
363 " }",
364 " DEF soxt->mouseOverFontZPos Font ",
365 " {",
366 " name defaultFont:Bold",
367 " size 16",
368 " }",
369 " DEF soxt->mouseOverTextZPos Text2 { } ",
370 " }",
371 "}", NULL
372 };
373
374 int i, bufsize;
375 for (i = bufsize = 0; superimposed[i]; i++)
376 bufsize += strlen(superimposed[i]) + 1;
377 char * buf = new char[bufsize + 1];
378 for (i = bufsize = 0; superimposed[i]; i++) {
379 strcpy(buf + bufsize, superimposed[i]);
380 bufsize += strlen(superimposed[i]);
381 buf[bufsize] = '\n';
382 bufsize++;
383 }
384 SoInput * input = new SoInput;
385 input->setBuffer(buf, bufsize);
386 SbBool ok = SoDB::read(input, superimposition);
387 (void)ok; // FWJ added to avoid compiler warning
388 assert(ok);
389 delete input;
390 delete[] buf;
391 superimposition->ref();
392
393 sscale = (SoScale *) getSuperimpositionNode(superimposition, "soxt->scale");
394 stranslation = (SoTranslation *) getSuperimpositionNode(superimposition, "soxt->translation");
395 sgeometry = (SoCoordinate3 *) getSuperimpositionNode(superimposition, "soxt->geometry");
396 axisSwitch = (SoSwitch *) getSuperimpositionNode(superimposition, "soxt->axisSwitch");
397 animSpeedOutlineSwitch = (SoSwitch *) getSuperimpositionNode(superimposition, "soxt->animSpeedOutlineSwitch");
398 animSpeedSwitch = (SoSwitch *) getSuperimpositionNode(superimposition, "soxt->animSpeedSwitch");
399 curInfoSwitch = (SoSwitch *) getSuperimpositionNode(superimposition, "soxt->curInfoSwitch");
400 curInfoTrans = (SoTranslation *) getSuperimpositionNode(superimposition, "soxt->curInfoTrans");
401 curInfoFont = (SoFont *) getSuperimpositionNode(superimposition, "soxt->curInfoFont");
402 curInfoText = (SoText2 *) getSuperimpositionNode(superimposition, "soxt->curInfoText");
403 mouseOverTransLogName = (SoTranslation*)getSuperimpositionNode(superimposition, "soxt->mouseOverTransLogName");
404 mouseOverFontLogName = (SoFont *) getSuperimpositionNode(superimposition, "soxt->mouseOverFontLogName");
405 mouseOverTextLogName = (SoText2 *) getSuperimpositionNode(superimposition, "soxt->mouseOverTextLogName");
406 mouseOverTransSolid = (SoTranslation *) getSuperimpositionNode(superimposition, "soxt->mouseOverTransSolid");
407 mouseOverFontSolid = (SoFont *) getSuperimpositionNode(superimposition, "soxt->mouseOverFontSolid");
408 mouseOverTextSolid = (SoText2 *) getSuperimpositionNode(superimposition, "soxt->mouseOverTextSolid");
409 mouseOverTransMaterial = (SoTranslation*)getSuperimpositionNode(superimposition, "soxt->mouseOverTransMaterial");
410 mouseOverFontMaterial = (SoFont *) getSuperimpositionNode(superimposition, "soxt->mouseOverFontMaterial");
411 mouseOverTextMaterial = (SoText2 *) getSuperimpositionNode(superimposition, "soxt->mouseOverTextMaterial");
412 mouseOverTransZPos = (SoTranslation *) getSuperimpositionNode(superimposition, "soxt->mouseOverTransZPos");
413 mouseOverFontZPos = (SoFont *) getSuperimpositionNode(superimposition, "soxt->mouseOverFontZPos");
414 mouseOverTextZPos = (SoText2 *) getSuperimpositionNode(superimposition, "soxt->mouseOverTextZPos");
415
416 SoCallback * cb = (SoCallback *) getSuperimpositionNode(superimposition, "soxt->callback");
417 cb->setCallback(superimpositionCB, this);
418
419 addSuperimposition(superimposition);
420 setSuperimpositionEnabled(superimposition, FALSE);
421 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
422 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
423 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
424
425 /////////////////////\SUPERIMPOSED SCENE///////////////////////////////////
426
427}
428
429
430// Adds a menu bar and menu items to the viewer.
432{
433 if (!parent)
434 SoDebugError::post("G4OpenInventorQtExaminerViewer::buildWidget",
435 "Error: Parent is null.");
436
437 // Common font for (almost) all widgets
438 font = new QFont;
439 font->setPointSize(12);
440 // This font setting does not propagate to added child widgets - Why?
441 parent->setFont(*font);
442 // This propagates everywhere but would affect UIQt!
443 // QApplication::setFont(*font);
444
445// MENU BAR
446
447 menubar = new QMenuBar(getRenderAreaWidget());
448 // FWJ DEBUG
449 // G4cout << "G4OpenInventorQtExaminerViewer: GOT A menubar=" <<
450 // menubar << G4endl;
451
452 filemenu = new QMenu("File");
453 menubar->addMenu(filemenu);
454
455 FileOpenBookmark = new QAction("Open Bookmark File", this);
456 FileOpenBookmark->setFont(*font);
457 connect(FileOpenBookmark, SIGNAL(triggered()), this,
458 SLOT(FileOpenBookmarkCB()));
459 filemenu->addAction(FileOpenBookmark);
460
461 FileNewBookmark = new QAction("New Bookmark File", this);
462 FileNewBookmark->setFont(*font);
463 connect(FileNewBookmark, SIGNAL(triggered()), this,
464 SLOT(FileNewBookmarkCB()));
465 filemenu->addAction(FileNewBookmark);
466
467 FileLoadRefPath = new QAction("Load Reference Path", this);
468 FileLoadRefPath->setFont(*font);
469 connect(FileLoadRefPath, SIGNAL(triggered()), this,
470 SLOT(FileLoadRefPathCB()));
471 filemenu->addAction(FileLoadRefPath);
472
473 FileSaveRefPath = new QAction("Save Reference Path", this);
474 FileSaveRefPath->setFont(*font);
475 connect(FileSaveRefPath, SIGNAL(triggered()), this,
476 SLOT(FileSaveRefPathCB()));
477 filemenu->addAction(FileSaveRefPath);
478
479 FileLoadSceneGraph = new QAction("Load scene graph", this);
480 FileLoadSceneGraph->setFont(*font);
481 connect(FileLoadSceneGraph, SIGNAL(triggered()), this,
482 SLOT(FileLoadSceneGraphCB()));
483 filemenu->addAction(FileLoadSceneGraph);
484
485 FileSaveSceneGraph = new QAction("Save scene graph", this);
486 FileSaveSceneGraph->setFont(*font);
487 connect(FileSaveSceneGraph, SIGNAL(triggered()), this,
488 SLOT(FileSaveSceneGraphCB()));
489 filemenu->addAction(FileSaveSceneGraph);
490
491 // Rest of File menu is done in G4OpenInventorQtViewer
492
493 toolsmenu = new QMenu("Tools");
494 menubar->addMenu(toolsmenu);
495
496 ToolsAnimateRefParticle = new QAction("Fly on Ref Path", this);
497 ToolsAnimateRefParticle->setFont(*font);
498 connect(ToolsAnimateRefParticle, SIGNAL(triggered()), this,
499 SLOT(ToolsAnimateRefParticleCB()));
500 toolsmenu->addAction(ToolsAnimateRefParticle);
501
502 ToolsRefPathStart = new QAction("Go to start of Ref Path", this);
503 ToolsRefPathStart->setFont(*font);
504 connect(ToolsRefPathStart, SIGNAL(triggered()), this,
505 SLOT(ToolsRefPathStartCB()));
506 toolsmenu->addAction(ToolsRefPathStart);
507
508 ToolsRefPathInvert = new QAction("Invert Ref Path", this);
509 ToolsRefPathInvert->setFont(*font);
510 connect(ToolsRefPathInvert, SIGNAL(triggered()), this,
511 SLOT(ToolsRefPathInvertCB()));
512 toolsmenu->addAction(ToolsRefPathInvert);
513
514 etcmenu = new QMenu("Etc");
515 menubar->addMenu(etcmenu);
516
517 // All Etc menu items are done in G4OpenInventorQtViewer
518
519 helpmenu = new QMenu("Help");
520 menubar->addMenu(helpmenu);
521
522 HelpControls = new QAction("Controls", this);
523 HelpControls->setFont(*font);
524 connect(HelpControls, SIGNAL(triggered()), this, SLOT(HelpControlsCB()));
525 helpmenu->addAction(HelpControls);
526
527 menubar->show();
528
529 // SoQtExaminerViewer::buildWidget(parent);
530
531 // APP VIEWER BUTTONS have their own box on upper left
532 // The built in viewer button list is PRIVATE
533
534 saveViewPtButton = new QPushButton;
535 saveViewPtButton->setIcon(QPixmap((const char **)saveViewPt_xpm));
536 saveViewPtButton->setIconSize(QSize(24,24));
537 saveViewPtButton->setToolTip("Bookmark this view");
538 connect(saveViewPtButton, SIGNAL(clicked()), this,
539 SLOT(SaveViewPtCB()));
540 addAppPushButton(saveViewPtButton);
541
542 nextViewPtButton = new QPushButton;
543 nextViewPtButton->setIconSize(QSize(24,24));
544 QCommonStyle style;
545 nextViewPtButton->setIcon(style.standardIcon(QStyle::SP_ArrowRight));
546 nextViewPtButton->setToolTip("Next bookmark");
547 connect(nextViewPtButton, SIGNAL(clicked()), this,
548 SLOT(NextViewPtCB()));
549 addAppPushButton(nextViewPtButton);
550
551 prevViewPtButton = new QPushButton;
552 prevViewPtButton->setIconSize(QSize(24,24));
553 prevViewPtButton->setIcon(style.standardIcon(QStyle::SP_ArrowLeft));
554 prevViewPtButton->setToolTip("Previous bookmark");
555 connect(prevViewPtButton, SIGNAL(clicked()), this,
556 SLOT(PrevViewPtCB()));
557 addAppPushButton(prevViewPtButton);
558
559 abbrOutputButton = new QPushButton;
560 abbrOutputButton->setCheckable(true);
561 abbrOutputButton->setIconSize(QSize(24,24));
562 abbrOutputButton->setIcon(QPixmap((const char **)pickext_xpm));
563 abbrOutputButton->setToolTip("Extended picking & readout");
564 connect(abbrOutputButton, SIGNAL(toggled(bool)), this,
565 SLOT(AbbrOutputCB(bool)));
566 addAppPushButton(abbrOutputButton);
567
568 pickRefPathButton = new QPushButton;
569 pickRefPathButton->setIconSize(QSize(24,24));
570 pickRefPathButton->setIcon(QPixmap((const char **)pickref_xpm));
571 pickRefPathButton->setToolTip("Pick ref trajectory");
572 connect(pickRefPathButton, SIGNAL(clicked()), this,
573 SLOT(PickRefPathCB()));
574 addAppPushButton(pickRefPathButton);
575
576 switchWireFrameButton = new QPushButton;
577 switchWireFrameButton->setCheckable(true);
578 switchWireFrameButton->setIconSize(QSize(24,24));
579 switchWireFrameButton->setIcon(QPixmap((const char **)wireframe_xpm));
580 switchWireFrameButton->setToolTip("Switch wireframe/solid");
581 connect(switchWireFrameButton, SIGNAL(toggled(bool)), this,
582 SLOT(SwitchWireFrameCB(bool)));
583 addAppPushButton(switchWireFrameButton);
584
585 switchAxesButton = new QPushButton;
586 switchAxesButton->setCheckable(true);
587 switchAxesButton->setText(QString("A"));
588 switchAxesButton->setToolTip("Axes on/off");
589 connect(switchAxesButton, SIGNAL(toggled(bool)), this,
590 SLOT(SwitchAxesCB(bool)));
591 addAppPushButton(switchAxesButton);
592
593 detachButton = new QPushButton;
594 detachButton->setIconSize(QSize(24,24));
595 detachButton->setIcon(style.standardIcon(QStyle::SP_CommandLink));
596 detachButton->setToolTip("Detach viewer window");
597 connect(detachButton, SIGNAL(clicked()), this,
598 SLOT(DetachCB()));
599 // Used for UIQt only so check and add later
600 // addAppPushButton(detachButton);
601
602 // HELP WINDOW
603
604 helpmsgbox = new QMessageBox(getParentWidget());
605 helpmsgbox->setWindowTitle("OIQt Controls");
606 helpmsgbox->setFont(*font);
607 QString messagetxt =
608"\nVIEWING mode (Hand cursor):\n\n\
609 Left-button + pointer move: rotate\n\
610 Shift+Left-button + pointer move: pan\n\
611 Middle-button + pointer move: pan\n\
612 Ctrl+Shift+Left-button + pointer move: zoom\n\
613 Mouse wheel: zoom\n\
614 Right-button: popup menu\n\n\
615PICKING mode (Arrow cursor):\n\n\
616 Click on a volume: geometry readout\n\
617 Click on a trajectory: particle & trajectory readout\n\
618 Ctrl + click on a volume: see daughters.\n\
619 Shift + click on a volume: see mother.\n\n\
620EXTENDED PICKING mode (Arrow+ viewer button):\n\n\
621 Hover the mouse over a volume or trajectory for\n\
622 overlayed readout.\n\n\
623ELEMENT NAVIGATION (requires Reference Path):\n\n\
624 Click on element in list: centers view on element\n\
625 Arrow keys: rotate in 90 degree steps around element \n\
626 Shift + Right Arrow: move to next element\n\
627 Shift + Left Arrow: move to previous element\n\n\
628FLY mode (requires Reference Path):\n\n\
629 Page Up: Increase speed\n\
630 Page Down: Decrease speed (& reverse if wanted)\n\
631 Up Arrow: raise camera above path\n\
632 Down Arror: lower camera below path\n\
633 Escape: Exit fly mode";
634 helpmsgbox->setText(messagetxt);
635 helpmsgbox->setModal(false);
636 // helpmsgbox->setWindowModality(Qt::NonModal);
637
638 // AUXILIARY LISTS WINDOW
639
640 // Bypass the namespace in order to make a persistent object
641 AuxWindowDialog = new Ui_Dialog;
642 AuxWindow = new QDialog(parent);
643 AuxWindowDialog->setupUi(AuxWindow);
644
645 // SIGNALS
646 connect(AuxWindowDialog->listWidget, SIGNAL(itemClicked(QListWidgetItem*)),
647 this, SLOT(LoadBookmarkCB(QListWidgetItem*)));
648 connect(AuxWindowDialog->listWidget1, SIGNAL(itemClicked(QListWidgetItem*)),
649 this, SLOT(LookAtSceneElementCB(QListWidgetItem*)));
650 connect(AuxWindowDialog->pushButton_2, SIGNAL(clicked()),
651 this, SLOT(DeleteBookmarkCB()));
652 connect(AuxWindowDialog->pushButton_3, SIGNAL(clicked()),
653 this, SLOT(RenameBookmarkCB()));
654 connect(AuxWindowDialog->pushButton, SIGNAL(clicked()),
655 this, SLOT(SortBookmarksCB()));
656
657 // FWJ Better to do this after viewer window is realized
658 // AuxWindow->show();
659 // AuxWindow->raise();
660 // AuxWindow->activateWindow();
661}
662
663
664// Called right after buttons and widgets get realized.
665// It sets the viewpoint last accessed.
667{
668 SoQtExaminerViewer::afterRealizeHook();
669
670 // Default height is used when selecting and viewing scene elements
671 // FWJ Added defaultHeight for Ortho camera
672 SoCamera *cam = getCamera();
673 if (cam) {
674 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
676 ((SoPerspectiveCamera *) cam)->heightAngle.getValue();
677 toggleCameraType();
679 ((SoOrthographicCamera *) cam)->height.getValue();
680 toggleCameraType();
681 } else {
683 ((SoOrthographicCamera *) cam)->height.getValue();
684 toggleCameraType();
685 cam = getCamera();
686 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
688 ((SoPerspectiveCamera *) cam)->heightAngle.getValue();
689 toggleCameraType();
690 }
691 }
692
693 // Open the default bookmark file
694 fileIn.open(fileName.c_str());
695 if (!fileIn.fail()) {
696 if (!loadViewPts()) {
697 QMessageBox msgbox;
698 msgbox.setFont(*font);
699 QString messagetxt = "Error reading bookmark file ";
700 messagetxt.append(QString(fileName.c_str()));
701 msgbox.setText(messagetxt);
702 msgbox.exec();
703 } else {
704 // Opens a file without erasing it
705 fileOut.open(fileName.c_str(), std::ios::in);
706 fileOut.seekp(0, std::ios::end); // For appending new data to the end
707 // FWJ DEBUG
708 // G4cout << "afterRealizeHook: opened EXISTING bookmark file"
709 // << G4endl;
710 if (viewPtList.size()) {
711 // FWJ disabled auto-selection of first viewpoint.
712 // Initial view should be user-controllable & not forced
713 // setViewPt();
715 }
716 }
717 fileIn.close();
718 } else {
719 // Creates a new default bookmark file
720 fileOut.open(fileName.c_str());
721 // FWJ DEBUG
722 // G4cout << "afterRealizeHook: Opened a NEW bookmark file" << G4endl;
723 }
724
725 fileIn.clear();
726
727 SoSeparator* root = (SoSeparator*) (getSceneManager()->getSceneGraph());
728 if (root == NULL)
729 SoDebugError::post("G4OpenInventorQtExaminerViewer::afterRealizeHook", "Root is null.");
730 else {
731 root->addChild(myCam); // For position/orientation calculation during animation
732 }
733
734 sceneChangeSensor = new SoNodeSensor;
735 sceneChangeSensor->setFunction(sceneChangeCB);
736 sceneChangeSensor->attach(root);
737 sceneChangeSensor->setData(this);
738
739 ///////////////////////////// MOUSEOVER & PICK /////////////////////
740
741 // Monitor mouseover events for displaying the name of scene elements
742 // An SoEventCallback is needed instead of using the default processSoEvent
743 // because that last one does not provide us with an SoPath to the object
744 // that was picked
745 SoEventCallback *moCB = new SoEventCallback;
746 moCB->addEventCallback(
747 SoLocation2Event::getClassTypeId(),
748 mouseoverCB, static_cast<void *>(this));
749 root->addChild(moCB);
750
751 // Override the default picking mechanism present in G4OpenInventorViewer
752 // because we want abbreviated output when picking a trajectory
753 SoEventCallback *pickCB = new SoEventCallback;
754 pickCB->addEventCallback(
755 SoMouseButtonEvent::getClassTypeId(),
756 pickingCB, static_cast<void *>(this));
757 root->addChild(pickCB);
758
759 ///////////////////////////// MOUSEOVER & PICK /////////////////////
760
761 AuxWindow->show();
762 AuxWindow->raise();
763 AuxWindow->activateWindow();
764
765 auto UI = G4UImanager::GetUIpointer();
766 uiQt = dynamic_cast<G4UIQt*>(UI->GetG4UIWindow());
767 // This explicitly sets the TabWidget as parent before addTab():
768 if (uiQt) {
769 viewerParent = getParentWidget();
770 viewerParent2 = viewerParent->parentWidget();
771 uiQt->AddTabWidget(getParentWidget(), *fName);
772 uiQtTabIndex = uiQt->GetViewerTabWidget()->currentIndex();
773 // attached = TRUE;
774 addAppPushButton(detachButton);
775 }
776}
777
778
779// This method locates a named node in the superimposed or original scene.
780// FWJ RENAME THIS
781SoNode*
783 const char* name)
784{
785 if (!searcher)
786 searcher = new SoSearchAction;
787 searcher->reset();
788 searcher->setName(SbName(name));
789 searcher->setInterest(SoSearchAction::FIRST);
790 searcher->setSearchingAll(TRUE);
791 searcher->apply(root);
792 assert(searcher->getPath());
793 return searcher->getPath()->getTail();
794}
795
796
797// FWJ don't know why userdata is called "closure"
798// It contains the this pointer!
800 SoAction * action)
801{
802 if (closure)
803 ((G4OpenInventorQtExaminerViewer*)closure)->superimpositionEvent(action);
804}
805
806
807// Renders and positions speed indicator and longitudinal
808// distance/viewpoint name on the drawing canvas
810{
811
812 if (!action->isOfType(SoGLRenderAction::getClassTypeId()))
813 return;
814 SbViewportRegion vpRegion =
815 ((SoGLRenderAction*)action)->getViewportRegion();
816 SbVec2s viewportSize = vpRegion.getViewportSizePixels();
817
818 // Aspect is WIDTH/HEIGHT
819 float aspect = float(viewportSize[0]) / float(viewportSize[1]);
820
821 // FWJ DEBUG
822 // G4cout << "SPEVENT X0 Y0 DX DY aspect: " << vpRegion.getViewportOrigin()[0] <<
823 // " " << vpRegion.getViewportOrigin()[1] <<
824 // " " << viewportSize[0] <<
825 // " " << viewportSize()[1] <<
826 // " " << aspect << G4endl;
827
828 // Translation and scale factor for animation speed indicator...
829
830 float factorx = 1.0f / float(viewportSize[1]) * 220.0f;
831 float factory = factorx;
832
833 if (aspect > 1.0f) {
834 stranslation->translation.setValue(SbVec3f(0.0f, -0.4f, 0.0f));
835 } else {
836 stranslation->translation.setValue(SbVec3f(0.0f, -0.4f / aspect, 0.0f));
837 factorx /= aspect;
838 factory /= aspect;
839 }
840 if (viewportSize[0] > 500)
841 factorx *= 500.0f / 400.0f;
842 else
843 factorx *= float(viewportSize[0]) / 400.0f;
844
845 sscale->scaleFactor.setValue(SbVec3f(factorx, factory, 1.0f));
846
847 // TEXT OVERLAY...
848
849 // FWJ Simplified and rewrote the following section to ease problems
850 // with the overlayed text after a viewer window resize.
851 // Result is now readable but needs further refinement of the scaling.
852
853 float xInfo, yInfo, xLogName, yLogName, xSolid, ySolid,
854 xMaterial, yMaterial, xZPos, yZPos;
855
856 // Base point for navigation distance or viewpoint name
857 // Origin is at center of render area.
858 xInfo = -.475;
859 yInfo = .475;
860 // Menu bar height in same coordinates:
861 float mbgap = 0.03;
862 if (aspect > 1.) xInfo = xInfo*aspect;
863 if (aspect < 1.) yInfo = yInfo/aspect;
864 yInfo = yInfo - mbgap*aspect;
865
866 // Following are relative to above base point
867 xLogName = 0.0;
868 yLogName = -.88 + mbgap*aspect;
869 xSolid = 0.0;
870 ySolid = -.91 + mbgap*aspect;
871 xMaterial = 0.0;
872 yMaterial = -.94 + mbgap*aspect;
873 xZPos = 0.0;
874 yZPos = -.97 + mbgap*aspect;
875
876 // Top line
877 curInfoTrans->translation.setValue(SbVec3f(xInfo, yInfo, 0.0));
878
879 // Bottom lines
880 mouseOverTransLogName->translation.setValue(SbVec3f(xLogName, yLogName, 0.0));
881 mouseOverTransSolid->translation.setValue(SbVec3f(xSolid, ySolid, 0.0));
882 mouseOverTransMaterial->translation.setValue(SbVec3f(xMaterial, yMaterial, 0.0));
883 mouseOverTransZPos->translation.setValue(SbVec3f(xZPos, yZPos, 0.0));
884
885 if (currentState == VIEWPOINT) { // Displaying viewpoint name
886 curInfoFont->size.setValue(15);
887 curInfoFont->name.setValue("defaultFont:Italic");
888 curInfoText->string.setValue(SbString(curViewPtName));
889 }
890 else if(currentState == GENERAL) { // Displaying longitudinal distance
891 curInfoFont->size.setValue(16);
892 curInfoFont->name.setValue("defaultFont:Bold");
893 curInfoText->string.setValue(SbString(""));
894 }
895 else {
896 if (refParticleIdx < (int) refParticleTrajectory.size() - 1) {
897 curInfoFont->size.setValue(16);
898 curInfoFont->name.setValue("defaultFont:Bold");
899 char zPos[20];
900 // FWJ need a better format here
901 snprintf(zPos, sizeof zPos, "%-7.2f [m]", refZPositions[refParticleIdx] / 1000);
902 curInfoText->string.setValue(SbString(zPos));
903 }
904 }
905}
906
907
908// Loads view point data from a file into a vector.
909
911{
912 bool error = false;
913 viewPtData tmp;
914 std::string token;
915 SbVec3f axis;
916 SbRotation orient;
917 float x{0.0}, y{0.0}, z{0.0}, angle{0.0};
918
919 // Gets the last view point accessed, stored in the first line of the data file.
920 fileIn >> token;
921 parseString<int>(viewPtIdx, token, error);
922 getline(fileIn, token); // Remove "\n"
923 // Converts data from string type into necessary types
924 while (getline(fileIn, token)) {
925
926 std::size_t end = token.find_last_not_of(' '); // Remove padded spaces
927 token = token.substr(0, end + 1);
928
929 char *vpName = new char[token.size() + 1];
930 strcpy(vpName, token.c_str());
931 tmp.viewPtName = vpName;
932 fileIn >> token;
933
934 parseString<float>(x, token, error);
935 fileIn >> token;
936 parseString<float>(y, token, error);
937 fileIn >> token;
938 parseString<float>(z, token, error);
939 fileIn >> token;
940 tmp.position = axis.setValue(x, y, z);
941
942 parseString<float>(x, token, error);
943 fileIn >> token;
944 parseString<float>(y, token, error);
945 fileIn >> token;
946 parseString<float>(z, token, error);
947 fileIn >> token;
948 parseString<float>(angle, token, error);
949 fileIn >> token;
950 orient.setValue(axis.setValue(x, y, z), angle);
951 tmp.orientation = orient.getValue();
952
953 int camType{0};
954 parseString<int>(camType, token, error);
955 fileIn >> token;
956 tmp.camType = (CameraType) camType;
957
958 parseString<float>(tmp.height, token, error);
959 fileIn >> token;
960 parseString<float>(tmp.focalDistance, token, error);
961 fileIn >> token;
962 parseString<float>(tmp.nearDistance, token, error);
963 fileIn >> token;
964 parseString<float>(tmp.farDistance, token, error);
965 fileIn >> token;
966 parseString<int>(tmp.viewportMapping, token, error);
967 fileIn >> token;
968 parseString<float>(tmp.aspectRatio, token, error);
969
970 getline(fileIn, token); // To remove "\n" characters
971 getline(fileIn, token);
972
973 if (error) {
974 viewPtIdx = 0;
975 viewPtList.clear();
976 return false;
977 }
978 viewPtList.push_back(tmp);
979 }
980
981 return true;
982}
983
984
985// Rotates camera 90 degrees around a scene element.
986// Rotation is animated for smoothness.
988{
989 SoCamera *cam = getCamera();
990
991 SbRotation rot(rotAxis, M_PI / (2 * ROT_CNT));
992 rot.multVec(camDir, camDir);
993 rot.multVec(camUpVec, camUpVec);
994
995 SbVec3f camPosNew = prevPt - (camDir*distance);
996 cam->position = camPosNew;
997 cam->pointAt(prevPt, camUpVec);
998 cam->focalDistance = (prevPt - camPosNew).length();
999
1000 rotCnt--;
1001
1002 if (animateSensorRotation->isScheduled()) {
1003 animateSensorRotation->unschedule();
1004 }
1005
1006 animateSensorRotation->setBaseTime(SbTime::getTimeOfDay());
1007 animateSensorRotation->setInterval(SbTime(0.02));
1008 animateSensorRotation->schedule();
1009
1010}
1011
1012
1013// Slides camera along the beamline.
1014void G4OpenInventorQtExaminerViewer::moveCamera(float dist, bool lookdown)
1015{
1016
1017 SoCamera *cam = getCamera();
1018 SbVec3f p1, p2; // The particle moves from p1 to p2
1019 SbVec3f particleDir; // Direction vector from p1 to p2
1020 SbVec3f camPosNew{0.0f, 0.0f, 0.0f}; // New position of the camera
1021
1022 if(refParticleTrajectory.size() == 0) {
1023 //refParticleTrajectory hasn't been set yet
1024 if(dist)
1025 distance = dist;
1026 else
1027 distance = (cam->position.getValue() - center).length();
1028
1029 cam->position.setValue(center + offsetFromCenter*distance);
1030 cam->focalDistance = (cam->position.getValue() - center).length();
1031 cam->pointAt(center, upVector);
1032 }
1033 else {
1034
1035 // If we move forward past the last trajectory point,
1036 // go back to the beginning
1037 if (refParticleIdx >= (int) refParticleTrajectory.size() - 1) {
1039 dist = (prevPt - cam->position.getValue()).length();
1040 refParticleIdx = 0;
1041 }
1042 // If we move backward past the beginning,
1043 // go to the last trajectory point
1044 if (refParticleIdx < 0) {
1046 dist = (prevPt - cam->position.getValue()).length();
1047 refParticleIdx = (int) refParticleTrajectory.size() - 2;
1048 }
1049
1050 // Set start and end points
1053
1054 // Get the direction from p1 to p2
1055 particleDir = p2 - p1;
1056 particleDir.normalize();
1057
1058 if(prevParticleDir == SbVec3f(0,0,0)) {
1059 // First time entering BEAMLINE mode, look at
1060 // the element from the front, with camera upright
1061 if(lookdown)
1062 camDir = SbVec3f(0,0,1);
1063 else
1064 camDir = SbVec3f(1,0,0);
1065 camUpVec = SbVec3f(0,1,0);
1066
1067 // In case the start of the goes in a
1068 // direction other than +z, rotate the camera accordingly
1069 SbRotation rot(SbVec3f(0,0,1), particleDir);
1070 rot.multVec(camDir, camDir);
1071 rot.multVec(camUpVec, camUpVec);
1072
1073 }
1074 else if(particleDir != prevParticleDir) {
1075 // The beamline has changed direction
1076
1077 SbRotation rot(prevParticleDir, particleDir);
1078 rot.multVec(camDir, camDir);
1079 rot.multVec(camUpVec, camUpVec);
1080
1081 }
1082
1083 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
1084 if (!dist)
1085 distance = (prevPt - cam->position.getValue()).length();
1086 else
1087 distance = dist;
1088 }
1089
1090 // FWJ distance not relevant -- use focalDistance
1091 // if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
1092 // if (!dist)
1093 // distance = (prevPt - cam->position.getValue()).length();
1094 // else
1095 // distance = dist;
1096 // }
1097
1098
1099 float x,y,z;
1100 prevPt.getValue(x,y,z);
1101
1102
1103 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
1104 camPosNew = p2 - (camDir*distance);
1105 }
1106 if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
1107 // FWJ maintain focal distance
1108 camPosNew = p2 - (camDir*cam->focalDistance.getValue());
1109 // camPosNew = p2 - (camDir);
1110 }
1111
1112 cam->position = camPosNew;
1113 cam->pointAt(p2, camUpVec);
1114 cam->focalDistance = (p2 - camPosNew).length();
1115
1116 p2.getValue(x,y,z);
1117 camPosNew.getValue(x,y,z);
1118
1119 prevParticleDir = particleDir;
1120 prevPt = p1; // For accurate distance calculation
1121
1122 }
1123
1124}
1125
1126
1128 SoEventCallback *eventCB)
1129{
1130 SoHandleEventAction* action = eventCB->getAction();
1131 const SoPickedPoint *pp = action->getPickedPoint();
1133
1134 if(pp != NULL) {
1135
1136 SoPath* path = pp->getPath();
1137 SoNode* node = ((SoFullPath*)path)->getTail();
1138
1139 if(node->getTypeId() == SoLineSet::getClassTypeId()) {
1140
1141 if(This->pickRefPathFlag) {
1142 This->pickRefPathFlag = false;
1143 if(This->viewingBeforePickRef != This->isViewing())
1144 This->setViewing(This->viewingBeforePickRef);
1145 else
1146 This->setComponentCursor(SoQtCursor(SoQtCursor::DEFAULT));
1147
1148 // The trajectory is a set of lines stored in a LineSet
1149 SoLineSet * trajectory = (SoLineSet *)node;
1150 // FWJ DEBUG
1151 // G4cout << "FOUND trajectory LineSet" << trajectory << G4endl;
1152
1153 // The set of all trajectories is stored in a Seperator group node
1154 // one level above the LineSet that was picked. The nodes under that
1155 // seperator are as follows (in this order): Material, LightModel,
1156 // ResetTransform, MatrixTransform, Coordinate3, DrawStyle, LineSet
1157 SoSeparator * grpNode =
1158 (SoSeparator*)(((SoFullPath*)path)->getNodeFromTail(1));
1159
1160 // The node that contains the coordinates for the trajectory is a
1161 // Coordinate3 node which occurs before the LineSet node. We iterate
1162 // back through the nodes in the group until we find the Coordinate3 node
1163 int nodeIndex = grpNode->findChild(trajectory);
1164 SoNode * tmpNode;
1165 // FWJ needs initialization
1166 SoCoordinate3 * coords = 0;
1167 // SoCoordinate3 * coords;
1168 // We allow only 100 iterations, in case the node isn't found
1169 // (should take only a few iterations)
1170 for(int i = 0; i < 100; ++i) {
1171 --nodeIndex;
1172
1173 tmpNode = grpNode->getChild(nodeIndex);
1174 if(tmpNode->getTypeId() == SoCoordinate3::getClassTypeId()) {
1175 //node found
1176 coords = (SoCoordinate3 *)tmpNode;
1177 break;
1178 }
1179 }
1180
1181 if(coords == NULL) {
1182 G4warn << "Could not find the coordinates node"
1183 " for the picked trajectory." << G4endl;
1184 G4warn << " Reference trajectory not set" << G4endl;
1185 return;
1186 }
1187 // FWJ DEBUG
1188 // G4cout << "FOUND SoCoordinate3 node " << coords << G4endl;
1189
1190
1191 if ((This->lshiftdown) || (This->rshiftdown))
1192 This->setReferencePath(trajectory, coords, true); //APPENDING
1193 else
1194 This->setReferencePath(trajectory, coords, false);
1195
1196 return;
1197
1198 }
1199 else if(This->abbrOutputFlag) {
1200
1201 G4AttHolder* attHolder = dynamic_cast<G4AttHolder*>(node);
1202 if(attHolder && attHolder->GetAttDefs().size()) {
1203
1204 std::string strTrajPoint = "G4TrajectoryPoint:";
1205 std::ostringstream oss;
1206 for (std::size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
1207 G4cout << G4AttCheck(attHolder->GetAttValues()[i],
1208 attHolder->GetAttDefs()[i]);
1209 oss << G4AttCheck(attHolder->GetAttValues()[i],
1210 attHolder->GetAttDefs()[i]);
1211 if(oss.str().find(strTrajPoint) != std::string::npos) {
1212
1213 // Last attribute displayed was a trajectory point. Since we
1214 // want abbreviated output, display the last one and exit
1215 // (unless we're already at the last (and only) trajectory point)
1216 if(i != attHolder->GetAttDefs().size()-1) {
1217 G4cout << G4AttCheck(
1218 attHolder->GetAttValues()[attHolder->GetAttDefs().size()-1],
1219 attHolder->GetAttDefs()[attHolder->GetAttDefs().size()-1]);
1220 }
1221 break;
1222 }
1223 }
1224 } else {
1225 G4String name((char*)node->getName().getString());
1226 G4String cls((char*)node->getTypeId().getName().getString());
1227 G4warn << "SoNode : " << node
1228 << " SoType : " << cls
1229 << " name : " << name
1230 << G4endl;
1231 G4warn << "No attributes attached." << G4endl;
1232 }
1233
1234 return;
1235 }
1236 else{
1237 //Go to default behavior
1238 }
1239 }
1240 else {
1241 //Go to default behavior
1242 }
1243
1244 // Default behavior in G4OpenInventorViewer::SelectionCB
1245 G4AttHolder* attHolder = dynamic_cast<G4AttHolder*>(node);
1246 if(attHolder && attHolder->GetAttDefs().size()) {
1247 for (std::size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
1248 G4cout << G4AttCheck(attHolder->GetAttValues()[i],
1249 attHolder->GetAttDefs()[i]);
1250 }
1251 } else {
1252 G4String name((char*)node->getName().getString());
1253 G4String cls((char*)node->getTypeId().getName().getString());
1254 G4warn << "SoNode : " << node
1255 << " SoType : " << cls
1256 << " name : " << name
1257 << G4endl;
1258 G4warn << "No attributes attached." << G4endl;
1259 }
1260
1261 //Suppress other event handlers
1262 eventCB->setHandled();
1263 }
1264}
1265
1266
1267void G4OpenInventorQtExaminerViewer::mouseoverCB(void *aThis, SoEventCallback *eventCB)
1268{
1269 SoHandleEventAction* action = eventCB->getAction();
1270 const SoPickedPoint* pp = action->getPickedPoint();
1272
1273 if(!This->abbrOutputFlag)
1274 return;
1275
1276 if(pp != NULL) {
1277
1278 const SbViewportRegion & viewportRegion = action->getViewportRegion();
1279
1280 std::string sLogName;
1281 float x,y,z;
1282 std::stringstream ssZPos;
1283 std::stringstream ssSolids;
1284 std::stringstream ssMaterials;
1285 SoPath * path = pp->getPath();
1286 SoNode* node = ((SoFullPath*)path)->getTail();
1287
1288 if(node->getTypeId() == Geant4_SoPolyhedron::getClassTypeId()) {
1289
1290 sLogName = "Logical Volume: ";
1291 sLogName += ((Geant4_SoPolyhedron *)node)->getName().getString();
1292
1293 SoGetBoundingBoxAction bAction(viewportRegion);
1294 bAction.apply((SoFullPath*)path);
1295 SbBox3f bBox = bAction.getBoundingBox();
1296 SbVec3f centr = bBox.getCenter();
1297 centr.getValue(x,y,z);
1298 ssZPos << "Pos: " << x << " " << y << " " << z;
1299
1300 G4AttHolder* attHolder = dynamic_cast<G4AttHolder*>(node);
1301 if(attHolder && attHolder->GetAttDefs().size()) {
1302
1303 std::vector<const std::map<G4String,G4AttDef>*> vecDefs =
1304 attHolder->GetAttDefs();
1305 std::vector<const std::vector<G4AttValue>*> vecVals =
1306 attHolder->GetAttValues();
1307 for (std::size_t i = 0; i < vecDefs.size(); ++i) {
1308 const std::vector<G4AttValue> * vals = vecVals[i];
1309
1310 std::vector<G4AttValue>::const_iterator iValue;
1311
1312 for (iValue = vals->begin(); iValue != vals->end(); ++iValue) {
1313 const G4String& valueName = iValue->GetName();
1314 const G4String& value = iValue->GetValue();
1315
1316 if(valueName == "Solid") {
1317 if(ssSolids.str() == "")
1318 ssSolids << "Solid Name: " << value;
1319 else
1320 ssSolids << ", " << value;
1321 }
1322
1323 if(valueName == "Material") {
1324 if(ssMaterials.str() == "")
1325 ssMaterials << "Material Name: " << value;
1326 else
1327 ssMaterials << ", " << value;
1328 }
1329 }
1330 }
1331 }
1332 }
1333 // FWJ Mouseover for trajectories
1334 else if(node->getTypeId() == SoLineSet::getClassTypeId()) {
1335 // G4cout << "Trajectory!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
1336 G4AttHolder* attHolder = dynamic_cast<G4AttHolder*>(node);
1337 if(attHolder && attHolder->GetAttDefs().size()) {
1338 std::string strTrajPoint = "G4TrajectoryPoint:";
1339 std::ostringstream oss;
1340 G4String t1, t1Ch, t2, t3, t4;
1341 for (std::size_t i = 0; i < attHolder->GetAttDefs().size(); ++i) {
1342 // G4cout << "Getting index " << i << " from attHolder" << G4endl;
1343 // No, returns a vector!
1344 // G4AttValue* attValue = attHolder->GetAttValues()[i];
1345 const std::vector<G4AttValue>* vals = attHolder->GetAttValues()[i];
1346 std::vector<G4AttValue>::const_iterator iValue;
1347 for (iValue = vals->begin(); iValue != vals->end(); ++iValue) {
1348 const G4String& valueName = iValue->GetName();
1349 const G4String& value = iValue->GetValue();
1350 // G4cout << " valueName = " << valueName << G4endl;
1351 // G4cout << " value = " << value << G4endl;
1352 // LINE 1
1353 if (valueName == "PN") t1 = value;
1354 if (valueName == "Ch") {
1355 if (atof(value.c_str()) > 0)
1356 t1Ch = " +";
1357 else
1358 t1Ch = " ";
1359 t1Ch += value;
1360 }
1361 if (valueName == "PDG") {
1362 t1 += " ";
1363 t1 += value;
1364 t1 += t1Ch;
1365 This->mouseOverTextLogName->string.setValue(t1);
1366 }
1367 // G4cout << " t1 = " << t1 << G4endl;
1368 // LINE 2
1369 if (valueName == "EventID") t2 = "Evt " + value;
1370 if (valueName == "ID") t2 += " Trk " + value;
1371 if (valueName == "PID") {
1372 t2 += " Prt " + value;
1373 This->mouseOverTextSolid->string.setValue(t2);
1374 }
1375 // LINE 3
1376 if (valueName == "IKE") t3 = "KE " + value;
1377 if (valueName == "IMom") {
1378 // Remove units
1379 std::size_t ipos = value.rfind(" ");
1380 G4String value1 = value;
1381 value1.erase(ipos);
1382 t3 += " P (" + value1 + ")";
1383 }
1384 if (valueName == "IMag") {
1385 t3 += " " + value + "/c";
1386 // t3 += " " + value;
1387 This->mouseOverTextMaterial->string.setValue(t3);
1388 }
1389 // LINE 4
1390 if (valueName == "NTP") {
1391 std::ostringstream t4oss;
1392 t4oss << "TrjPts " << value;
1393 t4oss << " Pos " << pp->getPoint()[0] << " " << pp->getPoint()[1] <<
1394 " " << pp->getPoint()[2];
1395 This->mouseOverTextZPos->string.setValue(SbString(t4oss.str().c_str()));
1396 }
1397 }
1398// G4cout << " NOW CALLING G4AttCheck" << G4endl;
1399// G4cout << G4AttCheck(attHolder->GetAttValues()[i],
1400// attHolder->GetAttDefs()[i]);
1401// oss << G4AttCheck(attHolder->GetAttValues()[i],
1402// attHolder->GetAttDefs()[i]);
1403// if(oss.str().find(strTrajPoint) != std::string::npos) {
1404// // Last attribute displayed was a trajectory point. Since we
1405// // want abbreviated output, display the last one and exit
1406// // (unless we're already at the last (and only) trajectory point)
1407// if(i != attHolder->GetAttDefs().size()-1) {
1408// G4cout << G4AttCheck(
1409// attHolder->GetAttValues()[attHolder->GetAttDefs().size()-1],
1410// attHolder->GetAttDefs()[attHolder->GetAttDefs().size()-1]);
1411// }
1412// break;
1413// }
1414 }
1415 }
1416 This->setSuperimpositionEnabled(This->superimposition, TRUE);
1417 This->scheduleRedraw();
1418 eventCB->setHandled();
1419 return;
1420 }
1421
1422 bool redraw = false;
1423 if(std::string(This->mouseOverTextLogName->string.getValues(0)->getString()) != sLogName) {
1424 This->mouseOverTextLogName->string.setValue(SbString(sLogName.c_str()));
1425 redraw = true;
1426 }
1427 if(std::string(This->mouseOverTextSolid->string.getValues(0)->getString()) != ssSolids.str()) {
1428 This->mouseOverTextSolid->string.setValue(SbString(ssSolids.str().c_str()));
1429 redraw = true;
1430 }
1431 if(std::string(This->mouseOverTextMaterial->string.getValues(0)->getString()) != ssMaterials.str()) {
1432 This->mouseOverTextMaterial->string.setValue(SbString(ssMaterials.str().c_str()));
1433 redraw = true;
1434 }
1435 if(std::string(This->mouseOverTextZPos->string.getValues(0)->getString()) != ssZPos.str()) {
1436 This->mouseOverTextZPos->string.setValue(SbString(ssZPos.str().c_str()));
1437 redraw = true;
1438 }
1439
1440 if(redraw) {
1441 This->setSuperimpositionEnabled(This->superimposition, TRUE);
1442 This->scheduleRedraw();
1443 }
1444
1445 eventCB->setHandled();
1446 }
1447 else {
1448 if(std::string(This->mouseOverTextLogName->string.getValues(0)->getString()) != "") {
1449 This->mouseOverTextLogName->string.setValue(SbString(""));
1450 This->scheduleRedraw();
1451 }
1452 if(std::string(This->mouseOverTextSolid->string.getValues(0)->getString()) != "") {
1453 This->mouseOverTextSolid->string.setValue(SbString(""));
1454 This->scheduleRedraw();
1455 }
1456 if(std::string(This->mouseOverTextMaterial->string.getValues(0)->getString()) != "") {
1457 This->mouseOverTextMaterial->string.setValue(SbString(""));
1458 This->scheduleRedraw();
1459 }
1460 if(std::string(This->mouseOverTextZPos->string.getValues(0)->getString()) != "") {
1461 This->mouseOverTextZPos->string.setValue(SbString(""));
1462 This->scheduleRedraw();
1463 }
1464 }
1465}
1466
1467
1468// Called by hitting PageUp during animation.
1470 if (std::ceil(animateBtwPtsPeriod * 100) >= 4) {
1471 if (speedStep > 0.08)
1472 speedStep -= 0.02;
1473 else
1474 speedStep = 0.02;
1476 } else
1477 animateBtwPtsPeriod = 0.0;
1478
1480 int lastIdx = (int) refParticleTrajectory.size() - 1;
1483 }
1484}
1485
1486// Called by hitting PageDown during animation.
1490 if (std::floor(animateBtwPtsPeriod * 100) == 12) { // Errors in double representation
1491 speedStep = 0.08;
1492 } else if (animateBtwPtsPeriod > 0.12)
1493 speedStep += 0.02;
1494 } else {
1497 maxSpeed = 0.0f;
1498 if (animateSensor->isScheduled())
1499 animateSensor->unschedule();
1500 }
1501}
1502
1503
1504// Based on the user's interaction the speed indicator bar needs to be adjusted
1505
1507{
1508 assert(this->sgeometry != NULL);
1509
1510 SbVec3f * points = this->sgeometry->point.startEditing();
1511
1512 if (points[10][0] == 0.0f)
1513 this->animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_ALL);
1514 if (points[14][0] == 0.0f)
1515 this->animSpeedSwitch->whichChild.setValue(SO_SWITCH_ALL);
1516 points[10][0] = this->maxSpeed;
1517 points[11][0] = this->maxSpeed;
1518 points[14][0] = this->maxSpeed;
1519 points[15][0] = this->maxSpeed;
1520 this->sgeometry->point.finishEditing();
1521
1522 if (this->maxSpeed == 0.0f) {
1523 this->animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
1524 this->animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
1525 }
1526}
1527
1528
1530 switch (currentState) {
1531 case ANIMATION:
1532 case REVERSED_ANIMATION:
1533 case PAUSED_ANIMATION:
1535 SoQtExaminerViewer::actualRedraw();
1536 break;
1537 default:
1538 SoQtExaminerViewer::actualRedraw();
1539 break;
1540 }
1541}
1542
1543
1545 SoCoordinate3 *coords, bool append)
1546{
1547 // TODO: Color the reference path
1548 // Disable the color stuff for now: changes all trajectories
1549 // FWJ See G4OpenInventorXtExaminerViewer.cc for test code
1550
1551 // The trajectory is composed of all the polyline segments in the
1552 // multiple value field (SoMFInt32) numVertices.
1553 // For each of the numVertices.getNum()* polyline segments,
1554 // retrieve the points from the SoCoordinate3 node
1555
1556 SbVec3f refParticlePt;
1557
1558 if(!append)
1559 refParticleTrajectory.clear();
1560
1561 for(int i = 0; i < lineset->numVertices.getNum(); ++i) {
1562 for(int j = 0; j < lineset->numVertices[i]; ++j) {
1563 refParticlePt = coords->point[j];
1564 refParticleTrajectory.push_back(refParticlePt);
1565 }
1566 }
1567 // Remove points that are too close to each other
1571 sortElements();
1572}
1573
1574
1576{
1577 refZPositions.clear();
1578 refZPositions.push_back(0);
1579 float dist;
1580 for(unsigned int i=0; i < refParticleTrajectory.size() - 1; ++i) {
1581 dist = (refParticleTrajectory[i] -
1582 refParticleTrajectory[i + 1]).length();
1583 refZPositions.push_back(refZPositions[i] + dist);
1584 }
1585}
1586
1587
1589{
1590 SoSearchAction action;
1591 action.setType(SoLineSet::getClassTypeId(),false);
1592 action.setInterest(SoSearchAction::ALL);
1593 action.apply(getSceneGraph());
1594
1595 SoPathList &pathList = action.getPaths();
1596
1597 if(pathList.getLength() != 0) {
1598
1599 SoCoordinate3 * coords = NULL;
1600 std::vector<SoCoordinate3 *> coordvec;
1601 std::vector<SoLineSet *> linevec;
1602
1603 bool refPathFound = false;
1604 for(int i = 0; i < pathList.getLength(); ++i) {
1605 SoFullPath *path = (SoFullPath *)pathList[i];
1606
1607 G4AttHolder* attHolder = dynamic_cast<G4AttHolder*>(path->getTail());
1608 if(attHolder != nullptr)
1609 {
1610 for (std::size_t j = 0; j < attHolder->GetAttDefs().size(); ++j) {
1611 std::ostringstream oss;
1612 oss << G4AttCheck(attHolder->GetAttValues()[j],
1613 attHolder->GetAttDefs()[j]);
1614
1615 std::string findStr = "Type of trajectory (Type): ";
1616 std::string compareValue = "REFERENCE";
1617 std::size_t idx = oss.str().find(findStr);
1618
1619 if(idx != std::string::npos) {
1620 if(oss.str().substr(idx + findStr.size(),
1621 compareValue.size()) == compareValue) {
1622 coords = getCoordsNode(path);
1623 if(coords != NULL) {
1624 refPathFound = true;
1625 coordvec.push_back(coords);
1626 linevec.push_back((SoLineSet *)path->getTail());
1627 }
1628 break;
1629 }
1630 }
1631
1632 findStr = "Track ID (ID): ";
1633 idx = oss.str().find(findStr);
1634 if(idx != std::string::npos) {
1635 //index all primary tracks
1636 std::string tmpstr = oss.str().substr(idx + findStr.size(),1);
1637 std::istringstream buffer(tmpstr);
1638 int num;
1639 buffer >> num;
1640 if(num == 1) {
1641
1642 // Check if next character is a number,
1643 // in which case we don't have Track ID 1
1644 // FWJ attempt to fix Coverity issue.
1645 char nextChar = oss.str().at(idx+findStr.size()+1);
1646 // const char * nextChar =
1647 // oss.str().substr(idx + findStr.size() + 1,1).c_str();
1648 if(std::isdigit(nextChar))
1649 break; //Not a primary track, continue with next track
1650
1651 coords = getCoordsNode(path);
1652 if(coords != NULL) {
1653 coordvec.push_back(coords);
1654 linevec.push_back((SoLineSet *)path->getTail());
1655 break; //Found coords node, continue with next track
1656 }
1657 }
1658 else
1659 break; //Not a primary track, continue with next track
1660 }
1661 else{
1662 //Not a Track ID attribute, fall through
1663 }
1664 }
1665 }
1666
1667 if(refPathFound)
1668 break;
1669 }
1670
1671 if(coordvec.empty())
1672 return; //No track with a Coordinate3 node found
1673
1674 if(refPathFound) {
1675 //set ref path to last traj, coord in the vecs
1676 setReferencePath(linevec.back(), coordvec.back());
1677 return;
1678 }
1679 //else
1680
1681 int longestIdx = 0;
1682 float longestLength = 0.0;
1683 // For all paths
1684 for(unsigned int i=0;i < linevec.size(); ++i) {
1685
1686 //First generate a vector with all the points in this lineset
1687 std::vector<SbVec3f> trajectory;
1688 // For all lines in the i path
1689 for(int j=0; j < linevec[i]->numVertices.getNum(); ++j) {
1690 // For all points in line j
1691 for(int k=0; k < linevec[i]->numVertices[j]; ++k) {
1692 trajectory.push_back(coordvec[i]->point[k]);
1693 }
1694 }
1695
1696 // Then calculate the total length
1697 float tmpLength=0.0;
1698 for(unsigned int j=0; j < trajectory.size() - 1; ++j) {
1699 tmpLength += (trajectory[j] - trajectory[j + 1]).length();
1700 }
1701
1702 if(tmpLength > longestLength) {
1703 longestIdx = i;
1704 longestLength = tmpLength;
1705 }
1706 }
1707
1708 // Set the longest path as the reference path
1709 setReferencePath(linevec[longestIdx], coordvec[longestIdx]);
1710 }
1711}
1712
1713
1714SoCoordinate3 * G4OpenInventorQtExaminerViewer::getCoordsNode(SoFullPath *path)
1715{
1716 SoLineSet *trajectory = (SoLineSet *)path->getTail();
1717 SoSeparator * grpNode = (SoSeparator*)(((SoFullPath*)path)->getNodeFromTail(1));
1718 int nodeIndex = grpNode->findChild(trajectory);
1719 SoNode * tmpNode;
1720
1721 // We allow only 100 iterations, in case the node isn't found
1722 // (should take only a few iterations)
1723 for (int i = 0; i < 100; ++i) {
1724 --nodeIndex;
1725
1726 tmpNode = grpNode->getChild(nodeIndex);
1727 if(tmpNode->getTypeId() == SoCoordinate3::getClassTypeId()) {
1728 //node found
1729 return (SoCoordinate3 *)tmpNode;
1730 }
1731 }
1732 return NULL; //coords node not found
1733}
1734
1735
1736// Displays scene elements on the right side of listsDialog.
1737// else: scene graph is searched for Geant4_SoPolyhedron type nodes
1739{
1740 std::string field, eltName;
1741
1742 std::map<std::string, int> duplicates;
1743 std::map<std::string, int> sceneElts;
1744 SoSearchAction search;
1745 Geant4_SoPolyhedron *node;
1746 SoGroup *root = (SoGroup *)getSceneManager()->getSceneGraph();
1747
1748 SoBaseKit::setSearchingChildren(TRUE);
1749
1750 search.reset();
1751 search.setSearchingAll(TRUE);
1752 search.setInterest(SoSearchAction::ALL);
1753 search.setType(Geant4_SoPolyhedron::getClassTypeId(), 0);
1754
1755 // FWJ DEBUG
1756 // G4cout << "Searching for elements....." << G4endl;
1757 search.apply(root);
1758
1759 SoPathList &pl = search.getPaths();
1760
1761
1762 // First find which names occur more than once so we can append a counter to them
1763 for (int i = 0; i < pl.getLength(); i++) {
1764 SoFullPath *path = (SoFullPath *)pl[i];
1765 node = (Geant4_SoPolyhedron *)path->getTail();
1766 eltName = node->getName();
1767 // G4cout << " FOUND " << i << " " << eltName << G4endl;
1768 if(duplicates.count(eltName))
1769 duplicates[eltName]++;
1770 else
1771 duplicates[eltName] = 1;
1772 }
1773
1774 for(int i = 0; i < pl.getLength(); i++) {
1775 float x,y,z;
1776 std::stringstream ssCount;
1777 SoFullPath *path = (SoFullPath *)pl[i];
1778 node = (Geant4_SoPolyhedron *)path->getTail();
1779 eltName = node->getName();
1780 field = eltName;
1781 if(duplicates[eltName] == 1)
1782 ssCount << "";//duplicates[field]
1783 else {
1784 if(sceneElts.count(eltName))
1785 sceneElts[eltName]++;
1786 else
1787 sceneElts[eltName] = 1;
1788
1789 ssCount << sceneElts[eltName];
1790 field += "_";
1791 }
1792
1793 field += ssCount.str();
1794
1795 SoGetBoundingBoxAction bAction(getViewportRegion());
1796 bAction.apply(path);
1797 SbBox3f bBox = bAction.getBoundingBox();
1798
1799 SbVec3f centr = bBox.getCenter();
1800 centr.getValue(x,y,z);
1801
1802 path->ref();
1803 sceneElement el = { field, path, centr, 0.0 };
1804 sceneElements.push_back(el);
1805 }
1806}
1807
1808
1810{
1811 float x,y,z;
1812 a.getValue(x,y,z);
1813 return x*x + y*y + z*z;
1814}
1815
1816
1818 float &dist,
1819 SbVec3f &closestPoint,
1820 int &index)
1821{
1822 // a : Previous point on trajectory
1823 // b : Next point on trajectory
1824 // q : the point in space
1825 // dab, daq, dbq: distance between a & b, a & q, b & q
1826 //
1827 // Theory: A point p on a line ab is defined as:
1828 //
1829 // p(t) = a+t?(b?a)
1830 //
1831 // note: All are vectors except the parameter t
1832 //
1833 // When t is between 0 and 1 the point p is situated between a and b on ab.
1834 // The point p is defined in terms of the parameter t, subsequently so does
1835 // the distance from the query point q to the point p. To find the minimum
1836 // of that distance we differentiate it and set equal to zero:
1837 //
1838 // diff(Norm(p(t)- q)) = 0
1839 //
1840 // note: diff means taking the derivative with regard to t
1841 //
1842 // The resulting t is given in the code below. The square of the distance
1843 // between p and q is given by:
1844 //
1845 // d^2 = (Norm(p(t)-q))^2
1846 //
1847 // The expression found is given in the code below (current_dist)
1848 //
1849 // Ref: http://programmizm.sourceforge.net/blog/2012/
1850 // distance-from-a-point-to-a-polyline
1851 //
1852 // --PLG
1853
1854 const std::size_t count = refParticleTrajectory.size();
1855 assert(count>0);
1856
1857 SbVec3f b = refParticleTrajectory[0];
1858 SbVec3f dbq = b - q;
1859 float sqrDist = sqrlen(dbq);
1860 closestPoint = b;
1861 index = 0;
1862 for (std::size_t i = 1; i < count; ++i) {
1863 const SbVec3f a = b;
1864 const SbVec3f daq = dbq;
1865 b = refParticleTrajectory[i];
1866 dbq = b - q;
1867 const SbVec3f dab = a - b;
1868
1869 float dab_x, dab_y, dab_z;
1870 dab.getValue(dab_x,dab_y,dab_z);
1871 float daq_x, daq_y, daq_z;
1872 daq.getValue(daq_x, daq_y, daq_z);
1873 float dbq_x, dbq_y, dbq_z;
1874 dbq.getValue(dbq_x, dbq_y, dbq_z);
1875
1876 const float inv_sqrlen = 1./sqrlen(dab);
1877 const float t = (dab_x*daq_x + dab_y*daq_y + dab_z*daq_z)*inv_sqrlen;
1878
1879 if (t<0.) {
1880 // The trajectory point occurs before point a
1881 // Go to the next point
1882 continue;
1883 }
1884 float current_dist;
1885 if (t<=1.) {
1886 // The trajectory point occurs between a and b.
1887 // Compute the distance to that point
1888 current_dist = daq_x*daq_x + daq_y*daq_y + daq_z*daq_z
1889 - t*(daq_x*dab_x + daq_y*dab_y + daq_z*dab_z)
1890 + t*t*(dab_x*dab_x + dab_y*dab_y + dab_z*dab_z);
1891 }
1892 else { //t>1.
1893 // The trajectory point occurs after b.
1894 // Get the distance to point b
1895 current_dist = sqrlen(dbq);
1896 }
1897
1898 if (current_dist < sqrDist) {
1899 sqrDist = current_dist;
1900 closestPoint = a + t*(b-a);
1901 index = (int) i;
1902 }
1903 }
1904
1905 dist = std::sqrt(sqrDist);
1906}
1907
1908
1910{
1911 if(refParticleTrajectory.empty())
1912 return;
1913
1914 float * trajLength = new float[refParticleTrajectory.size()];
1915 typedef std::map<elementForSorting, sceneElement> sortedMap;
1916 sortedMap sorted;
1917
1918 // For every point on the reference trajectory, compute
1919 // the total length from the start
1920 SbVec3f prevPoint;
1921 std::vector<SbVec3f>::iterator itRef = refParticleTrajectory.begin();
1922 int trajIndex = 0;
1923 prevPoint = *itRef;
1924 trajLength[trajIndex] = 0.0;
1925 ++itRef;
1926 ++trajIndex;
1927 for(; itRef != refParticleTrajectory.end(); ++itRef, ++trajIndex) {
1928 trajLength[trajIndex] = trajLength[trajIndex-1] + (*itRef - prevPoint).length();
1929 prevPoint = *itRef;
1930 }
1931
1932 // Compute the smallest distance between the element
1933 // and the reference trajectory (find the closest point),
1934 // then map the element to the trajectory length of that
1935 // point (calculated above)
1936 SoGetBoundingBoxAction bAction(getViewportRegion());
1937 SbVec3f elementCoord;
1938 std::vector<sceneElement>::iterator itEl;
1939 int elementIndex;
1941 for(itEl = sceneElements.begin(), elementIndex = 0;
1942 itEl != sceneElements.end(); ++itEl, ++elementIndex) {
1943 bAction.apply(itEl->path);
1944
1945 // FWJ sceneElement already has a center
1946 elementCoord = itEl->center;
1947 // ... and this sometimes returns an empty box!
1948 // elementCoord = bAction.getBoundingBox().getCenter();
1949 // if (bAction.getBoundingBox().isEmpty()) {
1950 // G4cout << "sortElements: Box is empty!" << G4endl;
1951 // G4cout << " element name=" << itEl->name << G4endl;
1952 // }
1953
1954 int index;
1955 distanceToTrajectory(elementCoord, el.smallestDistance, el.closestPoint, index);
1956 itEl->closestPointZCoord = el.closestPointZCoord = trajLength[index];
1957 el.distanceToBeamlineStart = (itEl->center - refParticleTrajectory[0]).length();
1958
1959 // This map of the scene elements (or their coordinates rather)
1960 // is automatically sorted by trajectory length (Z coord), then
1961 // by the distance between the element and the point in case the Z coord
1962 // is the same as another element. This is done by using as a key
1963 // an element structure which implements the operator for weak ordering
1964 sorted.insert(std::make_pair(el,*itEl));
1965 }
1966
1967 // store the sorted elements into the vector field
1968 sceneElements.clear();
1969
1970 sortedMap::iterator itSorted = sorted.begin();
1971 for(; itSorted != sorted.end(); itSorted++)
1972 sceneElements.push_back(itSorted->second);
1973
1974 zcoordSetFlag = true;
1975
1977
1978 delete[] trajLength;
1979}
1980
1981
1983{
1984 // FWJ DEBUG
1985 // G4cout << "Populating ELEMENT LIST..." << G4endl;
1986
1987 AuxWindowDialog->listWidget1->clear();
1988 // int size = sceneElements.size();
1989
1990 std::vector<sceneElement>::const_iterator it;
1991 std::stringstream ss;
1992
1993 for(it=sceneElements.begin(); it!=sceneElements.end(); ++it) {
1994 ss << it->name;
1995 if(zcoordSetFlag)
1996 ss << " [" << it->closestPointZCoord << "]";
1997
1998 new QListWidgetItem(ss.str().c_str(), AuxWindowDialog->listWidget1);
1999 ss.str("");
2000 }
2001}
2002
2003
2004// Called when user clicks a scene element in listsDialog.
2005// Zooms onto that element.
2006void
2007G4OpenInventorQtExaminerViewer::LookAtSceneElementCB(QListWidgetItem* item)
2008{
2009 // FWJ DEBUG
2010 // G4cout << "AuxWindow: listWidget1 select element CALLBACK" << G4endl;
2011
2012 SoCamera * cam = getCamera();
2013
2014 if (SoQtExaminerViewer::isAnimating())
2015 stopAnimating();
2016
2019 if (animateSensor->isScheduled())
2020 animateSensor->unschedule();
2021 setSuperimpositionEnabled(superimposition, FALSE);
2022 maxSpeed = 0.0f;
2023 scheduleRedraw();
2024 restoreCamera();
2026 } else if (currentState == VIEWPOINT)
2027 setSuperimpositionEnabled(superimposition, FALSE);
2028
2029 std::string elementField = qPrintable(item->text());
2030
2031 std::size_t idx = elementField.find_last_of("[");
2032 if(idx == std::string::npos)
2033 idx = elementField.size(); //if "[" not found for whatever reason (list not sorted)
2034 else
2035 idx--; // To get rid of the space that is between the name and '['
2036
2037 bool error = false;
2038 SoFullPath *path;
2039 SoSearchAction search;
2040 SoNode *root = getSceneManager()->getSceneGraph();
2041 int counter = 0;
2042 std::size_t idxUnderscore = elementField.find_last_of("_");
2043
2044 parseString<int>(counter,
2045 elementField.substr(idxUnderscore + 1, idx), error);
2046
2047 SoBaseKit::setSearchingChildren(TRUE);
2048 search.reset();
2049 search.setSearchingAll(TRUE);
2050
2051 // G4cout << " Starting search for elementField " << elementField
2052 // << G4endl;
2053
2054 if(error) { // No counter is present => element name was not modified
2055 curEltName = elementField.substr(0, idx);
2056 search.setName(curEltName.c_str());
2057 search.apply(root);
2058
2059 path = (SoFullPath *)search.getPath();
2060 }
2061 else {
2062 curEltName = elementField.substr(0, idxUnderscore);
2063 search.setInterest(SoSearchAction::ALL);
2064 search.setName(curEltName.c_str());
2065 search.apply(root);
2066
2067 SoPathList &pl = search.getPaths();
2068 path = (SoFullPath *)pl[counter - 1]; // Since counter starts at 1, not 0
2069 }
2070
2071 G4ThreeVector global;
2072
2073 // FWJ FLIP THIS
2074 if ((idx > 0) && (path)) {
2075
2076 if(!refParticleTrajectory.empty()) {
2077
2078 SoGetBoundingBoxAction bAction(getViewportRegion());
2079 bAction.apply(path);
2080 SbBox3f bBox = bAction.getBoundingBox();
2081 SbVec3f elementCoord = bBox.getCenter();
2082
2083 refParticleIdx = 0;
2084 SbVec3f p;
2085
2086 float absLengthNow, absLengthMin;
2087 int maxIdx = (int) refParticleTrajectory.size() - 2;
2088 int targetIdx = 0;
2089 SbVec3f dir;
2090
2092 absLengthMin = (p - elementCoord).length();
2094
2095 // Find a ref. particle's point closest to element's global coords
2096 while (refParticleIdx < maxIdx) {
2098 absLengthNow = (p - elementCoord).length();
2099
2100 if (absLengthNow < absLengthMin) {
2101 absLengthMin = absLengthNow;
2102 targetIdx = refParticleIdx;
2103 }
2105 }
2106
2107 if (currentState != BEAMLINE) { // Set up default zoom
2108 SbVec3f p1, pN;
2110 prevParticleDir = SbVec3f(0,0,0); //so that moveCamera() knows sets default parameters
2111
2112 p1 = prevPt = refParticleTrajectory[0];
2114 distance = (pN - p1).length() / 10;
2115
2116 // FWJ Rather than switching to a default height, it is more flexible
2117 // to keep the same height(magnification) while moving the camera.
2118 // if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2119 // ((SoOrthographicCamera *) cam)->height.setValue(defaultHeight);
2120 // // FWJ Restore the default height instead of hard-wired value
2121 // // ((SoOrthographicCamera *) cam)->height.setValue(10000.0f);
2122 // }
2123 // else if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2124
2125 // FWJ required to avoid extreme perspective after camera move:
2126 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2127 ((SoPerspectiveCamera*)cam)->heightAngle.setValue(defaultHeightAngle);
2128
2129 } else {
2130 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2131 distance = (prevPt - cam->position.getValue()).length();
2132 }
2133 refParticleIdx = targetIdx;
2134
2135 //////////////////////////////////////////////////////////////
2136 setSuperimpositionEnabled(superimposition, TRUE);
2137 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
2138 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
2139 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
2140 scheduleRedraw();
2141 //////////////////////////////////////////////////////////////
2142
2144
2145 }
2146
2147 else {
2148 offsetFromCenter.setValue(0, 0, 1);
2149 distance = 50;// small number since using viewAll() for default zoom
2150 upVector.setValue(0, 1, 0);
2152 cam->viewAll(path, getViewportRegion());
2153 }
2154 }
2155
2156}
2157
2158
2159void G4OpenInventorQtExaminerViewer::FileLoadRefPathCB()
2160{
2161 // G4cout << "File: Load Ref Path CALLBACK" << G4endl;
2162
2163 QFileDialog filedialog(getParentWidget(), tr("Load Reference Path"));
2164 filedialog.setFileMode(QFileDialog::AnyFile);
2165 filedialog.setFont(*font);
2166 if (!filedialog.exec()) return;
2167 QStringList filenameinlist = filedialog.selectedFiles();
2168 QString filenamein = filenameinlist[0];
2169
2170 std::ifstream ifs(qPrintable(filenamein));
2171 if(ifs.is_open()) {
2172 refParticleTrajectory.clear();
2173 float x,y,z;
2174 while(ifs >> x >> y >> z) {
2175 refParticleTrajectory.push_back(SbVec3f(x,y,z));
2176 }
2177 ifs.close();
2178 } else {
2179 QMessageBox msgbox;
2180 msgbox.setFont(*font);
2181 QString messagetxt = "Reference Path file not found: ";
2182 messagetxt.append(filenamein);
2183 msgbox.setText(messagetxt);
2184 msgbox.exec();
2185 return;
2186 }
2187 if (refParticleTrajectory.size() < 2) {
2188 QMessageBox msgbox;
2189 msgbox.setFont(*font);
2190 QString messagetxt = "Invalid Reference Path";
2191 msgbox.setText(messagetxt);
2192 msgbox.exec();
2193 return;
2194 }
2195 // Following setReferencePath() ...
2199 sortElements();
2200}
2201
2202
2203void G4OpenInventorQtExaminerViewer::FileSaveRefPathCB()
2204{
2205 // G4cout << "File: Save Ref Path CALLBACK" << G4endl;
2206
2207 QFileDialog filedialog(getParentWidget(), tr("Save Reference Path"));
2208 filedialog.setFileMode(QFileDialog::AnyFile);
2209 // To enable confirmation of overwriting
2210 filedialog.setAcceptMode(QFileDialog::AcceptSave);
2211 filedialog.setFont(*font);
2212 if (!filedialog.exec()) return;
2213 QStringList filenameinlist = filedialog.selectedFiles();
2214 QString filenamein = filenameinlist[0];
2215
2216 std::ofstream ofs(qPrintable(filenamein));
2217 if (ofs.is_open()) {
2218 float x,y,z;
2219 for (unsigned int i=0; i < refParticleTrajectory.size(); ++i) {
2220 refParticleTrajectory[i].getValue(x,y,z);
2221 ofs << x << " " << y << " " << z << "\n";
2222 }
2223 ofs.close();
2224 } else {
2225 QMessageBox msgbox;
2226 msgbox.setFont(*font);
2227 QString messagetxt = "Error opening file ";
2228 messagetxt.append(filenamein);
2229 msgbox.setText(messagetxt);
2230 msgbox.exec();
2231 }
2232
2233}
2234
2236{
2237 if(refParticleTrajectory.empty())
2238 return;
2239
2240 SbVec3f p1, p2, p3, dirNow, dirNxt, dir, p2_tmp, p_start, p_corner, p_nxt;
2241 float avgDistBtwPts = 0;
2242 float totalDistBtwPts = 0;
2243 std::vector<SbVec3f> newRefParticleTrajectory;
2244 SbVec3f refPoint;
2245 std::size_t size = refParticleTrajectory.size() - 1;
2246 int numOfPts = 0;
2247 for (std::size_t i = 0; i < size; ++i) {
2248 p1 = refParticleTrajectory[i];
2249 p2 = refParticleTrajectory[i + 1];
2250 if (p1 == p2)
2251 continue;
2252 numOfPts++;
2253 totalDistBtwPts += (p2 - p1).length();
2254 }
2255 // Nothing useful to do (and fix Coverity)
2256 if (numOfPts <= 2) return;
2257
2258 avgDistBtwPts = totalDistBtwPts / numOfPts;
2259 float minDistAllowed = 0.75 * avgDistBtwPts;
2260 // float maxDistAllowed = 1.25 * avgDistBtwPts; // Pts tend to be close not far
2261
2262 float x, y, z;
2263 std::size_t i = 0, j = 0;
2264 while (i < size) {
2265 p1 = refParticleTrajectory[i];
2266 p2 = refParticleTrajectory[i + 1];
2267
2268 refPoint = p1;
2269 p1.getValue(x, y, z);
2270
2271 newRefParticleTrajectory.push_back(refPoint);
2272
2273 j = i;
2274 while ((p2 - p1).length() < minDistAllowed && j < (size - 1)) {
2275 j++;
2276
2277 p1 = refParticleTrajectory[j];
2278 p2 = refParticleTrajectory[j + 1];
2279 }
2280 if (j != i)
2281 i = j + 1;
2282 else
2283 i++;
2284 }
2285
2286 refParticleTrajectory.clear();
2287 refParticleTrajectory = newRefParticleTrajectory;
2288}
2289
2290
2292{
2293 SoCamera *cam = getCamera();
2294 camB4Animation.viewportMapping = cam->viewportMapping.getValue();
2295 camB4Animation.position = cam->position.getValue();
2296 camB4Animation.orientation = cam->orientation.getValue();
2297 camB4Animation.aspectRatio = cam->aspectRatio.getValue();
2298 camB4Animation.nearDistance = cam->nearDistance.getValue();
2299 camB4Animation.farDistance = cam->farDistance.getValue();
2300 camB4Animation.focalDistance = cam->focalDistance.getValue();
2301
2302 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2304 ((SoPerspectiveCamera *) cam)->heightAngle.getValue();
2306 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2308 ((SoOrthographicCamera *) cam)->height.getValue();
2310 }
2311}
2312
2313
2315{
2316 SoCamera *cam = getCamera();
2317
2318 cam->viewportMapping = camB4Animation.viewportMapping;
2319 cam->position = camB4Animation.position;
2320 cam->orientation = camB4Animation.orientation;
2321 cam->aspectRatio = camB4Animation.aspectRatio;
2322 cam->nearDistance = camB4Animation.nearDistance;
2323 cam->farDistance = camB4Animation.farDistance;
2324 cam->focalDistance = camB4Animation.focalDistance;
2325
2326 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2328 toggleCameraType();
2329 cam = getCamera();
2330 ((SoOrthographicCamera *) cam)->height.setValue(
2332 } else
2333 ((SoPerspectiveCamera *) cam)->heightAngle.setValue(
2335 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2337 toggleCameraType();
2338 cam = getCamera();
2339 ((SoPerspectiveCamera *) cam)->heightAngle.setValue(
2341 } else
2342 ((SoOrthographicCamera *) cam)->height.setValue(
2344 }
2345}
2346
2347
2349 SoSensor *sensor)
2350{
2351 SbTime curTime = SbTime::getTimeOfDay();
2353
2354 SoTimerSensor* s = (SoTimerSensor*) sensor;
2355
2356 float t = float((curTime - s->getBaseTime()).getValue())
2357 / This->animateBtwPtsPeriod;
2358
2359 if ((t > 1.0f) || (t + s->getInterval().getValue() > 1.0f))
2360 t = 1.0f;
2361 SbBool end = (t == 1.0f);
2362
2363 if (end) {
2364 This->animateSensorRotation->unschedule();
2365 if(This->rotCnt) {
2366 // rotations left
2367 This->rotateCamera();
2368 }
2369 else {
2370 // rotation over
2371 This->currentState = This->prevState;
2372 return;
2373 }
2374 }
2375
2376}
2377
2378
2379// Called repeatedly during reference particle animation
2380
2382 SoSensor *sensor)
2383{
2384 SbTime curTime = SbTime::getTimeOfDay();
2386 SoCamera *cam = This->getCamera();
2387 SoTimerSensor* s = (SoTimerSensor*) sensor;
2388
2389 float t = float((curTime - s->getBaseTime()).getValue())
2390 / This->animateBtwPtsPeriod;
2391
2392 if ((t > 1.0f) || (t + s->getInterval().getValue() > 1.0f))
2393 t = 1.0f;
2394 SbBool end = (t == 1.0f);
2395
2396 cam->orientation = SbRotation::slerp(This->camStartOrient, This->camEndOrient, t);
2397 cam->position = This->camStartPos + (This->camEndPos - This->camStartPos) * t;
2398
2399 if (end) {
2400 This->animateSensor->unschedule();
2401
2402 if (This->currentState == ANIMATION) {
2403 if (This->refParticleIdx < (int) (This->refParticleTrajectory.size() - 1))
2404 This->animateRefParticle();
2405 else {
2407 This->speedStep = START_STEP;
2408 }
2409 }
2410 if (This->currentState == REVERSED_ANIMATION) {
2411 if (This->refParticleIdx >= 1)
2412 This->animateRefParticle();
2413 else {
2415 This->speedStep = START_STEP;
2416 }
2417 }
2418 }
2419}
2420
2421
2423{
2424 if (SoQtExaminerViewer::isAnimating())
2425 stopAnimating();
2426
2427 SbRotation rot;
2428 SbVec3f p1{0.0, 0.0, 0.0}, p2{0.0, 0.0, 0.0}, p2_tmp, camUpV, camD, camD_tmp, leftRightAxis;
2429 float x1, y1, z1, x2, y2, z2;
2430
2431 if (currentState == ANIMATION) {
2434 } else if (currentState == REVERSED_ANIMATION) {
2437 } else if (currentState == PAUSED_ANIMATION) {
2438 if (refParticleIdx < (int) refParticleTrajectory.size()) {
2441 } else {
2444 }
2445 }
2446 p1.getValue(x1, y1, z1);
2447 p2.getValue(x2, y2, z2);
2448
2449 camD = p2 - p1;
2450 camD.normalize();
2451
2452 p2_tmp.setValue(x2, y1, z2);
2453 camD_tmp = p2_tmp - p1;
2454 camD_tmp.normalize();
2455
2456 camUpV.setValue(0, 1, 0);
2457 rot.setValue(camD_tmp, camD);
2458 rot.multVec(camUpV, camUpV);
2459
2460 leftRightAxis = camD.cross(camUpV);
2461
2462 myCam->position = p1;
2463 myCam->pointAt(p2, camUpV);
2464
2465 // Update camera position
2466 p1 = p1 + (up_down * camUpV) + (left_right * leftRightAxis);
2467 myCam->position = p1;
2468 // FWJ Try look-ahead here
2469 int idx = refParticleIdx + pathLookahead;
2470 idx = std::min(idx, (int)refParticleTrajectory.size() - 1);
2471 myCam->pointAt(refParticleTrajectory[idx], camUpV);
2472 // myCam->pointAt(refParticleTrajectory[idx], camUpVec);
2473 myCam->focalDistance = 0.1f;
2474}
2475
2476
2478{
2479 G4OpenInventorQtExaminerViewer::ToolsRefPathStartCB();
2480}
2481
2482
2483void G4OpenInventorQtExaminerViewer::ToolsRefPathStartCB()
2484{
2485 if (!refParticleTrajectory.size()) {
2486 QMessageBox msgbox;
2487 msgbox.setFont(*font);
2488 QString messagetxt = "No current reference path";
2489 msgbox.setText(messagetxt);
2490 msgbox.exec();
2491 return;
2492 }
2493
2494 if (currentState == ROTATING)
2495 return;
2498 if (animateSensor->isScheduled())
2499 animateSensor->unschedule();
2500 setSuperimpositionEnabled(superimposition, FALSE);
2501 maxSpeed = 0.0f;
2502 scheduleRedraw();
2503 } else {
2504 saveCurCamera();
2507 }
2508
2509 if (SoQtExaminerViewer::isAnimating())
2510 stopAnimating();
2511
2512 up_down = 0;
2513 left_right = 0;
2514 step = 1;
2515
2516 refParticleIdx = 0;
2518 setSuperimpositionEnabled(superimposition, TRUE);
2519 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
2520 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
2521 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
2522 scheduleRedraw();
2523
2524 // FWJ Disabled: this is set in moveCamera()
2525 // Zoom in
2526 // SoCamera *cam = getCamera();
2527 // cam->focalDistance = 0.1f;
2528
2529 prevParticleDir = SbVec3f(0,0,0);
2530
2531 //Default zoom
2532 SbVec3f p1 = refParticleTrajectory[0];
2533 SbVec3f pN = refParticleTrajectory[refParticleTrajectory.size() - 1];
2534 distance = (pN - p1).length() / 10;
2535
2536 moveCamera(distance, true);
2537}
2538
2539
2540void G4OpenInventorQtExaminerViewer::ToolsRefPathInvertCB()
2541{
2542 invertRefPath();
2543}
2544
2545
2553
2554
2556{
2557 SoCamera *cam = getCamera();
2558
2559 camStartPos = cam->position.getValue();
2560 camStartOrient = cam->orientation.getValue();
2561
2562 if (currentState != BEAMLINE)
2564
2565 camEndPos = myCam->position.getValue();
2566 camEndOrient = myCam->orientation.getValue();
2567
2568 if (animateSensor->isScheduled())
2569 animateSensor->unschedule();
2570
2571 animateSensor->setBaseTime(SbTime::getTimeOfDay());
2572 animateSensor->setInterval(SbTime(0.02));
2573
2574 animateSensor->schedule();
2575}
2576
2577
2579{
2580 escapeCallback = callback;
2581}
2582
2583
2585{
2586 // FWJ DEBUG
2587 // G4cout << "SCENE CHANGE callback" << G4endl;
2588 // NOTE: could/should be disabled during animation
2589
2592 if(This->newEvents) {
2593 This->findAndSetRefPath();
2594 This->newEvents = false;
2595 }
2596}
2597
2598
2599//////////////////////////////////// BOOKMARKS ///////////////////////////
2600
2601// Adds bookmarks to listsDialog.
2602
2604{
2605 std::size_t size = viewPtList.size();
2606 if (!size) return;
2607
2608 for (std::size_t i = 0; i < size; ++i) {
2609 new QListWidgetItem(viewPtList[i].viewPtName,
2610 AuxWindowDialog->listWidget);
2611 }
2612}
2613
2614
2615// Converts a string type word into a float type.
2616
2617template<class T>
2618void G4OpenInventorQtExaminerViewer::parseString(T &t, const std::string &s,
2619 bool &error)
2620{
2621 std::istringstream str(s);
2622 if ((str >> t).fail())
2623 error = true;
2624}
2625
2626
2627void
2628G4OpenInventorQtExaminerViewer::FileOpenBookmarkCB()
2629{
2630 // FWJ DEBUG
2631 // G4cout << "File: Open Bookmark File CALLBACK" << G4endl;
2632 QFileDialog filedialog(getParentWidget(), tr("Open bookmark file"));
2633 filedialog.setFileMode(QFileDialog::ExistingFile);
2634 filedialog.setFont(*font);
2635 if (!filedialog.exec()) return;
2636 QStringList filenameinlist = filedialog.selectedFiles();
2637 QString filenamein = filenameinlist[0];
2638
2639 fileIn.close();
2640 fileIn.open(qPrintable(filenamein));
2641 if (fileIn.fail()) {
2642 QMessageBox msgbox;
2643 msgbox.setFont(*font);
2644 QString messagetxt = "Error opening file: ";
2645 messagetxt.append(filenamein);
2646 msgbox.setText(messagetxt);
2647 msgbox.exec();
2648 // G4cout << "ERROR opening file " << filename << G4endl;
2649 fileIn.clear();
2650 return;
2651 }
2652 // Opens a file without erasing it
2654
2655 if (!loadViewPts()) {
2656 QMessageBox msgbox;
2657 msgbox.setFont(*font);
2658 QString messagetxt = "Error reading bookmark file: ";
2659 messagetxt.append(filenamein);
2660 msgbox.setText(messagetxt);
2661 msgbox.exec();
2662 // G4cout << "ERROR reading bookmark file " << filename << G4endl;
2663 fileIn.clear();
2664 return;
2665 }
2666
2667 fileName = qPrintable(filenamein);
2668 fileOut.open(fileName.c_str(), std::ios::in);
2669 fileOut.seekp(0, std::ios::end);
2670
2671 addViewPoints();
2672
2673 // LATER: display filename in lists window
2674
2675 fileIn.close();
2676 fileIn.clear();
2677}
2678
2679// Called before loading a new viewpoint file.
2680// Resets member fields to default values.
2681
2683{
2684 viewPtIdx = -1;
2685 viewPtList.clear();
2686 // setSuperimpositionEnabled(superimposition, FALSE);
2687 // scheduleRedraw();
2689 if (fileOut.is_open()) fileOut.close();
2690
2691 AuxWindowDialog->listWidget->clear();
2692 AuxWindowDialog->lineEdit->setText(QString(""));
2693}
2694
2695
2696void
2697G4OpenInventorQtExaminerViewer::FileNewBookmarkCB()
2698{
2699 // G4cout << "File: Open New Bookmark File CALLBACK" << G4endl;
2700 QFileDialog filedialog(getParentWidget(), tr("Open new bookmark file"));
2701 filedialog.setFileMode(QFileDialog::AnyFile);
2702 // To enable confirmation of overwriting
2703 filedialog.setAcceptMode(QFileDialog::AcceptSave);
2704 // But change the "Save" button text
2705 filedialog.setLabelText(QFileDialog::Accept, QString("New"));
2706 filedialog.setFont(*font);
2707 if (!filedialog.exec()) return;
2708 QStringList filenameinlist = filedialog.selectedFiles();
2709 QString filenamein = filenameinlist[0];
2710
2712 fileName = qPrintable(filenamein);
2713 fileOut.open(fileName.c_str());
2714 if (fileOut.fail()) {
2715 QMessageBox msgbox;
2716 msgbox.setFont(*font);
2717 QString messagetxt = "Error opening new bookmark file: ";
2718 messagetxt.append(filenamein);
2719 msgbox.setText(messagetxt);
2720 msgbox.exec();
2721 // G4cout << "ERROR opening new bookmark file " << filename << G4endl;
2722 }
2723}
2724
2725
2726void
2727G4OpenInventorQtExaminerViewer::ToolsAnimateRefParticleCB()
2728{
2729 // G4cout << "Tools: Animate Ref Particle CALLBACK" << G4endl;
2730 if (!refParticleTrajectory.size()) {
2731 returnToAnim = true;
2732 G4warn << "No Reference Trajectory" << G4endl;
2733 return;
2734 }
2735
2736 ///////////////////////////////////////////////////////////////
2737 setSuperimpositionEnabled(superimposition, TRUE);
2739 axisSwitch->whichChild.setValue(SO_SWITCH_ALL);
2740 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_ALL);
2741 animSpeedSwitch->whichChild.setValue(SO_SWITCH_ALL);
2742 scheduleRedraw();
2743 ///////////////////////////////////////////////////////////////
2744
2745 SoCamera *cam = getCamera();
2746 // SbVec3f camDirOld, camDirNew, camDirNew_tmp, camUpVec, P0, P1, P1_tmp;
2747
2749 || currentState == ROTATING)
2750 return;
2751
2753
2754 saveCurCamera();
2757
2758 if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2759 toggleCameraType();
2760 cam = getCamera();
2761 }
2762
2763 refParticleIdx = 0; // Set the camera to the starting point of the animation
2766 left_right = up_down = 0;
2767
2768 cam->focalDistance = 0.1f;
2769 ((SoPerspectiveCamera *) cam)->heightAngle = 0.50f;
2770 }
2771
2774
2775 cam->position = (myCam)->position.getValue();
2776 cam->orientation = (myCam)->orientation.getValue();
2777 animateRefParticle(); // Animate the camera
2778}
2779
2780
2781void
2782G4OpenInventorQtExaminerViewer::SaveViewPtCB()
2783{
2784 // G4cout << "AppButton: Save Viewpoint CALLBACK" << G4endl;
2785 // First get viewpoint name ...
2786 // EMULATING getViewPtNameCB ...
2787 // bool ok;
2788 // Note QString() returns an empty string
2789
2790 // NONE OF THE FOLLOWING CHANGES THE FONT: FORGET IT FOR NOW
2791 QInputDialog* inputdialog = new QInputDialog(getParentWidget());
2792 inputdialog->setFont(*font);
2793 inputdialog->setWindowTitle(tr("Enter a name for the bookmark"));
2794 inputdialog->setLabelText("Bookmark name");
2795 // inputdialog->setTextEchoMode(QLineEdit::Normal);
2796 inputdialog->adjustSize();
2797 QString namein;
2798 if (inputdialog->exec() == QDialog::Accepted)
2799 namein=inputdialog->textValue().trimmed();
2800 else
2801 return;
2802 if (namein.isEmpty()) return;
2803
2804 // This easier approach failed: unable to set font size
2805 // QString namein =
2806 // QInputDialog::getText(getParentWidget(),
2807 // tr("Enter a name for the bookmark"),
2808 // tr("Bookmark name"), QLineEdit::Normal,
2809 // QString(), &ok);
2810
2811 namein.truncate(MAX_VP_NAME);
2812
2813 char* name = strdup(qPrintable(namein));
2814
2815 // FWJ DEBUG
2816 // G4cout << "QString is " << qPrintable(namein) << G4endl;
2817 // G4cout << "char[] is " << name << G4endl;
2818
2819 for (int i = 0; i < (int)viewPtList.size(); i++) {
2820 if (!strcmp(name, viewPtList[i].viewPtName)) {
2821 QMessageBox msgbox;
2822 msgbox.setText("Bookmark name is already in use");
2823 msgbox.exec();
2824 free(name);
2825 return;
2826 }
2827 }
2828
2829 if (viewPtIdx == -1) viewPtIdx = 0;
2830 saveViewPt(name);
2831
2832 saveViewPtItem = new QListWidgetItem(namein,
2833 AuxWindowDialog->listWidget);
2834 AuxWindowDialog->listWidget->setCurrentItem(saveViewPtItem);
2835 AuxWindowDialog->lineEdit->setText(namein);
2836 free(name);
2837}
2838
2839
2840// Saves current camera parameters to a viewpoint file.
2841
2843{
2844 SbVec3f axis;
2845 viewPtData tmp;
2846 float x, y, z, angle;
2847 SoCamera* camera = getCamera();
2848
2849 // NOTE: Xt VSN increments this at end of procedure
2850 // viewPtIdx++;
2851
2852 // FWJ DEBUG
2853 // G4cout << "saveViewPt: saving bookmark " << viewPtIdx << " " << name
2854 // << G4endl;
2855
2856 if (viewPtList.size() == 0) {
2858 }
2859
2860 tmp.viewPtName = name;
2861 tmp.viewportMapping = camera->viewportMapping.getValue();
2862 tmp.position = camera->position.getValue();
2863 tmp.orientation = camera->orientation.getValue();
2864 tmp.aspectRatio = camera->aspectRatio.getValue();
2865 tmp.nearDistance = camera->nearDistance.getValue();
2866 tmp.farDistance = camera->farDistance.getValue();
2867 tmp.focalDistance = camera->focalDistance.getValue();
2868
2869 // Save camera height (changed by zooming)
2870 if (camera->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2871 tmp.height = ((SoPerspectiveCamera *) camera)->heightAngle.getValue();
2872 tmp.camType = PERSPECTIVE;
2873 } else if (camera->isOfType(SoOrthographicCamera::getClassTypeId())) {
2874 tmp.height = ((SoOrthographicCamera *) camera)->height.getValue();
2875 tmp.camType = ORTHOGRAPHIC;
2876 } else {
2877 SoDebugError::post("G4OpenInventorQtExaminerViewer::saveViewPtCB",
2878 "Only Perspective and Orthographic cameras are supported.");
2879 return;
2880 }
2881
2882 viewPtList.push_back(tmp);
2883
2884 // Now save the view point to a .txt file
2885 // FWJ DEBUG
2886 // G4cout << "saveViewPt: writing to Bookmark file " << fileName << G4endl;
2887
2888 std::string vpName = name;
2889
2890 while ((int) vpName.size() <= MAX_VP_NAME)
2891 vpName += " ";
2892
2893 fileOut << vpName << std::endl;
2894 tmp.position.getValue(x, y, z);
2895 fileOut << x << " " << y << " " << z << std::endl;
2896
2897 // Reusing x, y and z for storing the axis
2898 tmp.orientation.getValue(axis, angle);
2899 axis.getValue(x, y, z);
2900 fileOut << x << " " << y << " " << z << " " << angle << std::endl;
2901
2902 fileOut << tmp.camType << " " << tmp.height << std::endl;
2903 fileOut << tmp.focalDistance << " ";
2904 fileOut << tmp.nearDistance << " ";
2905 fileOut << tmp.farDistance << std::endl;
2906 fileOut << tmp.viewportMapping << " ";
2907 fileOut << tmp.aspectRatio << "\n" << std::endl;
2908 fileOut.flush();
2909
2910 viewPtIdx++;
2911
2912 // FWJ DEBUG
2913 // G4cout << "saveViewPt: finished writing to file" << G4endl <<
2914 // " Next viewPtIdx is " << viewPtIdx << G4endl;
2915}
2916
2917
2918// Updates the viewPtIdx in a viewpoint file.
2919
2921{
2922 std::string idxStr;
2923 std::stringstream out;
2924
2925 out << viewPtIdx;
2926 idxStr = out.str();
2927 fileOut.seekp(0, std::ios::beg);
2928
2929 while ((int) idxStr.length() < MAX_VP_IDX) {
2930 idxStr += " ";
2931 }
2932
2933 // FWJ DEBUG
2934 // G4cout << "writeViewPtIdx: " << viewPtIdx << G4endl;
2935 fileOut << idxStr << "\n";
2936 fileOut.flush();
2937 fileOut.seekp(0, std::ios::end);
2938}
2939
2940
2941// Receives the name of the bookmark clicked and searches for it in viewPtList.
2942
2943void G4OpenInventorQtExaminerViewer::LoadBookmarkCB(QListWidgetItem* item)
2944{
2945 // FWJ DEBUG
2946 // G4cout << "AuxWindow: listWidget LoadBookmark CALLBACK" << G4endl;
2947
2948 for (int i = 0; i < (int)viewPtList.size(); i++) {
2949 if (!strcmp(viewPtList[i].viewPtName, qPrintable(item->text()))) {
2950 viewPtIdx = i;
2951 break;
2952 }
2953 }
2954 // G4cout << " FOUND viewPtIdx " << viewPtIdx << G4endl;
2955
2957 setViewPt();
2958 AuxWindowDialog->lineEdit->setText(item->text());
2959}
2960
2961
2962// Sets the viewpoint based on camera data that viewPtIdx is pointing to.
2963
2965{
2967 || currentState == ROTATING) {
2968 if (animateSensor->isScheduled()) animateSensor->unschedule();
2969 setSuperimpositionEnabled(superimposition, FALSE);
2970 maxSpeed = 0.0f;
2971 scheduleRedraw();
2972 }
2973
2974 SoCamera * camera = getCamera();
2975 if (camera == NULL) {
2976 G4warn << "setViewPt: Camera is null. Unable to set the viewpoint." <<
2977 G4endl;
2978 // String dialogName = (char *) "Missing Camera Node";
2979 // std::string msg = "Camera is null. Unable to set the viewpoint.";
2980 // warningMsgDialog(msg, dialogName, NULL);
2981 return;
2982 }
2983
2984 if (!viewPtList.size()) {
2985 G4warn << "setViewPt: There are no viewpoints to load." << G4endl;
2986 // String dialogName = (char *) "Missing Viewpoints";
2987 // std::string msg = "There are no viewpoints to load.";
2988 // warningMsgDialog(msg, dialogName, NULL);
2989 return;
2990 }
2991
2992 if (SoQtExaminerViewer::isAnimating()) stopAnimating();
2993
2994 if (currentState != VIEWPOINT) {
2996 //////////////////////////////////////////////////////////////
2997 setSuperimpositionEnabled(superimposition, TRUE);
2998 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
2999 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
3000 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
3001 scheduleRedraw();
3002 ///////////////////////////////////////////////////////////////
3003 }
3004
3005 curViewPtName = viewPtList[viewPtIdx].viewPtName;
3006 camera->viewportMapping = viewPtList[viewPtIdx].viewportMapping;
3007 camera->position = viewPtList[viewPtIdx].position;
3008 camera->orientation = viewPtList[viewPtIdx].orientation;
3009 camera->aspectRatio = viewPtList[viewPtIdx].aspectRatio;
3010 camera->nearDistance = viewPtList[viewPtIdx].nearDistance;
3011 camera->farDistance = viewPtList[viewPtIdx].farDistance;
3012 camera->focalDistance = viewPtList[viewPtIdx].focalDistance;
3013
3014 // Restore camera height (changed by zooming)
3015 if (camera->isOfType(SoPerspectiveCamera::getClassTypeId())) {
3016 if (viewPtList[viewPtIdx].camType == ORTHOGRAPHIC) {
3017 toggleCameraType();
3018 camera = getCamera();
3019 ((SoOrthographicCamera *) camera)->height.setValue(
3020 viewPtList[viewPtIdx].height);
3021 } else
3022 ((SoPerspectiveCamera *) camera)->heightAngle.setValue(
3023 viewPtList[viewPtIdx].height);
3024 } else if (camera->isOfType(SoOrthographicCamera::getClassTypeId())) {
3025 if (viewPtList[viewPtIdx].camType == PERSPECTIVE) {
3026 toggleCameraType();
3027 camera = getCamera();
3028 ((SoPerspectiveCamera *) camera)->heightAngle.setValue(
3029 viewPtList[viewPtIdx].height);
3030 } else
3031 ((SoOrthographicCamera *) camera)->height.setValue(
3032 viewPtList[viewPtIdx].height);
3033 } else {
3034 SoDebugError::post("G4OpenInventorQtExaminerViewer::setViewPt",
3035 "Only Perspective and Orthographic cameras are supported.");
3036 return;
3037 }
3038
3039}
3040
3041
3042void G4OpenInventorQtExaminerViewer::NextViewPtCB()
3043{
3044 // FWJ DEBUG
3045 // G4cout << "App Button: nextViewPt CALLBACK" << G4endl;
3046
3047 if (!viewPtList.size()) return;
3048 if (viewPtIdx >= (int)viewPtList.size() - 1)
3049 viewPtIdx = 0;
3050 else
3051 viewPtIdx++;
3052
3054 setViewPt();
3055 char* viewptname = viewPtList[viewPtIdx].viewPtName;
3056 AuxWindowDialog->lineEdit->setText(QString(viewptname));
3057}
3058
3059void G4OpenInventorQtExaminerViewer::PrevViewPtCB()
3060{
3061 // FWJ DEBUG
3062 // G4cout << "App Button: prevViewPt CALLBACK" << G4endl;
3063
3064 if (!viewPtList.size()) return;
3065 if (viewPtIdx == 0)
3066 viewPtIdx = (int) viewPtList.size() - 1;
3067 else
3068 viewPtIdx--;
3069
3071 setViewPt();
3072 char* viewptname = viewPtList[viewPtIdx].viewPtName;
3073 AuxWindowDialog->lineEdit->setText(QString(viewptname));
3074}
3075
3076
3077void G4OpenInventorQtExaminerViewer::AbbrOutputCB(bool checked)
3078{
3079 // FWJ DEBUG
3080 // G4cout << "App Button: abbrOutput CALLBACK" << G4endl;
3081
3082 abbrOutputFlag = checked;
3083}
3084
3085
3086void G4OpenInventorQtExaminerViewer::PickRefPathCB()
3087{
3088 // FWJ DEBUG
3089 // G4cout << "App Button: pickRefPath CALLBACK" << G4endl;
3090
3091 // Save viewing state and go to picking mode
3092 viewingBeforePickRef = isViewing();
3093 if(isViewing())
3094 setViewing(false);
3095 setComponentCursor(SoQtCursor(SoQtCursor::CROSSHAIR));
3096 pickRefPathFlag = true;
3097}
3098
3099
3100void G4OpenInventorQtExaminerViewer::SwitchWireFrameCB(bool checked)
3101{
3102 // FWJ DEBUG
3103 // G4cout << "App Button: switchWireFrame CALLBACK" << G4endl;
3104
3105 // if (switchWireFrameButton->isDown()) {
3106 if (checked) {
3107 setDrawStyle(SoQtViewer::STILL, SoQtViewer::VIEW_LINE);
3108 setDrawStyle(SoQtViewer::INTERACTIVE, SoQtViewer::VIEW_LINE);
3109 } else {
3110 setDrawStyle(SoQtViewer::STILL, SoQtViewer::VIEW_AS_IS);
3111 setDrawStyle(SoQtViewer::INTERACTIVE,
3112 SoQtViewer::VIEW_SAME_AS_STILL);
3113 }
3114}
3115
3116
3117void G4OpenInventorQtExaminerViewer::SwitchAxesCB(bool checked)
3118{
3119 // FWJ DEBUG
3120 // G4cout << "App Button: switchAxes CALLBACK" << G4endl;
3121 setFeedbackVisibility(checked);
3122 // if (checked) {
3123 // setFeedbackVisibility(TRUE);
3124 // } else {
3125 // setFeedbackVisibility(FALSE);
3126 // }
3127}
3128
3129
3130void G4OpenInventorQtExaminerViewer::DetachCB()
3131{
3132 // FWJ DEBUG
3133 // G4cout << "App Button: detach CALLBACK" << G4endl;
3134 uiQt->GetViewerTabWidget()->removeTab(uiQtTabIndex);
3135 viewerParent->setParent(viewerParent2);
3136 removeAppPushButton(detachButton);
3137 show();
3138}
3139
3140
3141void G4OpenInventorQtExaminerViewer::DeleteBookmarkCB()
3142{
3143 // FWJ DEBUG
3144 // G4cout << "Delete Button: DeleteBookmark CALLBACK" << G4endl;
3145
3146 // Maybe nothing selected yet
3147 QListWidgetItem* listitem = AuxWindowDialog->listWidget->currentItem();
3148 if (!listitem) return;
3149 if (!(listitem->isSelected())) return;
3150
3151 QString vpnamein = listitem->text();
3152
3153 char* vpName = strdup(qPrintable(vpnamein));
3154 // G4cout << "DELETING bookmark " << vpName << G4endl;
3155
3156 deleteViewPt(vpName);
3157 delete listitem;
3158 free(vpName);
3159}
3160
3161// Deletes current viewpoint the user is looking at.
3162// Updates the input file and bookmarks as well.
3163
3165{
3166 std::string line;
3167 std::size_t end;
3168 fileIn.open(fileName.c_str());
3169 std::ofstream out("temporaryFile.txt");
3170
3171 if (!vpName)
3172 vpName = viewPtList[viewPtIdx].viewPtName;
3173
3174 getline(fileIn, line); // Printing the viewpoint idx
3175 out << line << "\n";
3176
3177 while (getline(fileIn, line)) {
3178 end = line.find_last_not_of(' ');
3179 line = line.substr(0, end + 1);
3180 if (!strcmp(line.c_str(), vpName)) { // Equal
3181 while (line.size()) {
3182 getline(fileIn, line);
3183 }
3184
3185 while (getline(fileIn, line))
3186 out << line << "\n";
3187 } else {
3188 while (line.size()) {
3189 out << line << "\n";
3190 getline(fileIn, line);
3191 }
3192 out << "\n";
3193 }
3194 }
3195
3196 std::size_t idx = 0; // Remove viewpoint from the vector
3197 std::size_t size = viewPtList.size();
3198 while (idx < size) {
3199 if (!strcmp(viewPtList[idx].viewPtName, vpName)) {
3200 viewPtList.erase(viewPtList.begin() + idx);
3201 break;
3202 }
3203 idx++;
3204 }
3205
3206 out.close();
3207 fileOut.close();
3208 fileIn.clear();
3209 fileIn.close();
3210
3211 // FWJ check return status: error popups needed here
3212 int istat = remove(fileName.c_str());
3213 if (istat == -1) {
3214 QMessageBox msgbox;
3215 msgbox.setFont(*font);
3216 QString messagetxt = "Error removing bookmarks file";
3217 // messagetxt.append(filenamein);
3218 msgbox.setText(messagetxt);
3219 msgbox.exec();
3220 // G4cout << "Error removing bookmarks file" << G4endl;
3221 }
3222 istat = rename("temporaryFile.txt", fileName.c_str());
3223 if (istat == -1) {
3224 QMessageBox msgbox;
3225 msgbox.setFont(*font);
3226 QString messagetxt = "Error renaming bookmarks file";
3227 // messagetxt.append(filenamein);
3228 msgbox.setText(messagetxt);
3229 msgbox.exec();
3230 // G4cout << "Error renaming bookmarks file" << G4endl;
3231 }
3232 fileOut.open(fileName.c_str(), std::ios::in);
3233 fileOut.seekp(0, std::ios::end);
3234
3235 if (!viewPtList.size()) { // viewPtList is empty
3236 curViewPtName = (char *) empty.c_str();
3237 scheduleRedraw();
3238 } else {
3239 if (viewPtIdx >= (int) viewPtList.size())
3240 viewPtIdx--;
3242 setViewPt();
3243 }
3244}
3245
3246
3247void G4OpenInventorQtExaminerViewer::RenameBookmarkCB()
3248{
3249 // FWJ DEBUG
3250 // G4cout << "Rename Button: RenameBookmark CALLBACK" << G4endl;
3251 // Maybe nothing selected yet
3252 QListWidgetItem* listitem = AuxWindowDialog->listWidget->currentItem();
3253 if (!listitem) return;
3254 if (!(listitem->isSelected())) return;
3255
3256 QString vpnamein = listitem->text();
3257
3258 QInputDialog* inputdialog = new QInputDialog(getParentWidget());
3259 inputdialog->setFont(*font);
3260 inputdialog->setWindowTitle(tr("Enter"));
3261 inputdialog->setLabelText("New bookmark name");
3262 inputdialog->adjustSize();
3263 QString newnamein;
3264 if (inputdialog->exec() == QDialog::Accepted)
3265 newnamein=inputdialog->textValue().trimmed();
3266 else
3267 return;
3268 if (newnamein.isEmpty()) return;
3269
3270 char* newname = strdup(qPrintable(newnamein));
3271
3272 std::size_t size = viewPtList.size();
3273 for (std::size_t i = 0; i < size; ++i) {
3274 if (!strcmp(newname, viewPtList[i].viewPtName)) {
3275 QMessageBox msgbox;
3276 msgbox.setFont(*font);
3277 msgbox.setText("Bookmark name is already in use");
3278 msgbox.exec();
3279 }
3280 }
3281
3282 // G4cout << "RENAMING to " << newname << G4endl;
3283 renameViewPt(newname);
3284 listitem->setText(QString(newname));
3285 AuxWindowDialog->lineEdit->setText(newname);
3286 // if (currentState == VIEWPOINT)
3287 // scheduleRedraw();
3288
3289 free(newname);
3290}
3291
3292// Renames currently selected viewpoint.
3293
3295{
3296 std::size_t idx = 0, end, pos;
3297 std::size_t size = viewPtList.size();
3298 std::string line, newName;
3299 fileIn.open(fileName.c_str());
3300
3301 newName = vpName;
3302 while ((int) newName.size() < MAX_VP_NAME)
3303 newName += " ";
3304
3305 getline(fileIn, line);
3306 pos = fileIn.tellg();
3307 while (getline(fileIn, line)) {
3308 end = line.find_last_not_of(' ');
3309 line = line.substr(0, end + 1);
3310 if (!strcmp(line.c_str(), curViewPtName)) {
3311 fileOut.seekp(pos);
3312 fileOut << newName;
3313 fileOut.seekp(0, std::ios::end); // Set the file pointer to the end of the file
3314 break;
3315 }
3316 while (line.size())
3317 getline(fileIn, line);
3318 pos = fileIn.tellg();
3319 }
3320
3321 fileIn.close();
3322 fileIn.clear();
3323
3324 while (idx < size) {
3325 if (!strcmp(viewPtList[idx].viewPtName, curViewPtName)) {
3326 strcpy(viewPtList[idx].viewPtName, vpName);
3327 break;
3328 }
3329 idx++;
3330 }
3331}
3332
3333
3334void G4OpenInventorQtExaminerViewer::SortBookmarksCB()
3335{
3336 // FWJ NOTE: Xt version of this does not work (does nothing)
3337
3338 // G4cout << "Sort Button: SortBookmarks CALLBACK" << G4endl;
3339
3340 // FWJ List for sorting
3341 // The dialog list and bookmark file will be rewritten.
3342 // Simpler to populate this list from the data structure.
3343
3344 std::vector<std::string> charList;
3345
3346 if (viewPtList.size() < 2) return;
3347
3348 // Get current entries from the list
3349
3350 for (int i = 0; i < (int)viewPtList.size(); i++) {
3351
3352 charList.push_back(viewPtList[i].viewPtName);
3353 // G4cout << " Pushed " << i << " " << charList[i] << G4endl;
3354 }
3355
3356 std::sort(charList.begin(), charList.end());
3357
3358 // FWJ POPULATE the new dialog list
3359 // G4cout << " Populating Bookmark listWidget..." << G4endl;
3360 AuxWindowDialog->listWidget->clear();
3361
3362 for (int i = 0; i < (int)viewPtList.size(); i++) {
3363 // viewPtIdx has to be changed to account for a different order in viewPtList
3364 if (!strcmp(charList[i].c_str(), curViewPtName))
3365 viewPtIdx = i;
3366 new QListWidgetItem(charList[i].c_str(), AuxWindowDialog->listWidget);
3367
3368 }
3369
3370 sortViewPts(charList);
3371
3372}
3373
3374// Rewrites entire viewpoint file with sorted viewpoints.
3375
3376void G4OpenInventorQtExaminerViewer::sortViewPts(std::vector<std::string> sortedViewPts)
3377{
3378 SbVec3f axis;
3379 float x, y, z, angle;
3380 std::size_t sortIdx = 0, unsortIdx = 0;
3381
3382 if (fileOut.is_open())
3383 fileOut.close();
3384
3385 fileOut.open(fileName.c_str()); // Erase current viewpoint file
3386
3388
3389 std::size_t size = sortedViewPts.size();
3390 while (sortIdx < size) {
3391 while (strcmp(sortedViewPts[sortIdx].c_str(),
3392 viewPtList[unsortIdx].viewPtName))
3393 unsortIdx++;
3394
3395 std::string vpName = viewPtList[unsortIdx].viewPtName;
3396
3397 while ((int) vpName.size() < MAX_VP_NAME)
3398 vpName += " ";
3399 fileOut << vpName << std::endl;
3400 viewPtList[unsortIdx].position.getValue(x, y, z);
3401 fileOut << x << " " << y << " " << z << std::endl;
3402
3403 // Reusing x, y and z for storing the axis
3404 viewPtList[unsortIdx].orientation.getValue(axis, angle);
3405 axis.getValue(x, y, z);
3406 fileOut << x << " " << y << " " << z << " " << angle << std::endl;
3407
3408 fileOut << viewPtList[unsortIdx].camType << " "
3409 << viewPtList[unsortIdx].height << std::endl;
3410 fileOut << viewPtList[unsortIdx].focalDistance << " ";
3411
3412 fileOut << viewPtList[unsortIdx].nearDistance << " ";
3413
3414 fileOut << viewPtList[unsortIdx].farDistance << std::endl;
3415
3416 fileOut << viewPtList[unsortIdx].viewportMapping << " ";
3417 fileOut << viewPtList[unsortIdx].aspectRatio << "\n" << std::endl;
3418 fileOut.flush();
3419
3420 unsortIdx = 0;
3421 sortIdx++;
3422 }
3423}
3424
3425// Needed to implement mouse wheel zoom direction change.
3426// Does not work with MacOS trackpad: use Coin3d default handler.
3427// Emulating private method SoGuiFullViewerP::zoom()
3428#ifndef __APPLE__
3429void
3431{
3432 float multiplicator = float(std::exp(diffvalue));
3433 SoCamera *cam = getCamera();
3434
3435 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
3436 const float oldfocaldist = cam->focalDistance.getValue();
3437 const float newfocaldist = oldfocaldist * multiplicator;
3438
3439 SbVec3f direction;
3440 cam->orientation.getValue().multVec(SbVec3f(0, 0, -1), direction);
3441
3442 const SbVec3f oldpos = cam->position.getValue();
3443 const SbVec3f newpos = oldpos + (newfocaldist - oldfocaldist) * -direction;
3444 cam->position = newpos;
3445 cam->focalDistance = newfocaldist;
3446 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
3447 SoOrthographicCamera * oc = (SoOrthographicCamera *)cam;
3448 oc->height = oc->height.getValue() * multiplicator;
3449 }
3450}
3451#endif
3452
3453// Handling mouse and keyboard events
3454
3455SbBool
3457{
3458
3459 // FWJ DEBUG
3460 // G4cout << "processSoEvent ############" << ++processSoEventCount << G4endl;
3461
3462 SoCamera *cam = getCamera();
3463 const SoType type(ev->getTypeId());
3464
3465// Needed to implement mouse wheel zoom direction change.
3466// Does not work with MacOS trackpad: use Coin3d default handler.
3467#ifndef __APPLE__
3468 if (type.isDerivedFrom(SoMouseButtonEvent::getClassTypeId())) {
3469 SoMouseButtonEvent * me = (SoMouseButtonEvent *) ev;
3470
3471 // if (currentState == ANIMATION || currentState == REVERSED_ANIMATION
3472 // || currentState == PAUSED_ANIMATION) {
3473
3474 switch (me->getButton()) {
3475
3476 case SoMouseButtonEvent::BUTTON4: // Scroll wheel up
3477 if (me->getState() == SoButtonEvent::DOWN) {
3478 // G4cout << "SCROLL WHEEL UP" << G4endl;
3479 zoom(-0.1f);
3480 return TRUE;
3481 }
3482 break;
3483
3484 case SoMouseButtonEvent::BUTTON5: // Scroll wheel down
3485 if (me->getState() == SoButtonEvent::DOWN) {
3486 // G4cout << "SCROLL WHEEL DOWN" << G4endl;
3487 zoom(0.1f);
3488 return TRUE;
3489 }
3490 break;
3491
3492 default:
3493 break;
3494 }
3495 // }
3496 if (currentState == GENERAL) {
3497
3498 }
3499 }
3500#endif
3501
3502 if (type.isDerivedFrom(SoKeyboardEvent::getClassTypeId())) {
3503 SoKeyboardEvent* ke = (SoKeyboardEvent*)ev;
3504
3505 if (SoKeyboardEvent::isKeyPressEvent(ev, ke->getKey())) {
3506 switch (ke->getKey()) {
3507 case SoKeyboardEvent::E:
3508 if (externalQtApp) {
3509 // G4cout << "E KEY PRESSED" << G4endl;
3510 return TRUE;
3511 } else {
3512 G4cout <<
3513 "E KEY PRESSED, EXITING OIQT VIEWER SECONDARY LOOP" <<
3514 G4endl;
3515 SoQt::exitMainLoop();
3516 // escapeCallback();
3517 return TRUE;
3518 }
3519 case SoKeyboardEvent::LEFT_SHIFT:
3520 this->lshiftdown = true;
3521 return TRUE;
3522 case SoKeyboardEvent::RIGHT_SHIFT:
3523 this->rshiftdown = true;
3524 return TRUE;
3525 case SoKeyboardEvent::LEFT_CONTROL:
3526 this->lctrldown = true;
3527 return TRUE;
3528 case SoKeyboardEvent::RIGHT_CONTROL:
3529 this->rctrldown = true;
3530 return TRUE;
3531 case SoKeyboardEvent::SPACE:
3532 if (currentState == ANIMATION
3536 if (animateSensor->isScheduled())
3537 animateSensor->unschedule();
3538 return TRUE;
3539 } else if (currentState == PAUSED_ANIMATION) {
3540 if (maxSpeed) {
3541 if ((beforePausing == ANIMATION
3543 < (int) refParticleTrajectory.size() - 1)
3545 && refParticleIdx > 0)) {
3548 }
3549 }
3550 return TRUE;
3551 }
3552 break;
3553 case SoKeyboardEvent::ESCAPE:
3554 if (currentState == ANIMATION
3557
3558 if (animateSensor->isScheduled())
3559 animateSensor->unschedule();
3562 setSuperimpositionEnabled(superimposition, FALSE);
3563 maxSpeed = 0.0f;
3564 step = 1;
3565
3566 scheduleRedraw();
3567 if (currentState == VIEWPOINT) {
3568 setSuperimpositionEnabled(superimposition, TRUE);
3569 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
3570 animSpeedOutlineSwitch->whichChild.setValue(
3571 SO_SWITCH_NONE);
3572 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
3573
3574 scheduleRedraw();
3575 }
3576 restoreCamera();
3577 return TRUE;
3578 }
3579 break;
3580 case SoKeyboardEvent::DELETE:
3581 if (viewPtList.size()
3582 && (currentState != ANIMATION
3585 // FWJ IMPLEMENT LATER
3586 // String dialogName = (char *) "Delete Viewpoint";
3587 // std::string msg = "Are you sure you want to delete current viewpoint?";
3588 // warningMsgDialog(msg, dialogName, deleteViewPtCB);
3589 return TRUE;
3590 }
3591 break;
3592 case SoKeyboardEvent::LEFT_ARROW:
3593 switch (currentState) {
3594 case BEAMLINE:
3595 if ((this->lshiftdown) || (this->rshiftdown)) {
3597 moveCamera();
3598 }
3599 else if ((this->lctrldown) || (this->rctrldown)) {
3600 if (SoQtExaminerViewer::isAnimating())
3601 stopAnimating();
3604 animateBtwPtsPeriod = 0.08f;
3605
3606 SbVec3f tmp = camDir;
3607 tmp.negate();
3608 rotAxis = tmp;
3609
3610 rotCnt = ROT_CNT;
3611 moveCamera(); // To make sure camera is perpendicular to the beamline
3612 rotateCamera();
3613 }
3614 else {
3615 if (SoQtExaminerViewer::isAnimating())
3616 stopAnimating();
3619 animateBtwPtsPeriod = 0.08f;
3620
3621 SbVec3f tmp = camUpVec;
3622 tmp.negate();
3623 rotAxis = tmp;
3624
3625 rotCnt = ROT_CNT;
3626 moveCamera(); // To make sure camera is perpendicular to the beamline
3627 rotateCamera();
3628
3629 }
3630 return TRUE;
3631
3632 case ANIMATION:
3633 case REVERSED_ANIMATION:
3634 left_right -= 1.5f;
3635 return TRUE;
3636 case PAUSED_ANIMATION:
3637 left_right -= 1.5f;
3639 cam->position = myCam->position;
3640 return TRUE;
3641 case GENERAL:
3642 case VIEWPOINT:
3643 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3644 // Using this allows us to look around without
3645 // changing the camera parameters (camDir, camUpVec)
3646 this->bottomWheelMotion(
3647 this->getBottomWheelValue() + 0.1f);
3648
3649 return TRUE;
3650 }
3651 break;
3652 case ROTATING:
3653 // For this state, let the keyboard event
3654 // be handled by superclass
3655 break;
3656 default:
3657 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3658 "Unhandled viewer state");
3659 break;
3660 }
3661 break;
3662
3663 case SoKeyboardEvent::RIGHT_ARROW:
3664 switch(currentState) {
3665 case BEAMLINE:
3666 if ((this->lshiftdown) || (this->rshiftdown)) {
3668 moveCamera();
3669 }
3670 else if ((this->lctrldown) || (this->rctrldown)) {
3671 if (SoQtExaminerViewer::isAnimating())
3672 stopAnimating();
3675 animateBtwPtsPeriod = 0.08f;
3676
3677 rotAxis = camDir;
3678
3679 rotCnt = ROT_CNT;
3680 moveCamera(); // To make sure camera is perpendicular to the beamline
3681 rotateCamera();
3682 }
3683 else{
3684 if (SoQtExaminerViewer::isAnimating())
3685 stopAnimating();
3688 animateBtwPtsPeriod = 0.08f;
3689
3690 rotAxis = camUpVec;
3691
3692 rotCnt = ROT_CNT;
3693 moveCamera(); // To make sure camera is perpendicular to the beamline
3694 rotateCamera();
3695 }
3696 return TRUE;
3697
3698 case ANIMATION:
3699 case REVERSED_ANIMATION:
3700 left_right += 1.5f;
3701 return TRUE;
3702 case PAUSED_ANIMATION:
3703 left_right += 1.5f;
3705 cam->position = myCam->position;
3706 return TRUE;
3707 case GENERAL:
3708 case VIEWPOINT:
3709 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3710 // Using this allows us to look around without
3711 // changing the camera parameters (camDir, camUpVec)
3712 this->bottomWheelMotion(
3713 this->getBottomWheelValue() - 0.1f);
3714 return TRUE;
3715 }
3716 break;
3717 case ROTATING:
3718 // For this state, let the keyboard event
3719 // be handled by superclass
3720 break;
3721 default:
3722 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3723 "Unhandled viewer state");
3724 break;
3725 }
3726 break;
3727
3728 case SoKeyboardEvent::DOWN_ARROW:
3729 switch(currentState) {
3730 case BEAMLINE:
3731
3732 if ((this->lshiftdown) || (this->rshiftdown)) {
3734 moveCamera();
3735 }
3736 else{
3737 if (SoQtExaminerViewer::isAnimating())
3738 stopAnimating();
3741 animateBtwPtsPeriod = 0.08f;
3742
3743 rotAxis = camDir.cross(camUpVec);
3744
3745 rotCnt = ROT_CNT;
3746 moveCamera(); // To make sure camera is perpendicular to the beamline
3747 rotateCamera();
3748
3749 }
3750 return TRUE;
3751
3752 case ANIMATION:
3753 case REVERSED_ANIMATION:
3754 up_down -= 1.5f;
3755 return TRUE;
3756 case PAUSED_ANIMATION:
3757 up_down -= 1.5f;
3759 cam->position = myCam->position;
3760 return TRUE;
3761 case GENERAL:
3762 case VIEWPOINT:
3763 // Using this allows us to look around without
3764 // changing the camera parameters (camDir, camUpVec)
3765 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3766 this->leftWheelMotion(this->getLeftWheelValue() - 0.1f);
3767 return TRUE;
3768 }
3769 break;
3770 case ROTATING:
3771 // For this state, let the keyboard event
3772 // be handled by superclass
3773 break;
3774 default:
3775 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3776 "Unhandled viewer state");
3777 break;
3778 }
3779 break;
3780
3781 case SoKeyboardEvent::UP_ARROW:
3782 switch(currentState) {
3783 case BEAMLINE:
3784 if ((this->lshiftdown) || (this->rshiftdown)) {
3786 moveCamera();
3787 }
3788 else{
3789 if (SoQtExaminerViewer::isAnimating())
3790 stopAnimating();
3793 animateBtwPtsPeriod = 0.08f;
3794
3795 rotAxis = camUpVec.cross(camDir);
3796
3797 rotCnt = ROT_CNT;
3798 moveCamera();
3799
3800 rotateCamera();
3801
3802
3803 }
3804 return TRUE;
3805 case ANIMATION:
3806 case REVERSED_ANIMATION:
3807 up_down += 1.5f;
3808 return TRUE;
3809 case PAUSED_ANIMATION:
3810 up_down += 1.5f;
3812 cam->position = myCam->position;
3813 return TRUE;
3814 case GENERAL:
3815 case VIEWPOINT:
3816 // Using this allows us to look around without
3817 // changing the camera parameters (camDir, camUpVec)
3818 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3819 this->leftWheelMotion(this->getLeftWheelValue() + 0.1f);
3820 return TRUE;
3821 }
3822 break;
3823 case ROTATING:
3824 // For this state, let the keyboard event
3825 // be handled by superclass
3826 break;
3827 default:
3828 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3829 "Unhandled viewer state");
3830 break;
3831 }
3832 break;
3833
3834 case SoKeyboardEvent::PAGE_UP:
3835 switch(currentState) {
3836 case BEAMLINE:
3837 if (step < (int) refParticleTrajectory.size() / 5) // Magic number
3838 step++;
3839 return TRUE;
3840 case ANIMATION:
3841 incSpeed();
3843 if (maxSpeed > 0.8)
3845 scheduleRedraw();
3846
3847 return TRUE;
3848 case REVERSED_ANIMATION:
3849 if(!animateSensor->isScheduled()) {
3851 if (refParticleIdx
3852 < (int) refParticleTrajectory.size() - 1) {
3855 scheduleRedraw();
3857 }
3858 }
3859 else{
3861 decSpeed();
3862 scheduleRedraw();
3863 }
3864 return TRUE;
3865 case PAUSED_ANIMATION:
3867 if (maxSpeed > 0.8)
3869
3870 if (beforePausing == ANIMATION) {
3871 incSpeed();
3872 } else {
3873 decSpeed();
3876 }
3877
3878 scheduleRedraw();
3879 return TRUE;
3880 default: //fall through
3881 break;
3882 }
3883 break;
3884
3885 case SoKeyboardEvent::PAGE_DOWN:
3886 switch(currentState) {
3887 case BEAMLINE:
3888 if (step > 1)
3889 step--;
3890 return TRUE;
3891 case ANIMATION:
3892 if(!animateSensor->isScheduled()) {
3894 if (refParticleIdx > 1) {
3897 scheduleRedraw();
3899 }
3900 }
3901 else{
3903 decSpeed();
3904 scheduleRedraw();
3905 }
3906 return TRUE;
3907 case REVERSED_ANIMATION:
3908 incSpeed();
3910 if (maxSpeed < -0.8)
3912 scheduleRedraw();
3913 return TRUE;
3914 case PAUSED_ANIMATION:
3916 if (maxSpeed < -0.8)
3919 incSpeed();
3920 } else {
3921 decSpeed();
3924 }
3925 scheduleRedraw();
3926 return TRUE;
3927 default:
3928 //fall through
3929 break;
3930 }
3931 break;
3932
3933 // FROM XT VIEWER
3934 // case SoKeyboardEvent::E:
3935 // this->escapeCallback(this->examinerObject);
3936 // break;
3937
3938 default:
3939 break; // To get rid of compiler warnings
3940 }
3941 }
3942 if (SoKeyboardEvent::isKeyReleaseEvent(ev, ke->getKey())) {
3943 switch (ke->getKey()) {
3944 case SoKeyboardEvent::LEFT_SHIFT:
3945 this->lshiftdown = false;
3946 return TRUE;
3947 case SoKeyboardEvent::RIGHT_SHIFT:
3948 this->rshiftdown = false;
3949 return TRUE;
3950 case SoKeyboardEvent::LEFT_CONTROL:
3951 this->lctrldown = false;
3952 return TRUE;
3953 case SoKeyboardEvent::RIGHT_CONTROL:
3954 this->rctrldown = false;
3955 return TRUE;
3956 default:
3957 break;
3958 }
3959 }
3960 }
3961
3962 // Pass the event on to the viewer
3963 // Need some checks here as in Xt viewer?
3964
3966 || currentState == ROTATING)
3967 return FALSE;
3968 else
3969 return SoQtExaminerViewer::processSoEvent(ev);
3970}
3971
3972
3973// REMAINDER OF MENU BAR FUNCTIONS...
3974
3975
3976void G4OpenInventorQtExaminerViewer::FileLoadSceneGraphCB()
3977{
3978 // G4cout << "File: Load scene graph CALLBACK" << G4endl;
3979
3980 QFileDialog filedialog(getParentWidget(), tr("Load Scene Graph"));
3981 filedialog.setFileMode(QFileDialog::AnyFile);
3982 filedialog.setFont(*font);
3983 if (!filedialog.exec()) return;
3984 QStringList filenameinlist = filedialog.selectedFiles();
3985 QString filenamein = filenameinlist[0];
3986
3987 SoInput sceneInput;
3988
3989 if (sceneInput.openFile(qPrintable(filenamein))) {
3990 // Read the whole file into the database
3991 newSceneGraph = SoDB::readAll(&sceneInput);
3992 if (newSceneGraph == NULL) {
3993 QMessageBox msgbox;
3994 msgbox.setFont(*font);
3995 QString messagetxt = "Error reading scene graph file ";
3996 messagetxt.append(filenamein);
3997 msgbox.setText(messagetxt);
3998 msgbox.exec();
3999 sceneInput.closeFile();
4000 return;
4001 }
4002 } else {
4003 QMessageBox msgbox;
4004 msgbox.setFont(*font);
4005 QString messagetxt = "Error opening scene graph file ";
4006 messagetxt.append(filenamein);
4007 msgbox.setText(messagetxt);
4008 msgbox.exec();
4009 return;
4010 }
4011
4012 SoSeparator* root = (SoSeparator*)getSceneGraph();
4013 root->unref();
4014 newSceneGraph->ref();
4015 setSceneGraph(newSceneGraph);
4016}
4017
4018void G4OpenInventorQtExaminerViewer::FileSaveSceneGraphCB()
4019{
4020 // G4cout << "File: Save scene graph CALLBACK" << G4endl;
4021
4022 QFileDialog filedialog(getParentWidget(), tr("Save scene graph"));
4023 filedialog.setFileMode(QFileDialog::AnyFile);
4024 // To enable confirmation of overwriting
4025 filedialog.setAcceptMode(QFileDialog::AcceptSave);
4026 filedialog.setFont(*font);
4027 if (!filedialog.exec()) return;
4028 QStringList filenameinlist = filedialog.selectedFiles();
4029 QString filenamein = filenameinlist[0];
4030
4031 SoWriteAction writeAction;
4032 SoSeparator* root = (SoSeparator*)getSceneGraph();
4033
4034 SoOutput* out = writeAction.getOutput();
4035
4036 if (out->openFile(qPrintable(filenamein))) {
4037 out->setBinary(FALSE);
4038 writeAction.apply(root);
4039 out->closeFile();
4040 } else {
4041 QMessageBox msgbox;
4042 msgbox.setFont(*font);
4043 QString messagetxt = "Error opening file ";
4044 messagetxt.append(filenamein);
4045 msgbox.setText(messagetxt);
4046 msgbox.exec();
4047 }
4048}
4049
4050
4051void G4OpenInventorQtExaminerViewer::HelpControlsCB()
4052{
4053 // G4cout << "Help: Help Controls CALLBACK" << G4endl;
4054 helpmsgbox->show();
4055}
4056
4057
4062
4065
4067{
4068 if (requestedState == G4State_EventProc) viewer->newEvents = true;
4069 return true;
4070}
G4ApplicationState
@ G4State_EventProc
#define MAX_SPEED_INDICATOR
#define SPEED_INDICATOR_STEP
#define G4warn
Definition G4Scene.cc:41
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define M_PI
Definition SbMath.h:33
const std::vector< const std::vector< G4AttValue > * > & GetAttValues() const
const std::vector< const std::map< G4String, G4AttDef > * > & GetAttDefs() const
static void mouseoverCB(void *aThis, SoEventCallback *eventCB)
SoCoordinate3 * getCoordsNode(SoFullPath *path)
static void superimpositionCB(void *closure, SoAction *action)
static void animateSensorRotationCB(void *, SoSensor *)
void moveCamera(float dist=0, bool lookdown=false)
SoNode * getSuperimpositionNode(SoNode *, const char *name)
void parseString(T &t, const std::string &s, bool &error)
static void pickingCB(void *aThis, SoEventCallback *eventCB)
static void sceneChangeCB(void *, SoSensor *)
SbBool processSoEvent(const SoEvent *const event)
G4OpenInventorQtExaminerViewer(QWidget *parent=NULL, const char *name=NULL, SbBool embed=TRUE, SoQtFullViewer::BuildFlag flag=BUILD_ALL, SoQtViewer::Type type=BROWSER)
void sortViewPts(std::vector< std::string >)
void setReferencePath(SoLineSet *, SoCoordinate3 *, bool append=false)
void distanceToTrajectory(const SbVec3f &, float &, SbVec3f &, int &)
static void animateSensorCB(void *, SoSensor *)
G4bool AddTabWidget(QWidget *, QString)
Definition G4UIQt.cc:2283
QTabWidget * GetViewerTabWidget()
Definition G4UIQt.hh:171
static G4UImanager * GetUIpointer()
virtual G4bool Notify(G4ApplicationState requestedState)
HookEventProcState(G4OpenInventorQtExaminerViewer *)
QPushButton * pushButton_2
QPushButton * pushButton_3
void setupUi(QDialog *Dialog)
QLineEdit * lineEdit
QListWidget * listWidget
QListWidget * listWidget1
QPushButton * pushButton
#define TRUE
Definition globals.hh:41
#define FALSE
Definition globals.hh:38
const char * name(G4int ptype)