Geant4 11.1.1
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, y, z, angle;
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;
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; // 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;
1481 if (refParticleIdx < lastIdx && !animateSensor->isScheduled())
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 for (std::size_t j = 0; j < attHolder->GetAttDefs().size(); ++j) {
1609 std::ostringstream oss;
1610 oss << G4AttCheck(attHolder->GetAttValues()[j],
1611 attHolder->GetAttDefs()[j]);
1612
1613 std::string findStr = "Type of trajectory (Type): ";
1614 std::string compareValue = "REFERENCE";
1615 std::size_t idx = oss.str().find(findStr);
1616
1617 if(idx != std::string::npos) {
1618 if(oss.str().substr(idx + findStr.size(),
1619 compareValue.size()) == compareValue) {
1620 coords = getCoordsNode(path);
1621 if(coords != NULL) {
1622 refPathFound = true;
1623 coordvec.push_back(coords);
1624 linevec.push_back((SoLineSet *)path->getTail());
1625 }
1626 break;
1627 }
1628 }
1629
1630 findStr = "Track ID (ID): ";
1631 idx = oss.str().find(findStr);
1632 if(idx != std::string::npos) {
1633 //index all primary tracks
1634 std::string tmpstr = oss.str().substr(idx + findStr.size(),1);
1635 std::istringstream buffer(tmpstr);
1636 int num;
1637 buffer >> num;
1638 if(num == 1) {
1639
1640 // Check if next character is a number,
1641 // in which case we don't have Track ID 1
1642 // FWJ attempt to fix Coverity issue.
1643 char nextChar = oss.str().at(idx+findStr.size()+1);
1644 // const char * nextChar =
1645 // oss.str().substr(idx + findStr.size() + 1,1).c_str();
1646 if(std::isdigit(nextChar))
1647 break; //Not a primary track, continue with next track
1648
1649 coords = getCoordsNode(path);
1650 if(coords != NULL) {
1651 coordvec.push_back(coords);
1652 linevec.push_back((SoLineSet *)path->getTail());
1653 break; //Found coords node, continue with next track
1654 }
1655 }
1656 else
1657 break; //Not a primary track, continue with next track
1658 }
1659 else{
1660 //Not a Track ID attribute, fall through
1661 }
1662 }
1663
1664 if(refPathFound)
1665 break;
1666 }
1667
1668 if(coordvec.empty())
1669 return; //No track with a Coordinate3 node found
1670
1671 if(refPathFound) {
1672 //set ref path to last traj, coord in the vecs
1673 setReferencePath(linevec.back(), coordvec.back());
1674 return;
1675 }
1676 //else
1677
1678 int longestIdx = 0;
1679 float longestLength = 0.0;
1680 // For all paths
1681 for(unsigned int i=0;i < linevec.size(); ++i) {
1682
1683 //First generate a vector with all the points in this lineset
1684 std::vector<SbVec3f> trajectory;
1685 // For all lines in the i path
1686 for(int j=0; j < linevec[i]->numVertices.getNum(); ++j) {
1687 // For all points in line j
1688 for(int k=0; k < linevec[i]->numVertices[j]; ++k) {
1689 trajectory.push_back(coordvec[i]->point[k]);
1690 }
1691 }
1692
1693 // Then calculate the total length
1694 float tmpLength=0.0;
1695 for(unsigned int j=0; j < trajectory.size() - 1; ++j) {
1696 tmpLength += (trajectory[j] - trajectory[j + 1]).length();
1697 }
1698
1699 if(tmpLength > longestLength) {
1700 longestIdx = i;
1701 longestLength = tmpLength;
1702 }
1703 }
1704
1705 // Set the longest path as the reference path
1706 setReferencePath(linevec[longestIdx], coordvec[longestIdx]);
1707 }
1708}
1709
1710
1711SoCoordinate3 * G4OpenInventorQtExaminerViewer::getCoordsNode(SoFullPath *path)
1712{
1713 SoLineSet *trajectory = (SoLineSet *)path->getTail();
1714 SoSeparator * grpNode = (SoSeparator*)(((SoFullPath*)path)->getNodeFromTail(1));
1715 int nodeIndex = grpNode->findChild(trajectory);
1716 SoNode * tmpNode;
1717
1718 // We allow only 100 iterations, in case the node isn't found
1719 // (should take only a few iterations)
1720 for (int i = 0; i < 100; ++i) {
1721 --nodeIndex;
1722
1723 tmpNode = grpNode->getChild(nodeIndex);
1724 if(tmpNode->getTypeId() == SoCoordinate3::getClassTypeId()) {
1725 //node found
1726 return (SoCoordinate3 *)tmpNode;
1727 }
1728 }
1729 return NULL; //coords node not found
1730}
1731
1732
1733// Displays scene elements on the right side of listsDialog.
1734// else: scene graph is searched for Geant4_SoPolyhedron type nodes
1736{
1737 std::string field, eltName;
1738
1739 std::map<std::string, int> duplicates;
1740 std::map<std::string, int> sceneElts;
1741 SoSearchAction search;
1742 Geant4_SoPolyhedron *node;
1743 SoGroup *root = (SoGroup *)getSceneManager()->getSceneGraph();
1744
1745 SoBaseKit::setSearchingChildren(TRUE);
1746
1747 search.reset();
1748 search.setSearchingAll(TRUE);
1749 search.setInterest(SoSearchAction::ALL);
1750 search.setType(Geant4_SoPolyhedron::getClassTypeId(), 0);
1751
1752 // FWJ DEBUG
1753 // G4cout << "Searching for elements....." << G4endl;
1754 search.apply(root);
1755
1756 SoPathList &pl = search.getPaths();
1757
1758
1759 // First find which names occur more than once so we can append a counter to them
1760 for (int i = 0; i < pl.getLength(); i++) {
1761 SoFullPath *path = (SoFullPath *)pl[i];
1762 node = (Geant4_SoPolyhedron *)path->getTail();
1763 eltName = node->getName();
1764 // G4cout << " FOUND " << i << " " << eltName << G4endl;
1765 if(duplicates.count(eltName))
1766 duplicates[eltName]++;
1767 else
1768 duplicates[eltName] = 1;
1769 }
1770
1771 for(int i = 0; i < pl.getLength(); i++) {
1772 float x,y,z;
1773 std::stringstream ssCount;
1774 SoFullPath *path = (SoFullPath *)pl[i];
1775 node = (Geant4_SoPolyhedron *)path->getTail();
1776 eltName = node->getName();
1777 field = eltName;
1778 if(duplicates[eltName] == 1)
1779 ssCount << "";//duplicates[field]
1780 else {
1781 if(sceneElts.count(eltName))
1782 sceneElts[eltName]++;
1783 else
1784 sceneElts[eltName] = 1;
1785
1786 ssCount << sceneElts[eltName];
1787 field += "_";
1788 }
1789
1790 field += ssCount.str();
1791
1792 SoGetBoundingBoxAction bAction(getViewportRegion());
1793 bAction.apply(path);
1794 SbBox3f bBox = bAction.getBoundingBox();
1795
1796 SbVec3f centr = bBox.getCenter();
1797 centr.getValue(x,y,z);
1798
1799 path->ref();
1800 sceneElement el = { field, path, centr, 0.0 };
1801 sceneElements.push_back(el);
1802 }
1803}
1804
1805
1807{
1808 float x,y,z;
1809 a.getValue(x,y,z);
1810 return x*x + y*y + z*z;
1811}
1812
1813
1815 float &dist,
1816 SbVec3f &closestPoint,
1817 int &index)
1818{
1819 // a : Previous point on trajectory
1820 // b : Next point on trajectory
1821 // q : the point in space
1822 // dab, daq, dbq: distance between a & b, a & q, b & q
1823 //
1824 // Theory: A point p on a line ab is defined as:
1825 //
1826 // p(t) = a+t?(b?a)
1827 //
1828 // note: All are vectors except the parameter t
1829 //
1830 // When t is between 0 and 1 the point p is situated between a and b on ab.
1831 // The point p is defined in terms of the parameter t, subsequently so does
1832 // the distance from the query point q to the point p. To find the minimum
1833 // of that distance we differentiate it and set equal to zero:
1834 //
1835 // diff(Norm(p(t)- q)) = 0
1836 //
1837 // note: diff means taking the derivative with regard to t
1838 //
1839 // The resulting t is given in the code below. The square of the distance
1840 // between p and q is given by:
1841 //
1842 // d^2 = (Norm(p(t)-q))^2
1843 //
1844 // The expression found is given in the code below (current_dist)
1845 //
1846 // Ref: http://programmizm.sourceforge.net/blog/2012/
1847 // distance-from-a-point-to-a-polyline
1848 //
1849 // --PLG
1850
1851 const std::size_t count = refParticleTrajectory.size();
1852 assert(count>0);
1853
1854 SbVec3f b = refParticleTrajectory[0];
1855 SbVec3f dbq = b - q;
1856 float sqrDist = sqrlen(dbq);
1857 closestPoint = b;
1858 index = 0;
1859 for (std::size_t i = 1; i < count; ++i) {
1860 const SbVec3f a = b;
1861 const SbVec3f daq = dbq;
1862 b = refParticleTrajectory[i];
1863 dbq = b - q;
1864 const SbVec3f dab = a - b;
1865
1866 float dab_x, dab_y, dab_z;
1867 dab.getValue(dab_x,dab_y,dab_z);
1868 float daq_x, daq_y, daq_z;
1869 daq.getValue(daq_x, daq_y, daq_z);
1870 float dbq_x, dbq_y, dbq_z;
1871 dbq.getValue(dbq_x, dbq_y, dbq_z);
1872
1873 const float inv_sqrlen = 1./sqrlen(dab);
1874 const float t = (dab_x*daq_x + dab_y*daq_y + dab_z*daq_z)*inv_sqrlen;
1875
1876 if (t<0.) {
1877 // The trajectory point occurs before point a
1878 // Go to the next point
1879 continue;
1880 }
1881 float current_dist;
1882 if (t<=1.) {
1883 // The trajectory point occurs between a and b.
1884 // Compute the distance to that point
1885 current_dist = daq_x*daq_x + daq_y*daq_y + daq_z*daq_z
1886 - t*(daq_x*dab_x + daq_y*dab_y + daq_z*dab_z)
1887 + t*t*(dab_x*dab_x + dab_y*dab_y + dab_z*dab_z);
1888 }
1889 else { //t>1.
1890 // The trajectory point occurs after b.
1891 // Get the distance to point b
1892 current_dist = sqrlen(dbq);
1893 }
1894
1895 if (current_dist < sqrDist) {
1896 sqrDist = current_dist;
1897 closestPoint = a + t*(b-a);
1898 index = (int) i;
1899 }
1900 }
1901
1902 dist = std::sqrt(sqrDist);
1903}
1904
1905
1907{
1908 if(refParticleTrajectory.empty())
1909 return;
1910
1911 float * trajLength = new float[refParticleTrajectory.size()];
1912 typedef std::map<elementForSorting, sceneElement> sortedMap;
1913 sortedMap sorted;
1914
1915 // For every point on the reference trajectory, compute
1916 // the total length from the start
1917 SbVec3f prevPoint;
1918 std::vector<SbVec3f>::iterator itRef = refParticleTrajectory.begin();
1919 int trajIndex = 0;
1920 prevPoint = *itRef;
1921 trajLength[trajIndex] = 0.0;
1922 ++itRef;
1923 ++trajIndex;
1924 for(; itRef != refParticleTrajectory.end(); ++itRef, ++trajIndex) {
1925 trajLength[trajIndex] = trajLength[trajIndex-1] + (*itRef - prevPoint).length();
1926 prevPoint = *itRef;
1927 }
1928
1929 // Compute the smallest distance between the element
1930 // and the reference trajectory (find the closest point),
1931 // then map the element to the trajectory length of that
1932 // point (calculated above)
1933 SoGetBoundingBoxAction bAction(getViewportRegion());
1934 SbVec3f elementCoord;
1935 std::vector<sceneElement>::iterator itEl;
1936 int elementIndex;
1938 for(itEl = sceneElements.begin(), elementIndex = 0;
1939 itEl != sceneElements.end(); ++itEl, ++elementIndex) {
1940 bAction.apply(itEl->path);
1941
1942 // FWJ sceneElement already has a center
1943 elementCoord = itEl->center;
1944 // ... and this sometimes returns an empty box!
1945 // elementCoord = bAction.getBoundingBox().getCenter();
1946 // if (bAction.getBoundingBox().isEmpty()) {
1947 // G4cout << "sortElements: Box is empty!" << G4endl;
1948 // G4cout << " element name=" << itEl->name << G4endl;
1949 // }
1950
1951 int index;
1952 distanceToTrajectory(elementCoord, el.smallestDistance, el.closestPoint, index);
1953 itEl->closestPointZCoord = el.closestPointZCoord = trajLength[index];
1954 el.distanceToBeamlineStart = (itEl->center - refParticleTrajectory[0]).length();
1955
1956 // This map of the scene elements (or their coordinates rather)
1957 // is automatically sorted by trajectory length (Z coord), then
1958 // by the distance between the element and the point in case the Z coord
1959 // is the same as another element. This is done by using as a key
1960 // an element structure which implements the operator for weak ordering
1961 sorted.insert(std::make_pair(el,*itEl));
1962 }
1963
1964 // store the sorted elements into the vector field
1965 sceneElements.clear();
1966
1967 sortedMap::iterator itSorted = sorted.begin();
1968 for(; itSorted != sorted.end(); itSorted++)
1969 sceneElements.push_back(itSorted->second);
1970
1971 zcoordSetFlag = true;
1972
1974
1975 delete[] trajLength;
1976}
1977
1978
1980{
1981 // FWJ DEBUG
1982 // G4cout << "Populating ELEMENT LIST..." << G4endl;
1983
1984 AuxWindowDialog->listWidget1->clear();
1985 // int size = sceneElements.size();
1986
1987 std::vector<sceneElement>::const_iterator it;
1988 std::stringstream ss;
1989
1990 for(it=sceneElements.begin(); it!=sceneElements.end(); ++it) {
1991 ss << it->name;
1992 if(zcoordSetFlag)
1993 ss << " [" << it->closestPointZCoord << "]";
1994
1995 new QListWidgetItem(ss.str().c_str(), AuxWindowDialog->listWidget1);
1996 ss.str("");
1997 }
1998}
1999
2000
2001// Called when user clicks a scene element in listsDialog.
2002// Zooms onto that element.
2003void
2004G4OpenInventorQtExaminerViewer::LookAtSceneElementCB(QListWidgetItem* item)
2005{
2006 char* value;
2007 std::string elementField;
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
2017 value = strdup(qPrintable(item->text()));
2018 // G4cout << "LOOKING FOR BOOKMARK " << value << G4endl;
2019
2022 if (animateSensor->isScheduled())
2023 animateSensor->unschedule();
2024 setSuperimpositionEnabled(superimposition, FALSE);
2025 maxSpeed = 0.0f;
2026 scheduleRedraw();
2027 restoreCamera();
2029 } else if (currentState == VIEWPOINT)
2030 setSuperimpositionEnabled(superimposition, FALSE);
2031
2032 elementField = value;
2033
2034 std::size_t idx = elementField.find_last_of("[");
2035 if(idx == std::string::npos)
2036 idx = elementField.size(); //if "[" not found for whatever reason (list not sorted)
2037 else
2038 idx--; // To get rid of the space that is between the name and '['
2039
2040 bool error = false;
2041 SoFullPath *path;
2042 SoSearchAction search;
2043 SoNode *root = getSceneManager()->getSceneGraph();
2044 int counter;
2045 std::size_t idxUnderscore = elementField.find_last_of("_");
2046
2047 parseString<int>(counter,
2048 elementField.substr(idxUnderscore + 1, idx), error);
2049
2050 SoBaseKit::setSearchingChildren(TRUE);
2051 search.reset();
2052 search.setSearchingAll(TRUE);
2053
2054 // G4cout << " Starting search for elementField " << elementField
2055 // << G4endl;
2056
2057 if(error) { // No counter is present => element name was not modified
2058 curEltName = elementField.substr(0, idx);
2059 search.setName(curEltName.c_str());
2060 search.apply(root);
2061
2062 path = (SoFullPath *)search.getPath();
2063 }
2064 else {
2065 curEltName = elementField.substr(0, idxUnderscore);
2066 search.setInterest(SoSearchAction::ALL);
2067 search.setName(curEltName.c_str());
2068 search.apply(root);
2069
2070 SoPathList &pl = search.getPaths();
2071 path = (SoFullPath *)pl[counter - 1]; // Since counter starts at 1, not 0
2072 }
2073
2074 G4ThreeVector global;
2075
2076 // FWJ FLIP THIS
2077 if ((idx > 0) && (path)) {
2078
2079 if(!refParticleTrajectory.empty()) {
2080
2081 SoGetBoundingBoxAction bAction(getViewportRegion());
2082 bAction.apply(path);
2083 SbBox3f bBox = bAction.getBoundingBox();
2084 SbVec3f elementCoord = bBox.getCenter();
2085
2086 refParticleIdx = 0;
2087 SbVec3f p;
2088
2089 float absLengthNow, absLengthMin;
2090 int maxIdx = (int) refParticleTrajectory.size() - 2;
2091 int targetIdx = 0;
2092 SbVec3f dir;
2093
2095 absLengthMin = (p - elementCoord).length();
2097
2098 // Find a ref. particle's point closest to element's global coords
2099 while (refParticleIdx < maxIdx) {
2101 absLengthNow = (p - elementCoord).length();
2102
2103 if (absLengthNow < absLengthMin) {
2104 absLengthMin = absLengthNow;
2105 targetIdx = refParticleIdx;
2106 }
2108 }
2109
2110 if (currentState != BEAMLINE) { // Set up default zoom
2111 SbVec3f p1, pN;
2113 prevParticleDir = SbVec3f(0,0,0); //so that moveCamera() knows sets default parameters
2114
2115 p1 = prevPt = refParticleTrajectory[0];
2117 distance = (pN - p1).length() / 10;
2118
2119 // FWJ Rather than switching to a default height, it is more flexible
2120 // to keep the same height(magnification) while moving the camera.
2121 // if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2122 // ((SoOrthographicCamera *) cam)->height.setValue(defaultHeight);
2123 // // FWJ Restore the default height instead of hard-wired value
2124 // // ((SoOrthographicCamera *) cam)->height.setValue(10000.0f);
2125 // }
2126 // else if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2127
2128 // FWJ required to avoid extreme perspective after camera move:
2129 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2130 ((SoPerspectiveCamera*)cam)->heightAngle.setValue(defaultHeightAngle);
2131
2132 } else {
2133 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId()))
2134 distance = (prevPt - cam->position.getValue()).length();
2135 }
2136 refParticleIdx = targetIdx;
2137
2138 //////////////////////////////////////////////////////////////
2139 setSuperimpositionEnabled(superimposition, TRUE);
2140 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
2141 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
2142 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
2143 scheduleRedraw();
2144 //////////////////////////////////////////////////////////////
2145
2147
2148 }
2149
2150 else {
2151 offsetFromCenter.setValue(0, 0, 1);
2152 distance = 50;// small number since using viewAll() for default zoom
2153 upVector.setValue(0, 1, 0);
2155 cam->viewAll(path, getViewportRegion());
2156 }
2157 }
2158
2159}
2160
2161
2162void G4OpenInventorQtExaminerViewer::FileLoadRefPathCB()
2163{
2164 // G4cout << "File: Load Ref Path CALLBACK" << G4endl;
2165
2166 QFileDialog filedialog(getParentWidget(), tr("Load Reference Path"));
2167 filedialog.setFileMode(QFileDialog::AnyFile);
2168 filedialog.setFont(*font);
2169 if (!filedialog.exec()) return;
2170 QStringList filenameinlist = filedialog.selectedFiles();
2171 QString filenamein = filenameinlist[0];
2172
2173 // G4cout << "Input file name is " << qPrintable(filenamein) << G4endl;
2174
2175 char* filename = new char[filenamein.size()+1];
2176 filename = strdup(qPrintable(filenamein));
2177 // G4cout << "char[] file name is " << filename << G4endl;
2178
2179 std::ifstream ifs(filename);
2180 if(ifs.is_open()) {
2181 refParticleTrajectory.clear();
2182 float x,y,z;
2183 while(ifs >> x >> y >> z) {
2184 refParticleTrajectory.push_back(SbVec3f(x,y,z));
2185 }
2186 ifs.close();
2187 } else {
2188 QMessageBox msgbox;
2189 msgbox.setFont(*font);
2190 QString messagetxt = "Reference Path file not found: ";
2191 messagetxt.append(filenamein);
2192 msgbox.setText(messagetxt);
2193 msgbox.exec();
2194 return;
2195 }
2196 if (refParticleTrajectory.size() < 2) {
2197 QMessageBox msgbox;
2198 msgbox.setFont(*font);
2199 QString messagetxt = "Invalid Reference Path";
2200 msgbox.setText(messagetxt);
2201 msgbox.exec();
2202 return;
2203 }
2204 // Following setReferencePath() ...
2208 sortElements();
2209}
2210
2211
2212void G4OpenInventorQtExaminerViewer::FileSaveRefPathCB()
2213{
2214 // G4cout << "File: Save Ref Path CALLBACK" << G4endl;
2215
2216 QFileDialog filedialog(getParentWidget(), tr("Save Reference Path"));
2217 filedialog.setFileMode(QFileDialog::AnyFile);
2218 // To enable confirmation of overwriting
2219 filedialog.setAcceptMode(QFileDialog::AcceptSave);
2220 filedialog.setFont(*font);
2221 if (!filedialog.exec()) return;
2222 QStringList filenameinlist = filedialog.selectedFiles();
2223 QString filenamein = filenameinlist[0];
2224
2225 // G4cout << "Input file name is " << qPrintable(filenamein) << G4endl;
2226
2227 char* filename = new char[filenamein.size()+1];
2228 filename = strdup(qPrintable(filenamein));
2229 // G4cout << "char[] file name is " << filename << G4endl;
2230
2231 std::ofstream ofs(filename);
2232 if (ofs.is_open()) {
2233 float x,y,z;
2234 for (unsigned int i=0; i < refParticleTrajectory.size(); ++i) {
2235 refParticleTrajectory[i].getValue(x,y,z);
2236 ofs << x << " " << y << " " << z << "\n";
2237 }
2238 ofs.close();
2239 } else {
2240 QMessageBox msgbox;
2241 msgbox.setFont(*font);
2242 QString messagetxt = "Error opening file ";
2243 messagetxt.append(filenamein);
2244 msgbox.setText(messagetxt);
2245 msgbox.exec();
2246 }
2247
2248}
2249
2251{
2252 if(refParticleTrajectory.empty())
2253 return;
2254
2255 SbVec3f p1, p2, p3, dirNow, dirNxt, dir, p2_tmp, p_start, p_corner, p_nxt;
2256 float avgDistBtwPts = 0;
2257 float totalDistBtwPts = 0;
2258 std::vector<SbVec3f> newRefParticleTrajectory;
2259 SbVec3f refPoint;
2260 std::size_t size = refParticleTrajectory.size() - 1;
2261 int numOfPts = 0;
2262 for (std::size_t i = 0; i < size; ++i) {
2263 p1 = refParticleTrajectory[i];
2264 p2 = refParticleTrajectory[i + 1];
2265 if (p1 == p2)
2266 continue;
2267 numOfPts++;
2268 totalDistBtwPts += (p2 - p1).length();
2269 }
2270 // Nothing useful to do (and fix Coverity)
2271 if (numOfPts <= 2) return;
2272
2273 avgDistBtwPts = totalDistBtwPts / numOfPts;
2274 float minDistAllowed = 0.75 * avgDistBtwPts;
2275 // float maxDistAllowed = 1.25 * avgDistBtwPts; // Pts tend to be close not far
2276
2277 float x, y, z;
2278 std::size_t i = 0, j = 0;
2279 while (i < size) {
2280 p1 = refParticleTrajectory[i];
2281 p2 = refParticleTrajectory[i + 1];
2282
2283 refPoint = p1;
2284 p1.getValue(x, y, z);
2285
2286 newRefParticleTrajectory.push_back(refPoint);
2287
2288 j = i;
2289 while ((p2 - p1).length() < minDistAllowed && j < (size - 1)) {
2290 j++;
2291
2292 p1 = refParticleTrajectory[j];
2293 p2 = refParticleTrajectory[j + 1];
2294 }
2295 if (j != i)
2296 i = j + 1;
2297 else
2298 i++;
2299 }
2300
2301 refParticleTrajectory.clear();
2302 refParticleTrajectory = newRefParticleTrajectory;
2303}
2304
2305
2307{
2308 SoCamera *cam = getCamera();
2309 camB4Animation.viewportMapping = cam->viewportMapping.getValue();
2310 camB4Animation.position = cam->position.getValue();
2311 camB4Animation.orientation = cam->orientation.getValue();
2312 camB4Animation.aspectRatio = cam->aspectRatio.getValue();
2313 camB4Animation.nearDistance = cam->nearDistance.getValue();
2314 camB4Animation.farDistance = cam->farDistance.getValue();
2315 camB4Animation.focalDistance = cam->focalDistance.getValue();
2316
2317 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2319 ((SoPerspectiveCamera *) cam)->heightAngle.getValue();
2321 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2323 ((SoOrthographicCamera *) cam)->height.getValue();
2325 }
2326}
2327
2328
2330{
2331 SoCamera *cam = getCamera();
2332
2333 cam->viewportMapping = camB4Animation.viewportMapping;
2334 cam->position = camB4Animation.position;
2335 cam->orientation = camB4Animation.orientation;
2336 cam->aspectRatio = camB4Animation.aspectRatio;
2337 cam->nearDistance = camB4Animation.nearDistance;
2338 cam->farDistance = camB4Animation.farDistance;
2339 cam->focalDistance = camB4Animation.focalDistance;
2340
2341 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2343 toggleCameraType();
2344 cam = getCamera();
2345 ((SoOrthographicCamera *) cam)->height.setValue(
2347 } else
2348 ((SoPerspectiveCamera *) cam)->heightAngle.setValue(
2350 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2352 toggleCameraType();
2353 cam = getCamera();
2354 ((SoPerspectiveCamera *) cam)->heightAngle.setValue(
2356 } else
2357 ((SoOrthographicCamera *) cam)->height.setValue(
2359 }
2360}
2361
2362
2364 SoSensor *sensor)
2365{
2366 SbTime curTime = SbTime::getTimeOfDay();
2368
2369 SoTimerSensor* s = (SoTimerSensor*) sensor;
2370
2371 float t = float((curTime - s->getBaseTime()).getValue())
2372 / This->animateBtwPtsPeriod;
2373
2374 if ((t > 1.0f) || (t + s->getInterval().getValue() > 1.0f))
2375 t = 1.0f;
2376 SbBool end = (t == 1.0f);
2377
2378 if (end) {
2379 This->animateSensorRotation->unschedule();
2380 if(This->rotCnt) {
2381 // rotations left
2382 This->rotateCamera();
2383 }
2384 else {
2385 // rotation over
2386 This->currentState = This->prevState;
2387 return;
2388 }
2389 }
2390
2391}
2392
2393
2394// Called repeatedly during reference particle animation
2395
2397 SoSensor *sensor)
2398{
2399 SbTime curTime = SbTime::getTimeOfDay();
2401 SoCamera *cam = This->getCamera();
2402 SoTimerSensor* s = (SoTimerSensor*) sensor;
2403
2404 float t = float((curTime - s->getBaseTime()).getValue())
2405 / This->animateBtwPtsPeriod;
2406
2407 if ((t > 1.0f) || (t + s->getInterval().getValue() > 1.0f))
2408 t = 1.0f;
2409 SbBool end = (t == 1.0f);
2410
2411 cam->orientation = SbRotation::slerp(This->camStartOrient, This->camEndOrient, t);
2412 cam->position = This->camStartPos + (This->camEndPos - This->camStartPos) * t;
2413
2414 if (end) {
2415 This->animateSensor->unschedule();
2416
2417 if (This->currentState == ANIMATION) {
2418 if (This->refParticleIdx < (int) (This->refParticleTrajectory.size() - 1))
2419 This->animateRefParticle();
2420 else {
2422 This->speedStep = START_STEP;
2423 }
2424 }
2425 if (This->currentState == REVERSED_ANIMATION) {
2426 if (This->refParticleIdx >= 1)
2427 This->animateRefParticle();
2428 else {
2430 This->speedStep = START_STEP;
2431 }
2432 }
2433 }
2434}
2435
2436
2438{
2439 if (SoQtExaminerViewer::isAnimating())
2440 stopAnimating();
2441
2442 SbRotation rot;
2443 SbVec3f p1, p2, p2_tmp, camUpV, camD, camD_tmp, leftRightAxis;
2444 float x1, y1, z1, x2, y2, z2;
2445
2446 if (currentState == ANIMATION) {
2449 } else if (currentState == REVERSED_ANIMATION) {
2452 } else if (currentState == PAUSED_ANIMATION) {
2453 if (refParticleIdx < (int) refParticleTrajectory.size()) {
2456 } else {
2459 }
2460 }
2461 p1.getValue(x1, y1, z1);
2462 p2.getValue(x2, y2, z2);
2463
2464 camD = p2 - p1;
2465 camD.normalize();
2466
2467 p2_tmp.setValue(x2, y1, z2);
2468 camD_tmp = p2_tmp - p1;
2469 camD_tmp.normalize();
2470
2471 camUpV.setValue(0, 1, 0);
2472 rot.setValue(camD_tmp, camD);
2473 rot.multVec(camUpV, camUpV);
2474
2475 leftRightAxis = camD.cross(camUpV);
2476
2477 myCam->position = p1;
2478 myCam->pointAt(p2, camUpV);
2479
2480 // Update camera position
2481 p1 = p1 + (up_down * camUpV) + (left_right * leftRightAxis);
2482 myCam->position = p1;
2483 // FWJ Try look-ahead here
2484 int idx = refParticleIdx + pathLookahead;
2485 idx = std::min(idx, (int)refParticleTrajectory.size() - 1);
2486 myCam->pointAt(refParticleTrajectory[idx], camUpV);
2487 // myCam->pointAt(refParticleTrajectory[idx], camUpVec);
2488 myCam->focalDistance = 0.1f;
2489}
2490
2491
2493{
2494 G4OpenInventorQtExaminerViewer::ToolsRefPathStartCB();
2495}
2496
2497
2498void G4OpenInventorQtExaminerViewer::ToolsRefPathStartCB()
2499{
2500 if (!refParticleTrajectory.size()) {
2501 QMessageBox msgbox;
2502 msgbox.setFont(*font);
2503 QString messagetxt = "No current reference path";
2504 msgbox.setText(messagetxt);
2505 msgbox.exec();
2506 return;
2507 }
2508
2509 if (currentState == ROTATING)
2510 return;
2513 if (animateSensor->isScheduled())
2514 animateSensor->unschedule();
2515 setSuperimpositionEnabled(superimposition, FALSE);
2516 maxSpeed = 0.0f;
2517 scheduleRedraw();
2518 } else {
2519 saveCurCamera();
2522 }
2523
2524 if (SoQtExaminerViewer::isAnimating())
2525 stopAnimating();
2526
2527 up_down = 0;
2528 left_right = 0;
2529 step = 1;
2530
2531 refParticleIdx = 0;
2533 setSuperimpositionEnabled(superimposition, TRUE);
2534 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
2535 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
2536 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
2537 scheduleRedraw();
2538
2539 // FWJ Disabled: this is set in moveCamera()
2540 // Zoom in
2541 // SoCamera *cam = getCamera();
2542 // cam->focalDistance = 0.1f;
2543
2544 prevParticleDir = SbVec3f(0,0,0);
2545
2546 //Default zoom
2547 SbVec3f p1 = refParticleTrajectory[0];
2548 SbVec3f pN = refParticleTrajectory[refParticleTrajectory.size() - 1];
2549 distance = (pN - p1).length() / 10;
2550
2551 moveCamera(distance, true);
2552}
2553
2554
2555void G4OpenInventorQtExaminerViewer::ToolsRefPathInvertCB()
2556{
2557 invertRefPath();
2558}
2559
2560
2562{
2563 std::reverse(refParticleTrajectory.begin(),
2564 refParticleTrajectory.end());
2566 sortElements();
2567}
2568
2569
2571{
2572 SoCamera *cam = getCamera();
2573
2574 camStartPos = cam->position.getValue();
2575 camStartOrient = cam->orientation.getValue();
2576
2577 if (currentState != BEAMLINE)
2579
2580 camEndPos = myCam->position.getValue();
2581 camEndOrient = myCam->orientation.getValue();
2582
2583 if (animateSensor->isScheduled())
2584 animateSensor->unschedule();
2585
2586 animateSensor->setBaseTime(SbTime::getTimeOfDay());
2587 animateSensor->setInterval(SbTime(0.02));
2588
2589 animateSensor->schedule();
2590}
2591
2592
2594{
2595 escapeCallback = callback;
2596}
2597
2598
2600{
2601 // FWJ DEBUG
2602 // G4cout << "SCENE CHANGE callback" << G4endl;
2603 // NOTE: could/should be disabled during animation
2604
2607 if(This->newEvents) {
2608 This->findAndSetRefPath();
2609 This->newEvents = false;
2610 }
2611}
2612
2613
2614//////////////////////////////////// BOOKMARKS ///////////////////////////
2615
2616// Adds bookmarks to listsDialog.
2617
2619{
2620 std::size_t size = viewPtList.size();
2621 if (!size) return;
2622
2623 for (std::size_t i = 0; i < size; ++i) {
2624 new QListWidgetItem(viewPtList[i].viewPtName,
2625 AuxWindowDialog->listWidget);
2626 }
2627}
2628
2629
2630// Converts a string type word into a float type.
2631
2632template<class T>
2633void G4OpenInventorQtExaminerViewer::parseString(T &t, const std::string &s,
2634 bool &error)
2635{
2636 std::istringstream str(s);
2637 if ((str >> t).fail())
2638 error = true;
2639}
2640
2641
2642void
2643G4OpenInventorQtExaminerViewer::FileOpenBookmarkCB()
2644{
2645 // FWJ DEBUG
2646 // G4cout << "File: Open Bookmark File CALLBACK" << G4endl;
2647 QFileDialog filedialog(getParentWidget(), tr("Open bookmark file"));
2648 filedialog.setFileMode(QFileDialog::ExistingFile);
2649 filedialog.setFont(*font);
2650 if (!filedialog.exec()) return;
2651 QStringList filenameinlist = filedialog.selectedFiles();
2652 QString filenamein = filenameinlist[0];
2653
2654 char* filename = new char[filenamein.size()+1];
2655 filename = strdup(qPrintable(filenamein));
2656 // G4cout << "char[] file name is " << filename << G4endl;
2657
2658 fileIn.close();
2659 fileIn.open(filename);
2660 if (fileIn.fail()) {
2661 QMessageBox msgbox;
2662 msgbox.setFont(*font);
2663 QString messagetxt = "Error opening file: ";
2664 messagetxt.append(filenamein);
2665 msgbox.setText(messagetxt);
2666 msgbox.exec();
2667 // G4cout << "ERROR opening file " << filename << G4endl;
2668 fileIn.clear();
2669 return;
2670 }
2671 // Opens a file without erasing it
2673
2674 if (!loadViewPts()) {
2675 QMessageBox msgbox;
2676 msgbox.setFont(*font);
2677 QString messagetxt = "Error reading bookmark file: ";
2678 messagetxt.append(filenamein);
2679 msgbox.setText(messagetxt);
2680 msgbox.exec();
2681 // G4cout << "ERROR reading bookmark file " << filename << G4endl;
2682 fileIn.clear();
2683 return;
2684 }
2685
2686 fileName = filename;
2687 fileOut.open(fileName.c_str(), std::ios::in);
2688 fileOut.seekp(0, std::ios::end);
2689
2690 addViewPoints();
2691
2692 // LATER: display filename in lists window
2693
2694 fileIn.close();
2695 fileIn.clear();
2696}
2697
2698// Called before loading a new viewpoint file.
2699// Resets member fields to default values.
2700
2702{
2703 viewPtIdx = -1;
2704 viewPtList.clear();
2705 // setSuperimpositionEnabled(superimposition, FALSE);
2706 // scheduleRedraw();
2708 if (fileOut.is_open()) fileOut.close();
2709
2710 AuxWindowDialog->listWidget->clear();
2711 AuxWindowDialog->lineEdit->setText(QString(""));
2712}
2713
2714
2715void
2716G4OpenInventorQtExaminerViewer::FileNewBookmarkCB()
2717{
2718 // G4cout << "File: Open New Bookmark File CALLBACK" << G4endl;
2719 QFileDialog filedialog(getParentWidget(), tr("Open new bookmark file"));
2720 filedialog.setFileMode(QFileDialog::AnyFile);
2721 // To enable confirmation of overwriting
2722 filedialog.setAcceptMode(QFileDialog::AcceptSave);
2723 // But change the "Save" button text
2724 filedialog.setLabelText(QFileDialog::Accept, QString("New"));
2725 filedialog.setFont(*font);
2726 if (!filedialog.exec()) return;
2727 QStringList filenameinlist = filedialog.selectedFiles();
2728 QString filenamein = filenameinlist[0];
2729
2730 // G4cout << "Input file name is " << qPrintable(filenamein) << G4endl;
2731
2732 char* filename = new char[filenamein.size()+1];
2733 filename = strdup(qPrintable(filenamein));
2734 // G4cout << "char[] file name is " << filename << G4endl;
2735
2737 fileName = filename;
2738 fileOut.open(fileName.c_str());
2739 if (fileOut.fail()) {
2740 QMessageBox msgbox;
2741 msgbox.setFont(*font);
2742 QString messagetxt = "Error opening new bookmark file: ";
2743 messagetxt.append(filenamein);
2744 msgbox.setText(messagetxt);
2745 msgbox.exec();
2746 // G4cout << "ERROR opening new bookmark file " << filename << G4endl;
2747 }
2748}
2749
2750
2751void
2752G4OpenInventorQtExaminerViewer::ToolsAnimateRefParticleCB()
2753{
2754 // G4cout << "Tools: Animate Ref Particle CALLBACK" << G4endl;
2755 if (!refParticleTrajectory.size()) {
2756 returnToAnim = true;
2757 G4warn << "No Reference Trajectory" << G4endl;
2758 return;
2759 }
2760
2761 ///////////////////////////////////////////////////////////////
2762 setSuperimpositionEnabled(superimposition, TRUE);
2764 axisSwitch->whichChild.setValue(SO_SWITCH_ALL);
2765 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_ALL);
2766 animSpeedSwitch->whichChild.setValue(SO_SWITCH_ALL);
2767 scheduleRedraw();
2768 ///////////////////////////////////////////////////////////////
2769
2770 SoCamera *cam = getCamera();
2771 // SbVec3f camDirOld, camDirNew, camDirNew_tmp, camUpVec, P0, P1, P1_tmp;
2772
2774 || currentState == ROTATING)
2775 return;
2776
2778
2779 saveCurCamera();
2782
2783 if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
2784 toggleCameraType();
2785 cam = getCamera();
2786 }
2787
2788 refParticleIdx = 0; // Set the camera to the starting point of the animation
2791 left_right = up_down = 0;
2792
2793 cam->focalDistance = 0.1f;
2794 ((SoPerspectiveCamera *) cam)->heightAngle = 0.50f;
2795 }
2796
2799
2800 cam->position = (myCam)->position.getValue();
2801 cam->orientation = (myCam)->orientation.getValue();
2802 animateRefParticle(); // Animate the camera
2803}
2804
2805
2806void
2807G4OpenInventorQtExaminerViewer::SaveViewPtCB()
2808{
2809 // G4cout << "AppButton: Save Viewpoint CALLBACK" << G4endl;
2810 // First get viewpoint name ...
2811 // EMULATING getViewPtNameCB ...
2812 // bool ok;
2813 // Note QString() returns an empty string
2814
2815 // NONE OF THE FOLLOWING CHANGES THE FONT: FORGET IT FOR NOW
2816 QInputDialog* inputdialog = new QInputDialog(getParentWidget());
2817 inputdialog->setFont(*font);
2818 inputdialog->setWindowTitle(tr("Enter a name for the bookmark"));
2819 inputdialog->setLabelText("Bookmark name");
2820 // inputdialog->setTextEchoMode(QLineEdit::Normal);
2821 inputdialog->adjustSize();
2822 QString namein;
2823 if (inputdialog->exec() == QDialog::Accepted)
2824 namein=inputdialog->textValue().trimmed();
2825 else
2826 return;
2827 if (namein.isEmpty()) return;
2828
2829 // This easier approach failed: unable to set font size
2830 // QString namein =
2831 // QInputDialog::getText(getParentWidget(),
2832 // tr("Enter a name for the bookmark"),
2833 // tr("Bookmark name"), QLineEdit::Normal,
2834 // QString(), &ok);
2835
2836 const int nVPName = MAX_VP_NAME + 1;
2837 char* name = new char[nVPName];
2838 // strncpy(name, strName.c_str(), nVPName);
2839 namein.truncate(MAX_VP_NAME);
2840
2841 QByteArray ba = namein.toLocal8Bit();
2842 name = strdup(ba.constData());
2843 // name = strdup(qPrintable(namein))
2844
2845 // FWJ DEBUG
2846 // G4cout << "QString is " << qPrintable(namein) << G4endl;
2847 // G4cout << "char[] is " << name << G4endl;
2848
2849 for (int i = 0; i < (int)viewPtList.size(); i++) {
2850 if (!strcmp(name, viewPtList[i].viewPtName)) {
2851 QMessageBox msgbox;
2852 msgbox.setText("Bookmark name is already in use");
2853 msgbox.exec();
2854 return;
2855 }
2856 }
2857
2858 if (viewPtIdx == -1) viewPtIdx = 0;
2859 saveViewPt(name);
2860
2861 saveViewPtItem = new QListWidgetItem(namein,
2862 AuxWindowDialog->listWidget);
2863 AuxWindowDialog->listWidget->setCurrentItem(saveViewPtItem);
2864 AuxWindowDialog->lineEdit->setText(namein);
2865}
2866
2867
2868// Saves current camera parameters to a viewpoint file.
2869
2871{
2872 SbVec3f axis;
2873 viewPtData tmp;
2874 float x, y, z, angle;
2875 SoCamera* camera = getCamera();
2876
2877 // NOTE: Xt VSN increments this at end of procedure
2878 // viewPtIdx++;
2879
2880 // FWJ DEBUG
2881 // G4cout << "saveViewPt: saving bookmark " << viewPtIdx << " " << name
2882 // << G4endl;
2883
2884 if (viewPtList.size() == 0) {
2886 }
2887
2888 tmp.viewPtName = name;
2889 tmp.viewportMapping = camera->viewportMapping.getValue();
2890 tmp.position = camera->position.getValue();
2891 tmp.orientation = camera->orientation.getValue();
2892 tmp.aspectRatio = camera->aspectRatio.getValue();
2893 tmp.nearDistance = camera->nearDistance.getValue();
2894 tmp.farDistance = camera->farDistance.getValue();
2895 tmp.focalDistance = camera->focalDistance.getValue();
2896
2897 // Save camera height (changed by zooming)
2898 if (camera->isOfType(SoPerspectiveCamera::getClassTypeId())) {
2899 tmp.height = ((SoPerspectiveCamera *) camera)->heightAngle.getValue();
2900 tmp.camType = PERSPECTIVE;
2901 } else if (camera->isOfType(SoOrthographicCamera::getClassTypeId())) {
2902 tmp.height = ((SoOrthographicCamera *) camera)->height.getValue();
2903 tmp.camType = ORTHOGRAPHIC;
2904 } else {
2905 SoDebugError::post("G4OpenInventorQtExaminerViewer::saveViewPtCB",
2906 "Only Perspective and Orthographic cameras are supported.");
2907 return;
2908 }
2909
2910 viewPtList.push_back(tmp);
2911
2912 // Now save the view point to a .txt file
2913 // FWJ DEBUG
2914 // G4cout << "saveViewPt: writing to Bookmark file " << fileName << G4endl;
2915
2916 std::string vpName = name;
2917
2918 while ((int) vpName.size() <= MAX_VP_NAME)
2919 vpName += " ";
2920
2921 fileOut << vpName << std::endl;
2922 tmp.position.getValue(x, y, z);
2923 fileOut << x << " " << y << " " << z << std::endl;
2924
2925 // Reusing x, y and z for storing the axis
2926 tmp.orientation.getValue(axis, angle);
2927 axis.getValue(x, y, z);
2928 fileOut << x << " " << y << " " << z << " " << angle << std::endl;
2929
2930 fileOut << tmp.camType << " " << tmp.height << std::endl;
2931 fileOut << tmp.focalDistance << " ";
2932 fileOut << tmp.nearDistance << " ";
2933 fileOut << tmp.farDistance << std::endl;
2934 fileOut << tmp.viewportMapping << " ";
2935 fileOut << tmp.aspectRatio << "\n" << std::endl;
2936 fileOut.flush();
2937
2938 viewPtIdx++;
2939
2940 // FWJ DEBUG
2941 // G4cout << "saveViewPt: finished writing to file" << G4endl <<
2942 // " Next viewPtIdx is " << viewPtIdx << G4endl;
2943}
2944
2945
2946// Updates the viewPtIdx in a viewpoint file.
2947
2949{
2950 std::string idxStr;
2951 std::stringstream out;
2952
2953 out << viewPtIdx;
2954 idxStr = out.str();
2955 fileOut.seekp(0, std::ios::beg);
2956
2957 while ((int) idxStr.length() < MAX_VP_IDX) {
2958 idxStr += " ";
2959 }
2960
2961 // FWJ DEBUG
2962 // G4cout << "writeViewPtIdx: " << viewPtIdx << G4endl;
2963 fileOut << idxStr << "\n";
2964 fileOut.flush();
2965 fileOut.seekp(0, std::ios::end);
2966}
2967
2968
2969// Receives the name of the bookmark clicked and searches for it in viewPtList.
2970
2971void G4OpenInventorQtExaminerViewer::LoadBookmarkCB(QListWidgetItem* item)
2972{
2973 // FWJ DEBUG
2974 // G4cout << "AuxWindow: listWidget LoadBookmark CALLBACK" << G4endl;
2975
2976 const int nVPName = MAX_VP_NAME + 1;
2977 char* vpName = new char[nVPName];
2978
2979 vpName = strdup(qPrintable(item->text()));
2980 // G4cout << "LOOKING FOR BOOKMARK " << vpName << G4endl;
2981
2982 for (int i = 0; i < (int)viewPtList.size(); i++) {
2983 if (!strcmp(viewPtList[i].viewPtName, vpName)) {
2984 viewPtIdx = i;
2985 break;
2986 }
2987 }
2988 // G4cout << " FOUND viewPtIdx " << viewPtIdx << G4endl;
2989
2991 setViewPt();
2992 AuxWindowDialog->lineEdit->setText(item->text());
2993}
2994
2995
2996// Sets the viewpoint based on camera data that viewPtIdx is pointing to.
2997
2999{
3001 || currentState == ROTATING) {
3002 if (animateSensor->isScheduled()) animateSensor->unschedule();
3003 setSuperimpositionEnabled(superimposition, FALSE);
3004 maxSpeed = 0.0f;
3005 scheduleRedraw();
3006 }
3007
3008 SoCamera * camera = getCamera();
3009 if (camera == NULL) {
3010 G4warn << "setViewPt: Camera is null. Unable to set the viewpoint." <<
3011 G4endl;
3012 // String dialogName = (char *) "Missing Camera Node";
3013 // std::string msg = "Camera is null. Unable to set the viewpoint.";
3014 // warningMsgDialog(msg, dialogName, NULL);
3015 return;
3016 }
3017
3018 if (!viewPtList.size()) {
3019 G4warn << "setViewPt: There are no viewpoints to load." << G4endl;
3020 // String dialogName = (char *) "Missing Viewpoints";
3021 // std::string msg = "There are no viewpoints to load.";
3022 // warningMsgDialog(msg, dialogName, NULL);
3023 return;
3024 }
3025
3026 if (SoQtExaminerViewer::isAnimating()) stopAnimating();
3027
3028 if (currentState != VIEWPOINT) {
3030 //////////////////////////////////////////////////////////////
3031 setSuperimpositionEnabled(superimposition, TRUE);
3032 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
3033 animSpeedOutlineSwitch->whichChild.setValue(SO_SWITCH_NONE);
3034 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
3035 scheduleRedraw();
3036 ///////////////////////////////////////////////////////////////
3037 }
3038
3039 curViewPtName = viewPtList[viewPtIdx].viewPtName;
3040 camera->viewportMapping = viewPtList[viewPtIdx].viewportMapping;
3041 camera->position = viewPtList[viewPtIdx].position;
3042 camera->orientation = viewPtList[viewPtIdx].orientation;
3043 camera->aspectRatio = viewPtList[viewPtIdx].aspectRatio;
3044 camera->nearDistance = viewPtList[viewPtIdx].nearDistance;
3045 camera->farDistance = viewPtList[viewPtIdx].farDistance;
3046 camera->focalDistance = viewPtList[viewPtIdx].focalDistance;
3047
3048 // Restore camera height (changed by zooming)
3049 if (camera->isOfType(SoPerspectiveCamera::getClassTypeId())) {
3050 if (viewPtList[viewPtIdx].camType == ORTHOGRAPHIC) {
3051 toggleCameraType();
3052 camera = getCamera();
3053 ((SoOrthographicCamera *) camera)->height.setValue(
3054 viewPtList[viewPtIdx].height);
3055 } else
3056 ((SoPerspectiveCamera *) camera)->heightAngle.setValue(
3057 viewPtList[viewPtIdx].height);
3058 } else if (camera->isOfType(SoOrthographicCamera::getClassTypeId())) {
3059 if (viewPtList[viewPtIdx].camType == PERSPECTIVE) {
3060 toggleCameraType();
3061 camera = getCamera();
3062 ((SoPerspectiveCamera *) camera)->heightAngle.setValue(
3063 viewPtList[viewPtIdx].height);
3064 } else
3065 ((SoOrthographicCamera *) camera)->height.setValue(
3066 viewPtList[viewPtIdx].height);
3067 } else {
3068 SoDebugError::post("G4OpenInventorQtExaminerViewer::setViewPt",
3069 "Only Perspective and Orthographic cameras are supported.");
3070 return;
3071 }
3072
3073}
3074
3075
3076void G4OpenInventorQtExaminerViewer::NextViewPtCB()
3077{
3078 // FWJ DEBUG
3079 // G4cout << "App Button: nextViewPt CALLBACK" << G4endl;
3080
3081 if (!viewPtList.size()) return;
3082 if (viewPtIdx >= (int)viewPtList.size() - 1)
3083 viewPtIdx = 0;
3084 else
3085 viewPtIdx++;
3086
3088 setViewPt();
3089 char* viewptname = viewPtList[viewPtIdx].viewPtName;
3090 AuxWindowDialog->lineEdit->setText(QString(viewptname));
3091}
3092
3093void G4OpenInventorQtExaminerViewer::PrevViewPtCB()
3094{
3095 // FWJ DEBUG
3096 // G4cout << "App Button: prevViewPt CALLBACK" << G4endl;
3097
3098 if (!viewPtList.size()) return;
3099 if (viewPtIdx == 0)
3100 viewPtIdx = (int) viewPtList.size() - 1;
3101 else
3102 viewPtIdx--;
3103
3105 setViewPt();
3106 char* viewptname = viewPtList[viewPtIdx].viewPtName;
3107 AuxWindowDialog->lineEdit->setText(QString(viewptname));
3108}
3109
3110
3111void G4OpenInventorQtExaminerViewer::AbbrOutputCB(bool checked)
3112{
3113 // FWJ DEBUG
3114 // G4cout << "App Button: abbrOutput CALLBACK" << G4endl;
3115
3116 abbrOutputFlag = checked;
3117}
3118
3119
3120void G4OpenInventorQtExaminerViewer::PickRefPathCB()
3121{
3122 // FWJ DEBUG
3123 // G4cout << "App Button: pickRefPath CALLBACK" << G4endl;
3124
3125 // Save viewing state and go to picking mode
3126 viewingBeforePickRef = isViewing();
3127 if(isViewing())
3128 setViewing(false);
3129 setComponentCursor(SoQtCursor(SoQtCursor::CROSSHAIR));
3130 pickRefPathFlag = true;
3131}
3132
3133
3134void G4OpenInventorQtExaminerViewer::SwitchWireFrameCB(bool checked)
3135{
3136 // FWJ DEBUG
3137 // G4cout << "App Button: switchWireFrame CALLBACK" << G4endl;
3138
3139 // if (switchWireFrameButton->isDown()) {
3140 if (checked) {
3141 setDrawStyle(SoQtViewer::STILL, SoQtViewer::VIEW_LINE);
3142 setDrawStyle(SoQtViewer::INTERACTIVE, SoQtViewer::VIEW_LINE);
3143 } else {
3144 setDrawStyle(SoQtViewer::STILL, SoQtViewer::VIEW_AS_IS);
3145 setDrawStyle(SoQtViewer::INTERACTIVE,
3146 SoQtViewer::VIEW_SAME_AS_STILL);
3147 }
3148}
3149
3150
3151void G4OpenInventorQtExaminerViewer::SwitchAxesCB(bool checked)
3152{
3153 // FWJ DEBUG
3154 // G4cout << "App Button: switchAxes CALLBACK" << G4endl;
3155 setFeedbackVisibility(checked);
3156 // if (checked) {
3157 // setFeedbackVisibility(TRUE);
3158 // } else {
3159 // setFeedbackVisibility(FALSE);
3160 // }
3161}
3162
3163
3164void G4OpenInventorQtExaminerViewer::DetachCB()
3165{
3166 // FWJ DEBUG
3167 // G4cout << "App Button: detach CALLBACK" << G4endl;
3168 uiQt->GetViewerTabWidget()->removeTab(uiQtTabIndex);
3169 viewerParent->setParent(viewerParent2);
3170 removeAppPushButton(detachButton);
3171 show();
3172}
3173
3174
3175void G4OpenInventorQtExaminerViewer::DeleteBookmarkCB()
3176{
3177 // FWJ DEBUG
3178 // G4cout << "Delete Button: DeleteBookmark CALLBACK" << G4endl;
3179
3180 // Maybe nothing selected yet
3181 QListWidgetItem* listitem = AuxWindowDialog->listWidget->currentItem();
3182 if (!listitem) return;
3183 if (!(listitem->isSelected())) return;
3184
3185 QString vpnamein = listitem->text();
3186
3187 const int nVPName = MAX_VP_NAME + 1;
3188 char* vpName = new char[nVPName];
3189 vpName = strdup(qPrintable(vpnamein));
3190 // G4cout << "DELETING bookmark " << vpName << G4endl;
3191
3192 deleteViewPt(vpName);
3193 delete listitem;
3194}
3195
3196// Deletes current viewpoint the user is looking at.
3197// Updates the input file and bookmarks as well.
3198
3200{
3201 std::string line;
3202 std::size_t end;
3203 fileIn.open(fileName.c_str());
3204 std::ofstream out("temporaryFile.txt");
3205
3206 if (!vpName)
3207 vpName = viewPtList[viewPtIdx].viewPtName;
3208
3209 getline(fileIn, line); // Printing the viewpoint idx
3210 out << line << "\n";
3211
3212 while (getline(fileIn, line)) {
3213 end = line.find_last_not_of(' ');
3214 line = line.substr(0, end + 1);
3215 if (!strcmp(line.c_str(), vpName)) { // Equal
3216 while (line.size()) {
3217 getline(fileIn, line);
3218 }
3219
3220 while (getline(fileIn, line))
3221 out << line << "\n";
3222 } else {
3223 while (line.size()) {
3224 out << line << "\n";
3225 getline(fileIn, line);
3226 }
3227 out << "\n";
3228 }
3229 }
3230
3231 std::size_t idx = 0; // Remove viewpoint from the vector
3232 std::size_t size = viewPtList.size();
3233 while (idx < size) {
3234 if (!strcmp(viewPtList[idx].viewPtName, vpName)) {
3235 viewPtList.erase(viewPtList.begin() + idx);
3236 break;
3237 }
3238 idx++;
3239 }
3240
3241 out.close();
3242 fileOut.close();
3243 fileIn.clear();
3244 fileIn.close();
3245
3246 // FWJ check return status: error popups needed here
3247 int istat = remove(fileName.c_str());
3248 if (istat == -1) {
3249 QMessageBox msgbox;
3250 msgbox.setFont(*font);
3251 QString messagetxt = "Error removing bookmarks file";
3252 // messagetxt.append(filenamein);
3253 msgbox.setText(messagetxt);
3254 msgbox.exec();
3255 // G4cout << "Error removing bookmarks file" << G4endl;
3256 }
3257 istat = rename("temporaryFile.txt", fileName.c_str());
3258 if (istat == -1) {
3259 QMessageBox msgbox;
3260 msgbox.setFont(*font);
3261 QString messagetxt = "Error renaming bookmarks file";
3262 // messagetxt.append(filenamein);
3263 msgbox.setText(messagetxt);
3264 msgbox.exec();
3265 // G4cout << "Error renaming bookmarks file" << G4endl;
3266 }
3267 fileOut.open(fileName.c_str(), std::ios::in);
3268 fileOut.seekp(0, std::ios::end);
3269
3270 if (!viewPtList.size()) { // viewPtList is empty
3271 curViewPtName = (char *) empty.c_str();
3272 scheduleRedraw();
3273 } else {
3274 if (viewPtIdx >= (int) viewPtList.size())
3275 viewPtIdx--;
3277 setViewPt();
3278 }
3279}
3280
3281
3282void G4OpenInventorQtExaminerViewer::RenameBookmarkCB()
3283{
3284 // FWJ DEBUG
3285 // G4cout << "Rename Button: RenameBookmark CALLBACK" << G4endl;
3286 // Maybe nothing selected yet
3287 QListWidgetItem* listitem = AuxWindowDialog->listWidget->currentItem();
3288 if (!listitem) return;
3289 if (!(listitem->isSelected())) return;
3290
3291 QString vpnamein = listitem->text();
3292
3293 const int nVPName = MAX_VP_NAME + 1;
3294 // char* vpName = new char[nVPName];
3295 // vpName = strdup(qPrintable(vpnamein));
3296 // G4cout << "RENAMING bookmark " << vpName << G4endl;
3297
3298 QInputDialog* inputdialog = new QInputDialog(getParentWidget());
3299 inputdialog->setFont(*font);
3300 inputdialog->setWindowTitle(tr("Enter"));
3301 inputdialog->setLabelText("New bookmark name");
3302 inputdialog->adjustSize();
3303 QString newnamein;
3304 if (inputdialog->exec() == QDialog::Accepted)
3305 newnamein=inputdialog->textValue().trimmed();
3306 else
3307 return;
3308 if (newnamein.isEmpty()) return;
3309
3310 char* newname = new char[nVPName];
3311 newname = strdup(qPrintable(newnamein));
3312
3313 std::size_t size = viewPtList.size();
3314 for (std::size_t i = 0; i < size; ++i) {
3315 if (!strcmp(newname, viewPtList[i].viewPtName)) {
3316 QMessageBox msgbox;
3317 msgbox.setFont(*font);
3318 msgbox.setText("Bookmark name is already in use");
3319 msgbox.exec();
3320 }
3321 }
3322
3323 // G4cout << "RENAMING to " << newname << G4endl;
3324 renameViewPt(newname);
3325 listitem->setText(QString(newname));
3326 AuxWindowDialog->lineEdit->setText(newname);
3327 // if (currentState == VIEWPOINT)
3328 // scheduleRedraw();
3329
3330 delete[] newname;
3331}
3332
3333// Renames currently selected viewpoint.
3334
3336{
3337 std::size_t idx = 0, end, pos;
3338 std::size_t size = viewPtList.size();
3339 std::string line, newName;
3340 fileIn.open(fileName.c_str());
3341
3342 newName = vpName;
3343 while ((int) newName.size() < MAX_VP_NAME)
3344 newName += " ";
3345
3346 getline(fileIn, line);
3347 pos = fileIn.tellg();
3348 while (getline(fileIn, line)) {
3349 end = line.find_last_not_of(' ');
3350 line = line.substr(0, end + 1);
3351 if (!strcmp(line.c_str(), curViewPtName)) {
3352 fileOut.seekp(pos);
3353 fileOut << newName;
3354 fileOut.seekp(0, std::ios::end); // Set the file pointer to the end of the file
3355 break;
3356 }
3357 while (line.size())
3358 getline(fileIn, line);
3359 pos = fileIn.tellg();
3360 }
3361
3362 fileIn.close();
3363 fileIn.clear();
3364
3365 while (idx < size) {
3366 if (!strcmp(viewPtList[idx].viewPtName, curViewPtName)) {
3367 strcpy(viewPtList[idx].viewPtName, vpName);
3368 break;
3369 }
3370 idx++;
3371 }
3372}
3373
3374
3375void G4OpenInventorQtExaminerViewer::SortBookmarksCB()
3376{
3377 // FWJ NOTE: Xt version of this does not work (does nothing)
3378
3379 // G4cout << "Sort Button: SortBookmarks CALLBACK" << G4endl;
3380
3381 // FWJ List for sorting
3382 // The dialog list and bookmark file will be rewritten.
3383 // Simpler to populate this list from the data structure.
3384
3385 std::vector<std::string> charList;
3386
3387 if (viewPtList.size() < 2) return;
3388
3389 // Get current entries from the list
3390
3391 for (int i = 0; i < (int)viewPtList.size(); i++) {
3392
3393 charList.push_back(viewPtList[i].viewPtName);
3394 // G4cout << " Pushed " << i << " " << charList[i] << G4endl;
3395 }
3396
3397 std::sort(charList.begin(), charList.end());
3398
3399 // FWJ POPULATE the new dialog list
3400 // G4cout << " Populating Bookmark listWidget..." << G4endl;
3401 AuxWindowDialog->listWidget->clear();
3402
3403 for (int i = 0; i < (int)viewPtList.size(); i++) {
3404 // viewPtIdx has to be changed to account for a different order in viewPtList
3405 if (!strcmp(charList[i].c_str(), curViewPtName))
3406 viewPtIdx = i;
3407 new QListWidgetItem(charList[i].c_str(), AuxWindowDialog->listWidget);
3408
3409 }
3410
3411 sortViewPts(charList);
3412
3413}
3414
3415// Rewrites entire viewpoint file with sorted viewpoints.
3416
3417void G4OpenInventorQtExaminerViewer::sortViewPts(std::vector<std::string> sortedViewPts)
3418{
3419 SbVec3f axis;
3420 float x, y, z, angle;
3421 std::size_t sortIdx = 0, unsortIdx = 0;
3422
3423 if (fileOut.is_open())
3424 fileOut.close();
3425
3426 fileOut.open(fileName.c_str()); // Erase current viewpoint file
3427
3429
3430 std::size_t size = sortedViewPts.size();
3431 while (sortIdx < size) {
3432 while (strcmp(sortedViewPts[sortIdx].c_str(),
3433 viewPtList[unsortIdx].viewPtName))
3434 unsortIdx++;
3435
3436 std::string vpName = viewPtList[unsortIdx].viewPtName;
3437
3438 while ((int) vpName.size() < MAX_VP_NAME)
3439 vpName += " ";
3440 fileOut << vpName << std::endl;
3441 viewPtList[unsortIdx].position.getValue(x, y, z);
3442 fileOut << x << " " << y << " " << z << std::endl;
3443
3444 // Reusing x, y and z for storing the axis
3445 viewPtList[unsortIdx].orientation.getValue(axis, angle);
3446 axis.getValue(x, y, z);
3447 fileOut << x << " " << y << " " << z << " " << angle << std::endl;
3448
3449 fileOut << viewPtList[unsortIdx].camType << " "
3450 << viewPtList[unsortIdx].height << std::endl;
3451 fileOut << viewPtList[unsortIdx].focalDistance << " ";
3452
3453 fileOut << viewPtList[unsortIdx].nearDistance << " ";
3454
3455 fileOut << viewPtList[unsortIdx].farDistance << std::endl;
3456
3457 fileOut << viewPtList[unsortIdx].viewportMapping << " ";
3458 fileOut << viewPtList[unsortIdx].aspectRatio << "\n" << std::endl;
3459 fileOut.flush();
3460
3461 unsortIdx = 0;
3462 sortIdx++;
3463 }
3464}
3465
3466// Needed to implement mouse wheel zoom direction change.
3467// Does not work with MacOS trackpad: use Coin3d default handler.
3468// Emulating private method SoGuiFullViewerP::zoom()
3469#ifndef __APPLE__
3470void
3472{
3473 float multiplicator = float(std::exp(diffvalue));
3474 SoCamera *cam = getCamera();
3475
3476 if (cam->isOfType(SoPerspectiveCamera::getClassTypeId())) {
3477 const float oldfocaldist = cam->focalDistance.getValue();
3478 const float newfocaldist = oldfocaldist * multiplicator;
3479
3480 SbVec3f direction;
3481 cam->orientation.getValue().multVec(SbVec3f(0, 0, -1), direction);
3482
3483 const SbVec3f oldpos = cam->position.getValue();
3484 const SbVec3f newpos = oldpos + (newfocaldist - oldfocaldist) * -direction;
3485 cam->position = newpos;
3486 cam->focalDistance = newfocaldist;
3487 } else if (cam->isOfType(SoOrthographicCamera::getClassTypeId())) {
3488 SoOrthographicCamera * oc = (SoOrthographicCamera *)cam;
3489 oc->height = oc->height.getValue() * multiplicator;
3490 }
3491}
3492#endif
3493
3494// Handling mouse and keyboard events
3495
3496SbBool
3498{
3499
3500 // FWJ DEBUG
3501 // G4cout << "processSoEvent ############" << ++processSoEventCount << G4endl;
3502
3503 SoCamera *cam = getCamera();
3504 const SoType type(ev->getTypeId());
3505
3506// Needed to implement mouse wheel zoom direction change.
3507// Does not work with MacOS trackpad: use Coin3d default handler.
3508#ifndef __APPLE__
3509 if (type.isDerivedFrom(SoMouseButtonEvent::getClassTypeId())) {
3510 SoMouseButtonEvent * me = (SoMouseButtonEvent *) ev;
3511
3512 // if (currentState == ANIMATION || currentState == REVERSED_ANIMATION
3513 // || currentState == PAUSED_ANIMATION) {
3514
3515 switch (me->getButton()) {
3516
3517 case SoMouseButtonEvent::BUTTON4: // Scroll wheel up
3518 if (me->getState() == SoButtonEvent::DOWN) {
3519 // G4cout << "SCROLL WHEEL UP" << G4endl;
3520 zoom(-0.1f);
3521 return TRUE;
3522 }
3523 break;
3524
3525 case SoMouseButtonEvent::BUTTON5: // Scroll wheel down
3526 if (me->getState() == SoButtonEvent::DOWN) {
3527 // G4cout << "SCROLL WHEEL DOWN" << G4endl;
3528 zoom(0.1f);
3529 return TRUE;
3530 }
3531 break;
3532
3533 default:
3534 break;
3535 }
3536 // }
3537 if (currentState == GENERAL) {
3538
3539 }
3540 }
3541#endif
3542
3543 if (type.isDerivedFrom(SoKeyboardEvent::getClassTypeId())) {
3544 SoKeyboardEvent* ke = (SoKeyboardEvent*)ev;
3545
3546 if (SoKeyboardEvent::isKeyPressEvent(ev, ke->getKey())) {
3547 switch (ke->getKey()) {
3548 case SoKeyboardEvent::E:
3549 if (externalQtApp) {
3550 // G4cout << "E KEY PRESSED" << G4endl;
3551 return TRUE;
3552 } else {
3553 G4cout <<
3554 "E KEY PRESSED, EXITING OIQT VIEWER SECONDARY LOOP" <<
3555 G4endl;
3556 SoQt::exitMainLoop();
3557 // escapeCallback();
3558 return TRUE;
3559 }
3560 case SoKeyboardEvent::LEFT_SHIFT:
3561 this->lshiftdown = true;
3562 return TRUE;
3563 case SoKeyboardEvent::RIGHT_SHIFT:
3564 this->rshiftdown = true;
3565 return TRUE;
3566 case SoKeyboardEvent::LEFT_CONTROL:
3567 this->lctrldown = true;
3568 return TRUE;
3569 case SoKeyboardEvent::RIGHT_CONTROL:
3570 this->rctrldown = true;
3571 return TRUE;
3572 case SoKeyboardEvent::SPACE:
3573 if (currentState == ANIMATION
3577 if (animateSensor->isScheduled())
3578 animateSensor->unschedule();
3579 return TRUE;
3580 } else if (currentState == PAUSED_ANIMATION) {
3581 if (maxSpeed) {
3582 if ((beforePausing == ANIMATION
3584 < (int) refParticleTrajectory.size() - 1)
3586 && refParticleIdx > 0)) {
3589 }
3590 }
3591 return TRUE;
3592 }
3593 break;
3594 case SoKeyboardEvent::ESCAPE:
3595 if (currentState == ANIMATION
3598
3599 if (animateSensor->isScheduled())
3600 animateSensor->unschedule();
3603 setSuperimpositionEnabled(superimposition, FALSE);
3604 maxSpeed = 0.0f;
3605 step = 1;
3606
3607 scheduleRedraw();
3608 if (currentState == VIEWPOINT) {
3609 setSuperimpositionEnabled(superimposition, TRUE);
3610 axisSwitch->whichChild.setValue(SO_SWITCH_NONE);
3611 animSpeedOutlineSwitch->whichChild.setValue(
3612 SO_SWITCH_NONE);
3613 animSpeedSwitch->whichChild.setValue(SO_SWITCH_NONE);
3614
3615 scheduleRedraw();
3616 }
3617 restoreCamera();
3618 return TRUE;
3619 }
3620 break;
3621 case SoKeyboardEvent::DELETE:
3622 if (viewPtList.size()
3623 && (currentState != ANIMATION
3626 // FWJ IMPLEMENT LATER
3627 // String dialogName = (char *) "Delete Viewpoint";
3628 // std::string msg = "Are you sure you want to delete current viewpoint?";
3629 // warningMsgDialog(msg, dialogName, deleteViewPtCB);
3630 return TRUE;
3631 }
3632 break;
3633 case SoKeyboardEvent::LEFT_ARROW:
3634 switch (currentState) {
3635 case BEAMLINE:
3636 if ((this->lshiftdown) || (this->rshiftdown)) {
3638 moveCamera();
3639 }
3640 else if ((this->lctrldown) || (this->rctrldown)) {
3641 if (SoQtExaminerViewer::isAnimating())
3642 stopAnimating();
3645 animateBtwPtsPeriod = 0.08f;
3646
3647 SbVec3f tmp = camDir;
3648 tmp.negate();
3649 rotAxis = tmp;
3650
3651 rotCnt = ROT_CNT;
3652 moveCamera(); // To make sure camera is perpendicular to the beamline
3653 rotateCamera();
3654 }
3655 else {
3656 if (SoQtExaminerViewer::isAnimating())
3657 stopAnimating();
3660 animateBtwPtsPeriod = 0.08f;
3661
3662 SbVec3f tmp = camUpVec;
3663 tmp.negate();
3664 rotAxis = tmp;
3665
3666 rotCnt = ROT_CNT;
3667 moveCamera(); // To make sure camera is perpendicular to the beamline
3668 rotateCamera();
3669
3670 }
3671 return TRUE;
3672
3673 case ANIMATION:
3674 case REVERSED_ANIMATION:
3675 left_right -= 1.5f;
3676 return TRUE;
3677 case PAUSED_ANIMATION:
3678 left_right -= 1.5f;
3680 cam->position = myCam->position;
3681 return TRUE;
3682 case GENERAL:
3683 case VIEWPOINT:
3684 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3685 // Using this allows us to look around without
3686 // changing the camera parameters (camDir, camUpVec)
3687 this->bottomWheelMotion(
3688 this->getBottomWheelValue() + 0.1f);
3689
3690 return TRUE;
3691 }
3692 break;
3693 case ROTATING:
3694 // For this state, let the keyboard event
3695 // be handled by superclass
3696 break;
3697 default:
3698 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3699 "Unhandled viewer state");
3700 break;
3701 }
3702 break;
3703
3704 case SoKeyboardEvent::RIGHT_ARROW:
3705 switch(currentState) {
3706 case BEAMLINE:
3707 if ((this->lshiftdown) || (this->rshiftdown)) {
3709 moveCamera();
3710 }
3711 else if ((this->lctrldown) || (this->rctrldown)) {
3712 if (SoQtExaminerViewer::isAnimating())
3713 stopAnimating();
3716 animateBtwPtsPeriod = 0.08f;
3717
3718 rotAxis = camDir;
3719
3720 rotCnt = ROT_CNT;
3721 moveCamera(); // To make sure camera is perpendicular to the beamline
3722 rotateCamera();
3723 }
3724 else{
3725 if (SoQtExaminerViewer::isAnimating())
3726 stopAnimating();
3729 animateBtwPtsPeriod = 0.08f;
3730
3731 rotAxis = camUpVec;
3732
3733 rotCnt = ROT_CNT;
3734 moveCamera(); // To make sure camera is perpendicular to the beamline
3735 rotateCamera();
3736 }
3737 return TRUE;
3738
3739 case ANIMATION:
3740 case REVERSED_ANIMATION:
3741 left_right += 1.5f;
3742 return TRUE;
3743 case PAUSED_ANIMATION:
3744 left_right += 1.5f;
3746 cam->position = myCam->position;
3747 return TRUE;
3748 case GENERAL:
3749 case VIEWPOINT:
3750 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3751 // Using this allows us to look around without
3752 // changing the camera parameters (camDir, camUpVec)
3753 this->bottomWheelMotion(
3754 this->getBottomWheelValue() - 0.1f);
3755 return TRUE;
3756 }
3757 break;
3758 case ROTATING:
3759 // For this state, let the keyboard event
3760 // be handled by superclass
3761 break;
3762 default:
3763 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3764 "Unhandled viewer state");
3765 break;
3766 }
3767 break;
3768
3769 case SoKeyboardEvent::DOWN_ARROW:
3770 switch(currentState) {
3771 case BEAMLINE:
3772
3773 if ((this->lshiftdown) || (this->rshiftdown)) {
3775 moveCamera();
3776 }
3777 else{
3778 if (SoQtExaminerViewer::isAnimating())
3779 stopAnimating();
3782 animateBtwPtsPeriod = 0.08f;
3783
3784 rotAxis = camDir.cross(camUpVec);
3785
3786 rotCnt = ROT_CNT;
3787 moveCamera(); // To make sure camera is perpendicular to the beamline
3788 rotateCamera();
3789
3790 }
3791 return TRUE;
3792
3793 case ANIMATION:
3794 case REVERSED_ANIMATION:
3795 up_down -= 1.5f;
3796 return TRUE;
3797 case PAUSED_ANIMATION:
3798 up_down -= 1.5f;
3800 cam->position = myCam->position;
3801 return TRUE;
3802 case GENERAL:
3803 case VIEWPOINT:
3804 // Using this allows us to look around without
3805 // changing the camera parameters (camDir, camUpVec)
3806 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3807 this->leftWheelMotion(this->getLeftWheelValue() - 0.1f);
3808 return TRUE;
3809 }
3810 break;
3811 case ROTATING:
3812 // For this state, let the keyboard event
3813 // be handled by superclass
3814 break;
3815 default:
3816 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3817 "Unhandled viewer state");
3818 break;
3819 }
3820 break;
3821
3822 case SoKeyboardEvent::UP_ARROW:
3823 switch(currentState) {
3824 case BEAMLINE:
3825 if ((this->lshiftdown) || (this->rshiftdown)) {
3827 moveCamera();
3828 }
3829 else{
3830 if (SoQtExaminerViewer::isAnimating())
3831 stopAnimating();
3834 animateBtwPtsPeriod = 0.08f;
3835
3836 rotAxis = camUpVec.cross(camDir);
3837
3838 rotCnt = ROT_CNT;
3839 moveCamera();
3840
3841 rotateCamera();
3842
3843
3844 }
3845 return TRUE;
3846 case ANIMATION:
3847 case REVERSED_ANIMATION:
3848 up_down += 1.5f;
3849 return TRUE;
3850 case PAUSED_ANIMATION:
3851 up_down += 1.5f;
3853 cam->position = myCam->position;
3854 return TRUE;
3855 case GENERAL:
3856 case VIEWPOINT:
3857 // Using this allows us to look around without
3858 // changing the camera parameters (camDir, camUpVec)
3859 if ((!this->lshiftdown) && (!this->rshiftdown)) {
3860 this->leftWheelMotion(this->getLeftWheelValue() + 0.1f);
3861 return TRUE;
3862 }
3863 break;
3864 case ROTATING:
3865 // For this state, let the keyboard event
3866 // be handled by superclass
3867 break;
3868 default:
3869 SoDebugError::post("G4OpenInventorQtExaminerViewer::processSoEvent",
3870 "Unhandled viewer state");
3871 break;
3872 }
3873 break;
3874
3875 case SoKeyboardEvent::PAGE_UP:
3876 switch(currentState) {
3877 case BEAMLINE:
3878 if (step < (int) refParticleTrajectory.size() / 5) // Magic number
3879 step++;
3880 return TRUE;
3881 case ANIMATION:
3882 incSpeed();
3884 if (maxSpeed > 0.8)
3886 scheduleRedraw();
3887
3888 return TRUE;
3889 case REVERSED_ANIMATION:
3890 if(!animateSensor->isScheduled()) {
3892 if (refParticleIdx
3893 < (int) refParticleTrajectory.size() - 1) {
3896 scheduleRedraw();
3898 }
3899 }
3900 else{
3902 decSpeed();
3903 scheduleRedraw();
3904 }
3905 return TRUE;
3906 case PAUSED_ANIMATION:
3908 if (maxSpeed > 0.8)
3910
3911 if (beforePausing == ANIMATION) {
3912 incSpeed();
3913 } else {
3914 decSpeed();
3917 }
3918
3919 scheduleRedraw();
3920 return TRUE;
3921 default: //fall through
3922 break;
3923 }
3924 break;
3925
3926 case SoKeyboardEvent::PAGE_DOWN:
3927 switch(currentState) {
3928 case BEAMLINE:
3929 if (step > 1)
3930 step--;
3931 return TRUE;
3932 case ANIMATION:
3933 if(!animateSensor->isScheduled()) {
3935 if (refParticleIdx > 1) {
3938 scheduleRedraw();
3940 }
3941 }
3942 else{
3944 decSpeed();
3945 scheduleRedraw();
3946 }
3947 return TRUE;
3948 case REVERSED_ANIMATION:
3949 incSpeed();
3951 if (maxSpeed < -0.8)
3953 scheduleRedraw();
3954 return TRUE;
3955 case PAUSED_ANIMATION:
3957 if (maxSpeed < -0.8)
3960 incSpeed();
3961 } else {
3962 decSpeed();
3965 }
3966 scheduleRedraw();
3967 return TRUE;
3968 default:
3969 //fall through
3970 break;
3971 }
3972 break;
3973
3974 // FROM XT VIEWER
3975 // case SoKeyboardEvent::E:
3976 // this->escapeCallback(this->examinerObject);
3977 // break;
3978
3979 default:
3980 break; // To get rid of compiler warnings
3981 }
3982 }
3983 if (SoKeyboardEvent::isKeyReleaseEvent(ev, ke->getKey())) {
3984 switch (ke->getKey()) {
3985 case SoKeyboardEvent::LEFT_SHIFT:
3986 this->lshiftdown = false;
3987 return TRUE;
3988 case SoKeyboardEvent::RIGHT_SHIFT:
3989 this->rshiftdown = false;
3990 return TRUE;
3991 case SoKeyboardEvent::LEFT_CONTROL:
3992 this->lctrldown = false;
3993 return TRUE;
3994 case SoKeyboardEvent::RIGHT_CONTROL:
3995 this->rctrldown = false;
3996 return TRUE;
3997 default:
3998 break;
3999 }
4000 }
4001 }
4002
4003 // Pass the event on to the viewer
4004 // Need some checks here as in Xt viewer?
4005
4007 || currentState == ROTATING)
4008 return FALSE;
4009 else
4010 return SoQtExaminerViewer::processSoEvent(ev);
4011}
4012
4013
4014// REMAINDER OF MENU BAR FUNCTIONS...
4015
4016
4017void G4OpenInventorQtExaminerViewer::FileLoadSceneGraphCB()
4018{
4019 // G4cout << "File: Load scene graph CALLBACK" << G4endl;
4020
4021 QFileDialog filedialog(getParentWidget(), tr("Load Scene Graph"));
4022 filedialog.setFileMode(QFileDialog::AnyFile);
4023 filedialog.setFont(*font);
4024 if (!filedialog.exec()) return;
4025 QStringList filenameinlist = filedialog.selectedFiles();
4026 QString filenamein = filenameinlist[0];
4027
4028 // G4cout << "Entered file name is " << qPrintable(filenamein) << G4endl;
4029
4030 char* filename = new char[filenamein.size()+1];
4031 filename = strdup(qPrintable(filenamein));
4032 // G4cout << "char[] file name is " << filename << G4endl;
4033
4034 SoInput sceneInput;
4035
4036 if (sceneInput.openFile(filename)) {
4037 // Read the whole file into the database
4038 newSceneGraph = SoDB::readAll(&sceneInput);
4039 if (newSceneGraph == NULL) {
4040 QMessageBox msgbox;
4041 msgbox.setFont(*font);
4042 QString messagetxt = "Error reading scene graph file ";
4043 messagetxt.append(filenamein);
4044 msgbox.setText(messagetxt);
4045 msgbox.exec();
4046 sceneInput.closeFile();
4047 return;
4048 }
4049 } else {
4050 QMessageBox msgbox;
4051 msgbox.setFont(*font);
4052 QString messagetxt = "Error opening scene graph file ";
4053 messagetxt.append(filenamein);
4054 msgbox.setText(messagetxt);
4055 msgbox.exec();
4056 return;
4057 }
4058
4059 SoSeparator* root = (SoSeparator*)getSceneGraph();
4060 root->unref();
4061 newSceneGraph->ref();
4062 setSceneGraph(newSceneGraph);
4063}
4064
4065void G4OpenInventorQtExaminerViewer::FileSaveSceneGraphCB()
4066{
4067 // G4cout << "File: Save scene graph CALLBACK" << G4endl;
4068
4069 QFileDialog filedialog(getParentWidget(), tr("Save scene graph"));
4070 filedialog.setFileMode(QFileDialog::AnyFile);
4071 // To enable confirmation of overwriting
4072 filedialog.setAcceptMode(QFileDialog::AcceptSave);
4073 filedialog.setFont(*font);
4074 if (!filedialog.exec()) return;
4075 QStringList filenameinlist = filedialog.selectedFiles();
4076 QString filenamein = filenameinlist[0];
4077
4078 // G4cout << "Entered file name is " << qPrintable(filenamein) << G4endl;
4079
4080 char* filename = new char[filenamein.size()+1];
4081 filename = strdup(qPrintable(filenamein));
4082 // G4cout << "char[] file name is " << filename << G4endl;
4083
4084 SoWriteAction writeAction;
4085 SoSeparator* root = (SoSeparator*)getSceneGraph();
4086
4087 SoOutput* out = writeAction.getOutput();
4088
4089 if (out->openFile(filename)) {
4090 out->setBinary(FALSE);
4091 writeAction.apply(root);
4092 out->closeFile();
4093 } else {
4094 QMessageBox msgbox;
4095 msgbox.setFont(*font);
4096 QString messagetxt = "Error opening file ";
4097 messagetxt.append(filenamein);
4098 msgbox.setText(messagetxt);
4099 msgbox.exec();
4100 }
4101}
4102
4103
4104void G4OpenInventorQtExaminerViewer::HelpControlsCB()
4105{
4106 // G4cout << "Help: Help Controls CALLBACK" << G4endl;
4107 helpmsgbox->show();
4108}
4109
4110
4112{
4113 viewer = vwr;
4114}
4115
4117{;}
4118
4120{
4121 if (requestedState == G4State_EventProc) viewer->newEvents = true;
4122 return true;
4123}
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:57
G4GLOB_DLL std::ostream G4cout
#define M_PI
Definition: SbMath.h:33
const std::vector< const std::vector< G4AttValue > * > & GetAttValues() const
Definition: G4AttHolder.hh:61
const std::vector< const std::map< G4String, G4AttDef > * > & GetAttDefs() const
Definition: G4AttHolder.hh:63
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:1860
QTabWidget * GetViewerTabWidget()
Definition: G4UIQt.hh:167
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
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)