Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/visualization/management/src/G4VisManager.cc

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

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 // 
 28 // GEANT4 Visualization Manager - John Allison 02/Jan/1996.
 29 // Michael Kelsey  31 Jan 2019 -- Add new command for electric field
 30 
 31 #include "G4VisManager.hh"
 32 
 33 #include "G4VisCommands.hh"
 34 #include "G4VisCommandsCompound.hh"
 35 #include "G4VisCommandsGeometry.hh"
 36 #include "G4VisCommandsGeometrySet.hh"
 37 #include "G4VisCommandsMultithreading.hh"
 38 #include "G4VisCommandsSet.hh"
 39 #include "G4VisCommandsScene.hh"
 40 #include "G4VisCommandsSceneAdd.hh"
 41 #include "G4VisCommandsPlotter.hh"
 42 #include "G4VisCommandsSceneHandler.hh"
 43 #include "G4VisCommandsTouchable.hh"
 44 #include "G4VisCommandsTouchableSet.hh"
 45 #include "G4VisCommandsViewer.hh"
 46 #include "G4VisCommandsViewerDefault.hh"
 47 #include "G4VisCommandsViewerSet.hh"
 48 #include "G4UImanager.hh"
 49 #include "G4VisStateDependent.hh"
 50 #include "G4UIdirectory.hh"
 51 #include "G4VGraphicsSystem.hh"
 52 #include "G4VSceneHandler.hh"
 53 #include "G4VViewer.hh"
 54 #include "G4VPhysicalVolume.hh"
 55 #include "G4LogicalVolume.hh"
 56 #include "G4VSolid.hh"
 57 #include "G4Vector3D.hh"
 58 #include "G4Point3D.hh"
 59 #include "G4RotationMatrix.hh"
 60 #include "G4Polyline.hh"
 61 #include "G4Polyhedron.hh"
 62 #include "G4NullModel.hh"
 63 #include "G4ModelingParameters.hh"
 64 #include "G4TransportationManager.hh"
 65 #include "G4VisCommandModelCreate.hh"
 66 #include "G4VisCommandsListManager.hh"
 67 #include "G4VisModelManager.hh"
 68 #include "G4VModelFactory.hh"
 69 #include "G4VisFilterManager.hh"
 70 #include "G4VTrajectoryModel.hh"
 71 #include "G4TrajectoryDrawByCharge.hh"
 72 #include "Randomize.hh"
 73 #include "G4RunManager.hh"
 74 #include "G4RunManagerFactory.hh"
 75 #include "G4EventManager.hh"
 76 #include "G4Run.hh"
 77 #include "G4Event.hh"
 78 #include <map>
 79 #include <set>
 80 #include <vector>
 81 #include <sstream>
 82 
 83 #include "G4Threading.hh"
 84 #include "G4AutoLock.hh"
 85 #include "G4GeometryWorkspace.hh" // no_geant4_module_check(!G4MULTITHREADED)
 86 #include "G4SolidsWorkspace.hh"
 87 #include <deque>
 88 #include <typeinfo>
 89 #include <chrono>
 90 #include <thread>
 91 #include <utility>
 92 
 93 #define G4warn G4cout
 94 
 95 // Local threading related variables
 96 namespace {
 97   G4bool mtRunInProgress = false;
 98   std::deque<const G4Event*> mtVisEventQueue;
 99   G4Thread* mtVisSubThread = 0;
100   [[maybe_unused]] G4Mutex mtVisSubThreadMutex = G4MUTEX_INITIALIZER;
101   //  G4Mutex visBeginOfRunMutex = G4MUTEX_INITIALIZER;
102   //  G4Mutex visBeginOfEventMutex = G4MUTEX_INITIALIZER;
103   G4Mutex visEndOfEventMutex = G4MUTEX_INITIALIZER;
104   //  G4Mutex visEndOfRunMutex = G4MUTEX_INITIALIZER;
105   G4bool isSubEventRunManagerType = false;
106   G4bool isValidViewForRun = false;
107   G4bool isFakeRun = false;
108 }
109 
110 // Class statics
111 G4VisManager* G4VisManager::fpInstance = 0;
112 
113 G4VisManager::Verbosity G4VisManager::fVerbosity = G4VisManager::warnings;
114 
115 G4VisManager::G4VisManager (const G4String& verbosityString)
116 : fVerbose         (1)
117 , fDefaultGraphicsSystemName("OGL")          // Override in G4VisExecutive
118 , fDefaultXGeometryString   ("600x600-0+0")  // Override in G4VisExecutive
119 , fDefaultGraphicsSystemBasis ("G4VisManager initialisation")
120 , fDefaultXGeometryStringBasis("G4VisManager initialisation")
121 , fInitialised     (false)
122 , fpGraphicsSystem (0)
123 , fpScene          (0)
124 , fpSceneHandler   (0)
125 , fpViewer         (0)
126 , fpStateDependent (0)
127 , fEventRefreshing          (false)
128 , fTransientsDrawnThisRun   (false)
129 , fTransientsDrawnThisEvent (false)
130 , fNoOfEventsDrawnThisRun   (0)
131 , fNKeepForPostProcessingRequests (0)
132 , fNKeepTheEventRequests    (0)
133 , fEventKeepingSuspended    (false)
134 , fDrawEventOnlyIfToBeKept  (false)
135 , fpRequestedEvent          (0)
136 , fReviewingKeptEvents      (false)
137 , fAbortReviewKeptEvents    (false)
138 , fReviewingPlots           (false)
139 , fAbortReviewPlots         (false)
140 , fIsDrawGroup              (false)
141 , fDrawGroupNestingDepth    (0)
142 , fIgnoreStateChanges       (false)
143 , fMaxEventQueueSize        (100)
144 , fWaitOnEventQueueFull     (true)
145 // All other objects use default constructors.
146 {
147   fpTrajDrawModelMgr = new G4VisModelManager<G4VTrajectoryModel>("/vis/modeling/trajectories");
148   fpTrajFilterMgr = new G4VisFilterManager<G4VTrajectory>("/vis/filtering/trajectories");
149   fpHitFilterMgr = new G4VisFilterManager<G4VHit>("/vis/filtering/hits");
150   fpDigiFilterMgr = new G4VisFilterManager<G4VDigi>("/vis/filtering/digi");
151 
152   VerbosityGuidanceStrings.push_back
153     ("Simple graded message scheme - digit or string (1st character defines):");
154   VerbosityGuidanceStrings.push_back
155     ("  0) quiet,         // Nothing is printed.");
156   VerbosityGuidanceStrings.push_back
157     ("  1) startup,       // Startup and endup messages are printed...");
158   VerbosityGuidanceStrings.push_back
159     ("  2) errors,        // ...and errors...");
160   VerbosityGuidanceStrings.push_back
161     ("  3) warnings,      // ...and warnings...");
162   VerbosityGuidanceStrings.push_back
163     ("  4) confirmations, // ...and confirming messages...");
164   VerbosityGuidanceStrings.push_back
165     ("  5) parameters,    // ...and parameters of scenes and views...");
166   VerbosityGuidanceStrings.push_back
167     ("  6) all            // ...and everything available.");
168 
169   if (fpInstance) {
170     G4Exception
171       ("G4VisManager::G4VisManager",
172        "visman0001", FatalException,
173        "Attempt to Construct more than one VisManager");
174   }
175 
176   fpInstance = this;
177   SetConcreteInstance(this);
178 
179   fpStateDependent = new G4VisStateDependent (this);
180   // No need to delete this; G4StateManager does this.
181 
182   fVerbosity = GetVerbosityValue(verbosityString);
183   if (fVerbosity >= startup) {
184       G4cout
185   << "Visualization Manager instantiating with verbosity \""
186   << VerbosityString(fVerbosity)
187   << "\"..." << G4endl;
188   }
189 
190   // Note: The specific graphics systems must be instantiated in a
191   // higher level library to avoid circular dependencies.  Also,
192   // some specifically need additional external libararies that the
193   // user must supply.  Therefore we ask the user to implement
194   // RegisterGraphicsSystems() and RegisterModelFactories()
195   // in a subclass.  We have to wait for the subclass to instantiate
196   // so RegisterGraphicsSystems() cannot be called from this
197   // constructor; it is called from Initialise().  So we ask the
198   // user:
199   //   (a) to write a subclass and implement  RegisterGraphicsSystems()
200   //       and RegisterModelFactories().  See
201   //       visualization/include/G4VisExecutive.hh/icc as an example.
202   //   (b) instantiate the subclass.
203   //   (c) invoke the Initialise() method of the subclass.
204   // For example:
205   //   ...
206   //   // Instantiate and initialise Visualization Manager.
207   //   G4VisManager* visManager = new G4VisExecutive;
208   //   visManager -> SetVerboseLevel (Verbose);
209   //   visManager -> Initialise ();
210   //   // (Don't forget to delete visManager;)
211   //   ...
212 
213   // Make top level command directory...
214   // vis commands should *not* be broadcast to workers
215   G4bool propagateToWorkers;
216   auto directory = new G4UIdirectory ("/vis/",propagateToWorkers=false);
217   directory -> SetGuidance ("Visualization commands.");
218   // Request commands in name order
219   directory -> Sort();  // Ordering propagates to sub-directories
220   fDirectoryList.push_back (directory);
221 
222   // Instantiate *basic* top level commands so that they can be used
223   // immediately after instantiation of the vis manager.  Other top
224   // level and lower level commands are instantiated later in
225   // RegisterMessengers.
226   G4VVisCommand::SetVisManager (this);  // Sets shared pointer
227   RegisterMessenger(new G4VisCommandVerbose);
228   RegisterMessenger(new G4VisCommandInitialize);
229 }
230 
231 G4VisManager::~G4VisManager()
232 {
233   G4UImanager* UImanager = G4UImanager::GetUIpointer();
234   UImanager->SetCoutDestination(nullptr);
235   std::size_t i;
236   for (i = 0; i < fSceneList.size (); ++i) {
237     delete fSceneList[i];
238   }
239   for (i = 0; i < fAvailableSceneHandlers.size (); ++i) {
240   if (fAvailableSceneHandlers[i] != NULL) {
241     delete fAvailableSceneHandlers[i];
242   }
243   }
244   for (i = 0; i < fAvailableGraphicsSystems.size (); ++i) {
245     if (fAvailableGraphicsSystems[i]) {
246       delete fAvailableGraphicsSystems[i];
247     }
248   }
249   if (fVerbosity >= startup) {
250     G4cout << "Graphics systems deleted." << G4endl;
251     G4cout << "Visualization Manager deleting..." << G4endl;
252   }
253   for (i = 0; i < fMessengerList.size (); ++i) {
254     delete fMessengerList[i];
255   }
256   for (i = 0; i < fDirectoryList.size (); ++i) {
257     delete fDirectoryList[i];
258   }
259 
260   delete fpDigiFilterMgr;
261   delete fpHitFilterMgr;
262   delete fpTrajFilterMgr;
263   delete fpTrajDrawModelMgr;
264   fpInstance = 0;
265 }
266 
267 G4VisManager* G4VisManager::GetInstance () {
268   if (!fpInstance) {
269     G4Exception
270       ("G4VisManager::GetInstance",
271        "visman0002", FatalException, "VisManager not yet instantiated");
272   }
273   return fpInstance;
274 }
275 
276 void G4VisManager::Initialise () {
277 
278   if (fInitialised && fVerbosity >= warnings) {
279     G4warn << "WARNING: G4VisManager::Initialise: already initialised."
280      << G4endl;
281     return;
282   }
283 
284   if (fVerbosity >= startup) {
285     G4cout << "Visualization Manager initialising..." << G4endl;
286   }
287 
288   if (fVerbosity >= parameters) {
289     G4cout <<
290       "\nYou have instantiated your own Visualization Manager, inheriting"
291       "\n  G4VisManager and implementing RegisterGraphicsSystems(), in which"
292       "\n  you should, normally, instantiate drivers which do not need"
293       "\n  external packages or libraries, and, optionally, drivers under"
294       "\n  control of environment variables."
295       "\n  Also you should implement RegisterModelFactories()."
296       "\n  See visualization/management/include/G4VisExecutive.hh/icc, for example."
297       "\n  In your main() you will have something like:"
298       "\n    G4VisManager* visManager = new G4VisExecutive;"
299       "\n    visManager -> SetVerboseLevel (Verbose);"
300       "\n    visManager -> Initialize ();"
301       "\n  (Don't forget to delete visManager;)"
302       "\n"
303    << G4endl;
304   }
305 
306   if (fVerbosity >= startup) {
307     G4cout << "Registering graphics systems..." << G4endl;
308   }
309 
310   RegisterGraphicsSystems ();
311 
312   if (fVerbosity >= startup) {
313     G4cout <<
314       "\nYou have successfully registered the following graphics systems."
315    << G4endl;
316     PrintAvailableGraphicsSystems (fVerbosity);
317     G4cout << G4endl;
318   }
319 
320   // Make command directories for commands instantiated in the
321   // modeling subcategory...
322   G4UIcommand* directory;
323   directory = new G4UIdirectory ("/vis/modeling/");
324   directory -> SetGuidance ("Modeling commands.");
325   fDirectoryList.push_back (directory);
326   directory = new G4UIdirectory ("/vis/modeling/trajectories/");
327   directory -> SetGuidance ("Trajectory model commands.");
328   fDirectoryList.push_back (directory);
329   directory = new G4UIdirectory ("/vis/modeling/trajectories/create/");
330   directory -> SetGuidance ("Create trajectory models and messengers.");
331   fDirectoryList.push_back (directory);
332 
333   // Filtering command directory
334   directory = new G4UIdirectory ("/vis/filtering/");
335   directory -> SetGuidance ("Filtering commands.");
336   fDirectoryList.push_back (directory);
337   directory = new G4UIdirectory ("/vis/filtering/trajectories/");
338   directory -> SetGuidance ("Trajectory filtering commands.");
339   fDirectoryList.push_back (directory);
340   directory = new G4UIdirectory ("/vis/filtering/trajectories/create/");
341   directory -> SetGuidance ("Create trajectory filters and messengers.");
342   fDirectoryList.push_back (directory);
343   directory = new G4UIdirectory ("/vis/filtering/hits/");
344   directory -> SetGuidance ("Hit filtering commands.");
345   fDirectoryList.push_back (directory);
346   directory = new G4UIdirectory ("/vis/filtering/hits/create/");
347   directory -> SetGuidance ("Create hit filters and messengers.");
348   fDirectoryList.push_back (directory);
349   directory = new G4UIdirectory ("/vis/filtering/digi/");
350   directory -> SetGuidance ("Digi filtering commands.");
351   fDirectoryList.push_back (directory);
352   directory = new G4UIdirectory ("/vis/filtering/digi/create/");
353   directory -> SetGuidance ("Create digi filters and messengers.");
354   fDirectoryList.push_back (directory);
355 
356   RegisterMessengers ();
357 
358   if (fVerbosity >= startup) {
359     G4cout << "Registering model factories..." << G4endl;
360   }
361 
362   RegisterModelFactories();
363 
364   if (fVerbosity >= startup) {
365     G4cout <<
366       "\nYou have successfully registered the following model factories."
367      << G4endl;
368     PrintAvailableModels (fVerbosity);
369     G4cout << G4endl;
370   }
371 
372   if (fVerbosity >= startup) {
373     PrintAvailableUserVisActions (fVerbosity);
374     G4cout << G4endl;
375   }
376 
377   InitialiseG4ColourMap();
378 
379   if (fVerbosity >= startup) {
380     G4cout <<
381     "Some /vis commands (optionally) take a string to specify colour."
382     "\n\"/vis/list\" to see available colours."
383     << G4endl;
384   }
385 
386   fInitialised = true;
387 }
388 
389 void G4VisManager::InitialiseG4ColourMap() const
390 {
391   G4Colour::InitialiseColourMap();  // Initialises (if not already initialised)
392 
393   // our forever 65 named colors taken long time ago from X11.
394   // Extracted from g4tools/include/tools/colors
395   // Copyright (C) 2010, Guy Barrand. All rights reserved.
396   // See the file tools.license for terms.
397 
398 #define TOOLS_COLORS_STAT(name,r,g,b) \
399 G4Colour::AddToMap(#name, G4Colour(r,g,b));
400 
401   //0-9
402   TOOLS_COLORS_STAT(aquamarine,0.496101F,0.996109F,0.828138F)
403   TOOLS_COLORS_STAT(mediumaquamarine,0.398444F,0.800793F,0.664073F)
404   //  TOOLS_COLORS_STAT(black,0,0,0) (already defined)
405   //  TOOLS_COLORS_STAT(blue,0,0,1) (already defined)
406   TOOLS_COLORS_STAT(cadetblue,0.371099F,0.617197F,0.62501F)
407   TOOLS_COLORS_STAT(cornflowerblue,0.390631F,0.58204F,0.925795F)
408   TOOLS_COLORS_STAT(darkslateblue,0.281254F,0.238285F,0.542977F)
409   TOOLS_COLORS_STAT(lightblue,0.675792F,0.843763F,0.898451F)
410   TOOLS_COLORS_STAT(lightsteelblue,0.68751F,0.765637F,0.867201F)
411   TOOLS_COLORS_STAT(mediumblue,0,0,0.800793F)
412 
413   //10-19
414   TOOLS_COLORS_STAT(mediumslateblue,0.480476F,0.406256F,0.929702F)
415   TOOLS_COLORS_STAT(midnightblue,0.0976577F,0.0976577F,0.437507F)
416   TOOLS_COLORS_STAT(navyblue,0,0,0.500008F)
417   TOOLS_COLORS_STAT(navy,0,0,0.500008F)
418   TOOLS_COLORS_STAT(skyblue,0.527352F,0.8047F,0.917983F)
419   TOOLS_COLORS_STAT(slateblue,0.414069F,0.351568F,0.800793F)
420   TOOLS_COLORS_STAT(steelblue,0.273442F,0.50782F,0.703136F)
421   TOOLS_COLORS_STAT(coral,0.996109F,0.496101F,0.312505F)
422   //  TOOLS_COLORS_STAT(cyan,0,1,1) (already defined)
423   TOOLS_COLORS_STAT(firebrick,0.695323F,0.132815F,0.132815F)
424 
425   //20-29
426   //  TOOLS_COLORS_STAT(brown,0.644541F,0.164065F,0.164065F) (already defined)
427   TOOLS_COLORS_STAT(gold,0.996109F,0.839857F,0)
428   TOOLS_COLORS_STAT(goldenrod,0.851575F,0.644541F,0.125002F)
429   //  TOOLS_COLORS_STAT(green,0,1,0) (already defined)
430   TOOLS_COLORS_STAT(darkgreen,0,0.390631F,0)
431   TOOLS_COLORS_STAT(darkolivegreen,0.332036F,0.417975F,0.183597F)
432   TOOLS_COLORS_STAT(forestgreen,0.132815F,0.542977F,0.132815F)
433   TOOLS_COLORS_STAT(limegreen,0.195315F,0.800793F,0.195315F)
434   TOOLS_COLORS_STAT(mediumseagreen,0.234379F,0.699229F,0.441413F)
435   TOOLS_COLORS_STAT(mediumspringgreen,0,0.976577F,0.601572F)
436 
437   //30-39
438   TOOLS_COLORS_STAT(palegreen,0.593759F,0.980484F,0.593759F)
439   TOOLS_COLORS_STAT(seagreen,0.17969F,0.542977F,0.339849F)
440   TOOLS_COLORS_STAT(springgreen,0,0.996109F,0.496101F)
441   TOOLS_COLORS_STAT(yellowgreen,0.601572F,0.800793F,0.195315F)
442   TOOLS_COLORS_STAT(darkslategrey,0.183597F,0.308598F,0.308598F)
443   TOOLS_COLORS_STAT(dimgrey,0.410163F,0.410163F,0.410163F)
444   TOOLS_COLORS_STAT(lightgrey,0.824231F,0.824231F,0.824231F)
445   //  TOOLS_COLORS_STAT(grey,0.750011F,0.750011F,0.750011F) (already defined)
446   TOOLS_COLORS_STAT(khaki,0.937514F,0.898451F,0.546883F)
447   //  TOOLS_COLORS_STAT(magenta,1,0,1) (already defined)
448 
449   //40-49
450   TOOLS_COLORS_STAT(maroon,0.68751F,0.187503F,0.375006F)
451   TOOLS_COLORS_STAT(orange,0.996109F,0.644541F,0)
452   TOOLS_COLORS_STAT(orchid,0.851575F,0.437507F,0.83595F)
453   TOOLS_COLORS_STAT(darkorchid,0.597665F,0.195315F,0.796887F)
454   TOOLS_COLORS_STAT(mediumorchid,0.726574F,0.332036F,0.824231F)
455   TOOLS_COLORS_STAT(pink,0.996109F,0.750011F,0.792981F)
456   TOOLS_COLORS_STAT(plum,0.863294F,0.62501F,0.863294F)
457   //  TOOLS_COLORS_STAT(red,1,0,0) (already defined)
458   TOOLS_COLORS_STAT(indianred,0.800793F,0.35938F,0.35938F)
459   TOOLS_COLORS_STAT(mediumvioletred,0.777356F,0.0820325F,0.519539F)
460 
461   //50-59
462   TOOLS_COLORS_STAT(orangered,0.996109F,0.269535F,0)
463   TOOLS_COLORS_STAT(violetred,0.812512F,0.125002F,0.562509F)
464   TOOLS_COLORS_STAT(salmon,0.976577F,0.500008F,0.445319F)
465   TOOLS_COLORS_STAT(sienna,0.62501F,0.320317F,0.175784F)
466   TOOLS_COLORS_STAT(tan,0.820325F,0.703136F,0.546883F)
467   TOOLS_COLORS_STAT(thistle,0.843763F,0.746105F,0.843763F)
468   TOOLS_COLORS_STAT(turquoise,0.250004F,0.875013F,0.812512F)
469   TOOLS_COLORS_STAT(darkturquoise,0,0.8047F,0.816419F)
470   TOOLS_COLORS_STAT(mediumturquoise,0.281254F,0.816419F,0.796887F)
471   TOOLS_COLORS_STAT(violet,0.929702F,0.50782F,0.929702F)
472 
473   //60-64
474   TOOLS_COLORS_STAT(blueviolet,0.539071F,0.167971F,0.882826F)
475   TOOLS_COLORS_STAT(wheat,0.957046F,0.867201F,0.699229F)
476   //  TOOLS_COLORS_STAT(white,1,1,1) (already defined)
477   //  TOOLS_COLORS_STAT(yellow,1,1,0) (already defined)
478   TOOLS_COLORS_STAT(greenyellow,0.675792F,0.996109F,0.18359F)
479 
480 #undef TOOLS_COLORS_STAT
481 }
482 
483 void G4VisManager::RegisterMessengers () {
484   
485   // Instantiate individual messengers/commands (often - but not
486   // always - one command per messenger).
487   
488   G4UIcommand* directory;
489   
490   directory = new G4UIdirectory ("/vis/geometry/");
491   directory -> SetGuidance("Operations on vis attributes of Geant4 geometry.");
492   fDirectoryList.push_back (directory);
493   RegisterMessenger(new G4VisCommandGeometryList);
494   RegisterMessenger(new G4VisCommandGeometryRestore);
495   
496   directory = new G4UIdirectory ("/vis/geometry/set/");
497   directory -> SetGuidance("Set vis attributes of Geant4 geometry.");
498   fDirectoryList.push_back (directory);
499   RegisterMessenger(new G4VisCommandGeometrySetColour);
500   RegisterMessenger(new G4VisCommandGeometrySetDaughtersInvisible);
501   RegisterMessenger(new G4VisCommandGeometrySetLineStyle);
502   RegisterMessenger(new G4VisCommandGeometrySetLineWidth);
503   RegisterMessenger(new G4VisCommandGeometrySetForceAuxEdgeVisible);
504   RegisterMessenger(new G4VisCommandGeometrySetForceCloud);
505   RegisterMessenger(new G4VisCommandGeometrySetForceLineSegmentsPerCircle);
506   RegisterMessenger(new G4VisCommandGeometrySetForceSolid);
507   RegisterMessenger(new G4VisCommandGeometrySetForceWireframe);
508   RegisterMessenger(new G4VisCommandGeometrySetVisibility);
509   
510   directory = new G4UIdirectory ("/vis/multithreading/");
511   directory -> SetGuidance("Commands unique to multithreading mode.");
512   fDirectoryList.push_back (directory);
513   RegisterMessenger(new G4VisCommandMultithreadingActionOnEventQueueFull);
514   RegisterMessenger(new G4VisCommandMultithreadingMaxEventQueueSize);
515 
516   directory = new G4UIdirectory ("/vis/set/");
517   directory -> SetGuidance
518     ("Set quantities for use in future commands where appropriate.");
519   fDirectoryList.push_back (directory);
520   RegisterMessenger(new G4VisCommandSetArrow3DLineSegmentsPerCircle);
521   RegisterMessenger(new G4VisCommandSetColour);
522   RegisterMessenger(new G4VisCommandSetExtentForField);
523   RegisterMessenger(new G4VisCommandSetLineWidth);
524   RegisterMessenger(new G4VisCommandSetTextColour);
525   RegisterMessenger(new G4VisCommandSetTextLayout);
526   RegisterMessenger(new G4VisCommandSetTextSize);
527   RegisterMessenger(new G4VisCommandSetTouchable);
528   RegisterMessenger(new G4VisCommandSetVolumeForField);
529 
530   directory = new G4UIdirectory ("/vis/scene/");
531   directory -> SetGuidance ("Operations on Geant4 scenes.");
532   fDirectoryList.push_back (directory);
533   RegisterMessenger(new G4VisCommandSceneActivateModel);
534   RegisterMessenger(new G4VisCommandSceneCreate);
535   RegisterMessenger(new G4VisCommandSceneEndOfEventAction);
536   RegisterMessenger(new G4VisCommandSceneEndOfRunAction);
537   RegisterMessenger(new G4VisCommandSceneList);
538   RegisterMessenger(new G4VisCommandSceneNotifyHandlers);
539   RegisterMessenger(new G4VisCommandSceneRemoveModel);
540   RegisterMessenger(new G4VisCommandSceneSelect);
541   RegisterMessenger(new G4VisCommandSceneShowExtents);
542 
543   directory = new G4UIdirectory ("/vis/scene/add/");
544   directory -> SetGuidance ("Add model to current scene.");
545   fDirectoryList.push_back (directory);
546   RegisterMessenger(new G4VisCommandSceneAddArrow);
547   RegisterMessenger(new G4VisCommandSceneAddArrow2D);
548   RegisterMessenger(new G4VisCommandSceneAddAxes);
549   RegisterMessenger(new G4VisCommandSceneAddDate);
550   RegisterMessenger(new G4VisCommandSceneAddDigis);
551   RegisterMessenger(new G4VisCommandSceneAddEventID);
552   RegisterMessenger(new G4VisCommandSceneAddExtent);
553   RegisterMessenger(new G4VisCommandSceneAddElectricField);
554   RegisterMessenger(new G4VisCommandSceneAddFrame);
555   RegisterMessenger(new G4VisCommandSceneAddGPS);
556   RegisterMessenger(new G4VisCommandSceneAddHits);
557   RegisterMessenger(new G4VisCommandSceneAddLine);
558   RegisterMessenger(new G4VisCommandSceneAddLine2D);
559   RegisterMessenger(new G4VisCommandSceneAddLocalAxes);
560   RegisterMessenger(new G4VisCommandSceneAddLogicalVolume);
561   RegisterMessenger(new G4VisCommandSceneAddLogo);
562   RegisterMessenger(new G4VisCommandSceneAddLogo2D);
563   RegisterMessenger(new G4VisCommandSceneAddMagneticField);
564   RegisterMessenger(new G4VisCommandSceneAddPlotter);
565   RegisterMessenger(new G4VisCommandSceneAddPSHits);
566   RegisterMessenger(new G4VisCommandSceneAddScale);
567   RegisterMessenger(new G4VisCommandSceneAddText);
568   RegisterMessenger(new G4VisCommandSceneAddText2D);
569   RegisterMessenger(new G4VisCommandSceneAddTrajectories);
570   RegisterMessenger(new G4VisCommandSceneAddUserAction);
571   RegisterMessenger(new G4VisCommandSceneAddVolume);
572   
573   RegisterMessenger(new G4VisCommandPlotterCreate);
574   RegisterMessenger(new G4VisCommandPlotterSetLayout);
575   RegisterMessenger(new G4VisCommandPlotterAddStyle);
576   RegisterMessenger(new G4VisCommandPlotterAddRegionStyle);
577   RegisterMessenger(new G4VisCommandPlotterAddRegionParameter);
578   RegisterMessenger(new G4VisCommandPlotterClear);
579   RegisterMessenger(new G4VisCommandPlotterClearRegion);
580   RegisterMessenger(new G4VisCommandPlotterList);
581   RegisterMessenger(new G4VisCommandPlotterAddRegionH1);
582   RegisterMessenger(new G4VisCommandPlotterAddRegionH2);
583   
584   directory = new G4UIdirectory ("/vis/sceneHandler/");
585   directory -> SetGuidance ("Operations on Geant4 scene handlers.");
586   fDirectoryList.push_back (directory);
587   RegisterMessenger(new G4VisCommandSceneHandlerAttach);
588   RegisterMessenger(new G4VisCommandSceneHandlerCreate);
589   RegisterMessenger(new G4VisCommandSceneHandlerList);
590   RegisterMessenger(new G4VisCommandSceneHandlerSelect);
591   
592   directory = new G4UIdirectory ("/vis/touchable/");
593   directory -> SetGuidance ("Operations on touchables.");
594   fDirectoryList.push_back (directory);
595   RegisterMessenger(new G4VisCommandsTouchable);
596   
597   directory = new G4UIdirectory ("/vis/touchable/set/");
598   directory -> SetGuidance ("Set vis attributes of current touchable.");
599   fDirectoryList.push_back (directory);
600   RegisterMessenger(new G4VisCommandsTouchableSet);
601 
602   directory = new G4UIdirectory ("/vis/viewer/");
603   directory -> SetGuidance ("Operations on Geant4 viewers.");
604   fDirectoryList.push_back (directory);
605   RegisterMessenger(new G4VisCommandViewerAddCutawayPlane);
606   RegisterMessenger(new G4VisCommandViewerCentreOn);
607   RegisterMessenger(new G4VisCommandViewerChangeCutawayPlane);
608   RegisterMessenger(new G4VisCommandViewerClear);
609   RegisterMessenger(new G4VisCommandViewerClearCutawayPlanes);
610   RegisterMessenger(new G4VisCommandViewerClearTransients);
611   RegisterMessenger(new G4VisCommandViewerClearVisAttributesModifiers);
612   RegisterMessenger(new G4VisCommandViewerClone);
613   RegisterMessenger(new G4VisCommandViewerColourByDensity);
614   RegisterMessenger(new G4VisCommandViewerCopyViewFrom);
615   RegisterMessenger(new G4VisCommandViewerCreate);
616   RegisterMessenger(new G4VisCommandViewerDolly);
617   RegisterMessenger(new G4VisCommandViewerFlush);
618   RegisterMessenger(new G4VisCommandViewerInterpolate);
619   RegisterMessenger(new G4VisCommandViewerList);
620   RegisterMessenger(new G4VisCommandViewerPan);
621   RegisterMessenger(new G4VisCommandViewerRebuild);
622   RegisterMessenger(new G4VisCommandViewerRefresh);
623   RegisterMessenger(new G4VisCommandViewerReset);
624   RegisterMessenger(new G4VisCommandViewerResetCameraParameters);
625   RegisterMessenger(new G4VisCommandViewerSave);
626   RegisterMessenger(new G4VisCommandViewerScale);
627   RegisterMessenger(new G4VisCommandViewerSelect);
628   RegisterMessenger(new G4VisCommandViewerUpdate);
629   RegisterMessenger(new G4VisCommandViewerZoom);
630   
631   directory = new G4UIdirectory ("/vis/viewer/default/");
632   directory -> SetGuidance("Set default values for future viewers.");
633   fDirectoryList.push_back (directory);
634   RegisterMessenger(new G4VisCommandViewerDefaultHiddenEdge);
635   RegisterMessenger(new G4VisCommandViewerDefaultStyle);
636   
637   directory = new G4UIdirectory ("/vis/viewer/set/");
638   directory -> SetGuidance ("Set view parameters of current viewer.");
639   fDirectoryList.push_back (directory);
640   RegisterMessenger(new G4VisCommandsViewerSet);
641   
642   // *Basic* top level commands were instantiated in the constructor
643   // so that they can be used immediately after instantiation of the
644   // vis manager.  Other top level commands, including "compound commands"
645   // (i.e., commands that invoke other commands) are instantiated here.
646 
647   RegisterMessenger(new G4VisCommandAbortReviewKeptEvents);
648   RegisterMessenger(new G4VisCommandAbortReviewPlots);
649   RegisterMessenger(new G4VisCommandDrawOnlyToBeKeptEvents);
650   RegisterMessenger(new G4VisCommandDrawTree);
651   RegisterMessenger(new G4VisCommandDrawView);
652   RegisterMessenger(new G4VisCommandDrawLogicalVolume);
653   RegisterMessenger(new G4VisCommandDrawVolume);
654   RegisterMessenger(new G4VisCommandEnable);
655   RegisterMessenger(new G4VisCommandList);
656   RegisterMessenger(new G4VisCommandOpen);
657   RegisterMessenger(new G4VisCommandPlot);
658   RegisterMessenger(new G4VisCommandReviewKeptEvents);
659   RegisterMessenger(new G4VisCommandReviewPlots);
660   RegisterMessenger(new G4VisCommandSpecify);
661 
662   // List manager commands
663   RegisterMessenger(new G4VisCommandListManagerList< G4VisModelManager<G4VTrajectoryModel> >
664         (fpTrajDrawModelMgr, fpTrajDrawModelMgr->Placement()));
665   RegisterMessenger(new G4VisCommandListManagerSelect< G4VisModelManager<G4VTrajectoryModel> >
666         (fpTrajDrawModelMgr, fpTrajDrawModelMgr->Placement()));
667   
668   // Trajectory filter manager commands
669   RegisterMessenger(new G4VisCommandListManagerList< G4VisFilterManager<G4VTrajectory> >
670                     (fpTrajFilterMgr, fpTrajFilterMgr->Placement()));
671   RegisterMessenger(new G4VisCommandManagerMode< G4VisFilterManager<G4VTrajectory> >
672                     (fpTrajFilterMgr, fpTrajFilterMgr->Placement()));
673   
674   // Hit filter manager commands
675   RegisterMessenger(new G4VisCommandListManagerList< G4VisFilterManager<G4VHit> >
676                     (fpHitFilterMgr, fpHitFilterMgr->Placement()));
677   RegisterMessenger(new G4VisCommandManagerMode< G4VisFilterManager<G4VHit> >
678                     (fpHitFilterMgr, fpHitFilterMgr->Placement()));
679   
680   // Digi filter manager commands
681   RegisterMessenger(new G4VisCommandListManagerList< G4VisFilterManager<G4VDigi> >
682                     (fpDigiFilterMgr, fpDigiFilterMgr->Placement()));
683   RegisterMessenger(new G4VisCommandManagerMode< G4VisFilterManager<G4VDigi> >
684                     (fpDigiFilterMgr, fpDigiFilterMgr->Placement()));
685 }
686 
687 #include <tools/histo/h1d>
688 #include <tools/histo/h2d>
689 
690 namespace {
691   struct PlotResults {
692     std::size_t fNumberOfPlots = 0;
693     std::size_t fTotalEntries = 0;
694   };
695   template <typename HT>  // tools::histo::h1d, etc
696   PlotResults ResultsOfHnPlots(const G4String& plotType) {  // h1, etc.
697     PlotResults plotResults;
698     auto ui = G4UImanager::GetUIpointer();
699     auto keepControlVerbose = ui->GetVerboseLevel();
700     ui->SetVerboseLevel(0);
701     auto status = ui->ApplyCommand("/analysis/" + plotType + "/getVector");
702     ui->SetVerboseLevel(keepControlVerbose);
703     if(status==G4UIcommandStatus::fCommandSucceeded) {
704       G4String hexString = ui->GetCurrentValues(G4String("/analysis/" + plotType + "/getVector"));
705       if(hexString.size()) {
706         void* ptr; std::istringstream is(hexString); is >> ptr;
707         auto vectorOfPlots = (const std::vector<HT*>*)ptr;
708         for (std::size_t i = 0; i < vectorOfPlots->size(); ++i) {
709           auto plot = (*vectorOfPlots)[i];
710           if (plot == nullptr) continue;  // Ignore deleted plots
711           ++plotResults.fNumberOfPlots;
712           plotResults.fTotalEntries += plot->entries();
713         }
714       }
715     }
716     return plotResults;
717   }
718   void PrintListOfPlots() {
719     std::size_t numberOfPlots = 0;
720     std::size_t numberOfEntries = 0;
721     PlotResults h1results = ResultsOfHnPlots<tools::histo::h1d>("h1");
722     numberOfPlots += h1results.fNumberOfPlots;
723     numberOfEntries += h1results.fTotalEntries;
724     PlotResults h2results = ResultsOfHnPlots<tools::histo::h2d>("h2");
725     numberOfPlots += h2results.fNumberOfPlots;
726     numberOfEntries += h2results.fTotalEntries;
727     if (numberOfPlots > 0) {
728       G4warn << "There are histograms that can be viewed with visualization:";
729       if (h1results.fNumberOfPlots > 0) {
730         G4warn << "\n  " << h1results.fNumberOfPlots << " h1 histograms(s)";
731       }
732       if (h2results.fNumberOfPlots > 0) {
733         G4warn << "\n  " << h2results.fNumberOfPlots << " h2 histograms(s)";
734       }
735       G4warn
736       << "\n  List them with \"/analysis/list\"."
737       << "\n  View them immediately with \"/vis/plot\" or \"/vis/reviewPlots\"."
738       << G4endl;
739       if (numberOfEntries == 0) {
740         G4warn <<
741         "  But...there are no entries. To make your histograms available for"
742         "\n  plotting in this UI session, use CloseFile(false) in your"
743         "\n  EndOfRunAction and Reset() in your BeginOfRunAction."
744         << G4endl;
745       }
746     }
747   }
748 }
749 
750 void G4VisManager::Enable() {
751   if (IsValidView ()) {
752     SetConcreteInstance(this);
753     if (fVerbosity >= confirmations) {
754       G4cout << "G4VisManager::Enable: visualization enabled." << G4endl;
755     }
756     if (fVerbosity >= warnings) {
757       std::size_t nKeptEvents = 0;
758       const G4Run* run = G4RunManager::GetRunManager()->GetCurrentRun();
759       if (run) nKeptEvents = run->GetNumberOfKeptEvents();
760       G4String isare("are"),plural("s");
761       if (nKeptEvents == 1) {isare = "is"; plural = "";}
762       G4cout <<
763       "There " << isare << ' ' << nKeptEvents << " kept event" << plural << '.'
764       << G4endl;
765       if (nKeptEvents > 0) {
766         G4cout <<
767   "  \"/vis/reviewKeptEvents\" to review one by one."
768   "\n  To see accumulated, \"/vis/enable\", then \"/vis/viewer/flush\" or \"/vis/viewer/rebuild\"."
769         << G4endl;
770       }
771       PrintListOfPlots();
772     }
773   }
774   else {
775     if (fVerbosity >= warnings) {
776       G4warn <<
777   "G4VisManager::Enable: WARNING: visualization remains disabled for"
778   "\n  above reasons.  Rectifying with valid vis commands will"
779   "\n  automatically enable."
780        << G4endl;
781     }
782   }
783 }
784 
785 void G4VisManager::Disable() {
786   SetConcreteInstance(0);
787   if (fVerbosity >= confirmations) {
788     G4cout <<
789     "G4VisManager::Disable: visualization disabled."
790     "\n  The pointer returned by GetConcreteInstance will be zero."
791     "\n  Note that it will become enabled after some valid vis commands."
792      << G4endl;
793   }
794   if (fVerbosity >= warnings) {
795     G4int currentTrajectoryType =
796     G4RunManagerKernel::GetRunManagerKernel()->GetTrackingManager()->GetStoreTrajectory();
797     if (currentTrajectoryType > 0) {
798     G4warn <<
799       "You may wish to disable trajectory production too:"
800       "\n  \"/tracking/storeTrajectory 0\""
801       "\nbut don't forget to re-enable with"
802       "\n  \"/vis/enable\""
803       "\n  \"/tracking/storeTrajectory " << currentTrajectoryType
804       << "\"\n  and maybe \"/vis/viewer/rebuild\""
805       << G4endl;
806     }
807   }
808 }
809 
810 const G4GraphicsSystemList& G4VisManager::GetAvailableGraphicsSystems () {
811   std::size_t nSystems = fAvailableGraphicsSystems.size ();
812   if (nSystems == 0) {
813     if (fVerbosity >= warnings) {
814       G4warn << "G4VisManager::GetAvailableGraphicsSystems: WARNING: no"
815   "\n graphics system available!"
816   "\n  1) Did you have environment variables G4VIS_BUILD_xxxx_DRIVER set"
817   "\n     when you compiled/built the visualization code?"
818   "\n  2) Did you instantiate your own Visualization Manager and forget"
819   "\n     to implement RegisterGraphicsSystems correctly?"
820   "\n  3) You can register your own graphics system, e.g.,"
821   "\n     visManager->RegisterGraphicsSystem(new MyGraphicsSystem);)"
822   "\n     after instantiating your vis manager and before"
823   "\n     visManager->Initialize()."
824        << G4endl;
825     }
826   }
827   return fAvailableGraphicsSystems;
828 }
829 
830 G4bool G4VisManager::RegisterGraphicsSystem (G4VGraphicsSystem* pSystem) {
831   G4bool happy = true;
832   if (pSystem) {
833     fAvailableGraphicsSystems.push_back (pSystem);
834     if (fVerbosity >= confirmations) {
835       G4cout << "G4VisManager::RegisterGraphicsSystem: "
836        << pSystem -> GetName ();
837       if (pSystem -> GetNickname () != "") {
838   G4cout << " (" << pSystem -> GetNickname () << ")";
839       }
840       G4cout << " registered." << G4endl;
841     }
842   }
843   else {
844     if (fVerbosity >= errors) {
845       G4warn << "G4VisManager::RegisterGraphicsSystem: null pointer!"
846        << G4endl;
847     }
848     happy=false;
849   }
850   return happy;
851 }
852 
853 const G4VTrajectoryModel*
854 G4VisManager::CurrentTrajDrawModel() const
855 {
856   assert (0 != fpTrajDrawModelMgr);
857 
858   const G4VTrajectoryModel* model = fpTrajDrawModelMgr->Current();
859 
860   if (0 == model) {
861     // No model was registered with the trajectory model manager.
862     // Use G4TrajectoryDrawByCharge as a fallback.
863     fpTrajDrawModelMgr->Register(new G4TrajectoryDrawByCharge("DefaultModel"));
864     if (fVerbosity >= warnings) {
865       G4warn<<"G4VisManager: Using G4TrajectoryDrawByCharge as fallback trajectory model."<<G4endl;
866       G4warn<<"See commands in /vis/modeling/trajectories/ for other options."<<G4endl;
867     }
868   }
869 
870   model = fpTrajDrawModelMgr->Current();
871   assert (0 != model); // Should definitely exist now
872 
873   return model;
874 }
875 
876 void G4VisManager::RegisterModel(G4VTrajectoryModel* model)
877 {
878   fpTrajDrawModelMgr->Register(model);
879 }
880 
881 void
882 G4VisManager::RegisterModelFactory(G4TrajDrawModelFactory* factory) 
883 {
884   fpTrajDrawModelMgr->Register(factory);
885 }
886 
887 void G4VisManager::RegisterModel(G4VFilter<G4VTrajectory>* model)
888 {
889   fpTrajFilterMgr->Register(model);
890 }
891 
892 void
893 G4VisManager::RegisterModelFactory(G4TrajFilterFactory* factory)
894 {
895   fpTrajFilterMgr->Register(factory);
896 }
897 
898 void G4VisManager::RegisterModel(G4VFilter<G4VHit>* model)
899 {
900   fpHitFilterMgr->Register(model);
901 }
902 
903 void
904 G4VisManager::RegisterModelFactory(G4HitFilterFactory* factory)
905 {
906   fpHitFilterMgr->Register(factory);
907 }
908 
909 void G4VisManager::RegisterModel(G4VFilter<G4VDigi>* model)
910 {
911   fpDigiFilterMgr->Register(model);
912 }
913 
914 void
915 G4VisManager::RegisterModelFactory(G4DigiFilterFactory* factory)
916 {
917   fpDigiFilterMgr->Register(factory);
918 }
919 
920 void G4VisManager::SelectTrajectoryModel(const G4String& model) 
921 {
922    fpTrajDrawModelMgr->SetCurrent(model);
923 }
924 
925 void G4VisManager::BeginDraw (const G4Transform3D& objectTransform)
926 {
927   if (G4Threading::IsWorkerThread()) return;
928   
929   fDrawGroupNestingDepth++;
930   if (fDrawGroupNestingDepth > 1) {
931     G4Exception
932       ("G4VisManager::BeginDraw",
933        "visman0008", JustWarning,
934        "Nesting detected. It is illegal to nest Begin/EndDraw."
935        "\n Ignored");
936     return;
937   }
938   if (IsValidView ()) {
939     ClearTransientStoreIfMarked();
940     fpSceneHandler -> BeginPrimitives (objectTransform);
941     fIsDrawGroup = true;
942   }
943 }
944 
945 void G4VisManager::EndDraw ()
946 {
947   if (G4Threading::IsWorkerThread()) return;
948   
949   fDrawGroupNestingDepth--;
950   if (fDrawGroupNestingDepth != 0) {
951     if (fDrawGroupNestingDepth < 0) fDrawGroupNestingDepth = 0;
952     return;
953   }
954   if (IsValidView ()) {
955     fpSceneHandler -> EndPrimitives ();
956   }
957   fIsDrawGroup = false;
958 }
959 
960 void G4VisManager::BeginDraw2D (const G4Transform3D& objectTransform)
961 {
962   if (G4Threading::IsWorkerThread()) return;
963   
964   fDrawGroupNestingDepth++;
965   if (fDrawGroupNestingDepth > 1) {
966     G4Exception
967       ("G4VisManager::BeginDraw2D",
968        "visman0009", JustWarning,
969        "Nesting detected. It is illegal to nest Begin/EndDraw2D."
970        "\n Ignored");
971     return;
972   }
973   if (IsValidView ()) {
974     ClearTransientStoreIfMarked();
975     fpSceneHandler -> BeginPrimitives2D (objectTransform);
976     fIsDrawGroup = true;
977   }
978 }
979 
980 void G4VisManager::EndDraw2D ()
981 {
982   if (G4Threading::IsWorkerThread()) return;
983   
984   fDrawGroupNestingDepth--;
985   if (fDrawGroupNestingDepth != 0) {
986     if (fDrawGroupNestingDepth < 0) fDrawGroupNestingDepth = 0;
987     return;
988   }
989   if (IsValidView ()) {
990     fpSceneHandler -> EndPrimitives2D ();
991   }
992   fIsDrawGroup = false;
993 }
994 
995 template <class T> void G4VisManager::DrawT
996 (const T& graphics_primitive, const G4Transform3D& objectTransform) {
997   if (G4Threading::IsWorkerThread()) return;
998   
999   if (fIsDrawGroup) {
1000     if (objectTransform != fpSceneHandler->GetObjectTransformation()) {
1001       G4Exception
1002   ("G4VSceneHandler::DrawT",
1003    "visman0010", FatalException,
1004    "Different transform detected in Begin/EndDraw group.");
1005     }
1006     fpSceneHandler -> AddPrimitive (graphics_primitive);
1007   } else {
1008     if (IsValidView ()) {
1009       ClearTransientStoreIfMarked();
1010       fpSceneHandler -> BeginPrimitives (objectTransform);
1011       fpSceneHandler -> AddPrimitive (graphics_primitive);
1012       fpSceneHandler -> EndPrimitives ();
1013     }
1014   }
1015 }
1016 
1017 template <class T> void G4VisManager::DrawT2D
1018 (const T& graphics_primitive, const G4Transform3D& objectTransform) {
1019   if (G4Threading::IsWorkerThread()) return;
1020   
1021   if (fIsDrawGroup) {
1022     if (objectTransform != fpSceneHandler->GetObjectTransformation()) {
1023       G4Exception
1024   ("G4VSceneHandler::DrawT",
1025    "visman0011", FatalException,
1026    "Different transform detected in Begin/EndDraw2D group.");
1027     }
1028     fpSceneHandler -> AddPrimitive (graphics_primitive);
1029   } else {
1030     if (IsValidView ()) {
1031       ClearTransientStoreIfMarked();
1032       fpSceneHandler -> BeginPrimitives2D (objectTransform);
1033       fpSceneHandler -> AddPrimitive (graphics_primitive);
1034       fpSceneHandler -> EndPrimitives2D ();
1035     }
1036   }
1037 }
1038 
1039 void G4VisManager::Draw (const G4Circle& circle,
1040        const G4Transform3D& objectTransform)
1041 {
1042   DrawT (circle, objectTransform);
1043 }
1044 
1045 void G4VisManager::Draw (const G4Polyhedron& polyhedron,
1046        const G4Transform3D& objectTransform)
1047 {
1048   DrawT (polyhedron, objectTransform);
1049 }
1050 
1051 void G4VisManager::Draw (const G4Polyline& line,
1052        const G4Transform3D& objectTransform)
1053 {
1054   DrawT (line, objectTransform);
1055 }
1056 
1057 void G4VisManager::Draw (const G4Polymarker& polymarker,
1058        const G4Transform3D& objectTransform)
1059 {
1060   DrawT (polymarker, objectTransform);
1061 }
1062 
1063 void G4VisManager::Draw (const G4Square& square,
1064        const G4Transform3D& objectTransform)
1065 {
1066   DrawT (square, objectTransform);
1067 }
1068 
1069 void G4VisManager::Draw (const G4Text& text,
1070        const G4Transform3D& objectTransform)
1071 {
1072   DrawT (text, objectTransform);
1073 }
1074 
1075 void G4VisManager::Draw2D (const G4Circle& circle,
1076          const G4Transform3D& objectTransform)
1077 {
1078   DrawT2D (circle, objectTransform);
1079 }
1080 
1081 void G4VisManager::Draw2D (const G4Polyhedron& polyhedron,
1082          const G4Transform3D& objectTransform)
1083 {
1084   DrawT2D (polyhedron, objectTransform);
1085 }
1086 
1087 void G4VisManager::Draw2D (const G4Polyline& line,
1088          const G4Transform3D& objectTransform)
1089 {
1090   DrawT2D (line, objectTransform);
1091 }
1092 
1093 void G4VisManager::Draw2D (const G4Polymarker& polymarker,
1094          const G4Transform3D& objectTransform)
1095 {
1096   DrawT2D (polymarker, objectTransform);
1097 }
1098 
1099 void G4VisManager::Draw2D (const G4Square& square,
1100          const G4Transform3D& objectTransform)
1101 {
1102   DrawT2D (square, objectTransform);
1103 }
1104 
1105 void G4VisManager::Draw2D (const G4Text& text,
1106          const G4Transform3D& objectTransform)
1107 {
1108   DrawT2D (text, objectTransform);
1109 }
1110 
1111 void G4VisManager::Draw (const G4VHit& hit) {
1112   if (G4Threading::IsWorkerThread()) return;
1113   
1114   if (fIsDrawGroup) {
1115     fpSceneHandler -> AddCompound (hit);
1116   } else {
1117     if (IsValidView ()) {
1118       ClearTransientStoreIfMarked();
1119       fpSceneHandler -> AddCompound (hit);
1120     }
1121   }
1122 }
1123 
1124 void G4VisManager::Draw (const G4VDigi& digi) {
1125   if (G4Threading::IsWorkerThread()) return;
1126   
1127   if (fIsDrawGroup) {
1128     fpSceneHandler -> AddCompound (digi);
1129   } else {
1130     if (IsValidView ()) {
1131       ClearTransientStoreIfMarked();
1132       fpSceneHandler -> AddCompound (digi);
1133     }
1134   }
1135 }
1136 
1137 void G4VisManager::Draw (const G4VTrajectory& traj) {
1138   if (G4Threading::IsWorkerThread()) return;
1139   
1140   // A trajectory needs a trajectories model to provide G4Atts, etc.
1141   static G4TrajectoriesModel trajectoriesModel;
1142   trajectoriesModel.SetCurrentTrajectory(&traj);
1143   G4RunManager* runManager = G4RunManagerFactory::GetMasterRunManager();
1144   const G4Run* currentRun  = runManager->GetCurrentRun();
1145   if (currentRun) {
1146     trajectoriesModel.SetRunID(currentRun->GetRunID());
1147   }
1148   const G4Event* currentEvent =
1149   G4EventManager::GetEventManager()->GetConstCurrentEvent();
1150   if (currentEvent) {
1151     trajectoriesModel.SetEventID(currentEvent->GetEventID());
1152   }
1153   if (fIsDrawGroup) {
1154     fpSceneHandler -> SetModel (&trajectoriesModel);
1155     fpSceneHandler -> AddCompound (traj);
1156     fpSceneHandler -> SetModel (0);
1157   } else {
1158     if (IsValidView ()) {
1159       ClearTransientStoreIfMarked();
1160       fpSceneHandler -> SetModel (&trajectoriesModel);
1161       fpSceneHandler -> AddCompound (traj);
1162       fpSceneHandler -> SetModel (0);
1163     }
1164   }
1165 }
1166 
1167 void G4VisManager::Draw (const G4LogicalVolume& logicalVol,
1168        const G4VisAttributes& attribs,
1169        const G4Transform3D& objectTransform) {
1170   if (G4Threading::IsWorkerThread()) return;
1171   
1172   // Find corresponding solid.
1173   G4VSolid* pSol = logicalVol.GetSolid ();
1174   Draw (*pSol, attribs, objectTransform);
1175 }
1176 
1177 void G4VisManager::Draw (const G4VSolid& solid,
1178        const G4VisAttributes& attribs,
1179        const G4Transform3D& objectTransform) {
1180   if (G4Threading::IsWorkerThread()) return;
1181   
1182   if (fIsDrawGroup) {
1183     fpSceneHandler -> PreAddSolid (objectTransform, attribs);
1184     solid.DescribeYourselfTo (*fpSceneHandler);
1185     fpSceneHandler -> PostAddSolid ();
1186   } else {
1187     if (IsValidView ()) {
1188       ClearTransientStoreIfMarked();
1189       fpSceneHandler -> PreAddSolid (objectTransform, attribs);
1190       solid.DescribeYourselfTo (*fpSceneHandler);
1191       fpSceneHandler -> PostAddSolid ();
1192     }
1193   }
1194 }
1195 
1196 void G4VisManager::Draw (const G4VPhysicalVolume& physicalVol,
1197        const G4VisAttributes& attribs,
1198        const G4Transform3D& objectTransform) {
1199   if (G4Threading::IsWorkerThread()) return;
1200   
1201   // Note: It is tempting to use a temporary model here, as for
1202   // trajectories, in order to get at the G4Atts of the physical
1203   // volume.  I tried it (JA).  But it's not easy to pass the
1204   // vis attributes.  Also other aspects of the model seem not to
1205   // be properly set up.  So, the idea has been abandoned for the time
1206   // being.  The model pointer will be null.  So when picking there
1207   // will be no G4Atts from this physical volume.
1208   //
1209   // If this is called from DrawHit, for example, the user may G4Atts to the
1210   // hit and these will be available with "/vis/scene/add/hits".
1211   //
1212   // Find corresponding logical volume and solid.
1213   G4LogicalVolume* pLV  = physicalVol.GetLogicalVolume ();
1214   G4VSolid*        pSol = pLV -> GetSolid ();
1215   Draw (*pSol, attribs, objectTransform);
1216 }
1217 
1218 void G4VisManager::DrawGeometry
1219 (G4VPhysicalVolume* v, const G4Transform3D& t)
1220 // Draws a geometry tree starting at the specified physical volume.
1221 {
1222   auto modelingParameters = fpSceneHandler->CreateModelingParameters();
1223   auto depth = G4PhysicalVolumeModel::UNLIMITED;
1224   const G4bool useFullExtent = true;
1225   G4PhysicalVolumeModel aPVModel(v,depth,t,modelingParameters,useFullExtent);
1226   aPVModel.DescribeYourselfTo(*fpSceneHandler);
1227   delete modelingParameters;
1228 }
1229 
1230 void G4VisManager::CreateSceneHandler (const G4String& name) {
1231   if (!fInitialised) Initialise ();
1232   if (fpGraphicsSystem) {
1233     G4VSceneHandler* pSceneHandler =
1234       fpGraphicsSystem -> CreateSceneHandler (name);
1235     if (pSceneHandler) {
1236       fAvailableSceneHandlers.push_back (pSceneHandler);
1237       fpSceneHandler = pSceneHandler;                         // Make current.
1238     }
1239     else {
1240       if (fVerbosity >= errors) {
1241   G4warn << "ERROR in G4VisManager::CreateSceneHandler during "
1242          << fpGraphicsSystem -> GetName ()
1243          << " scene handler creation.\n  No action taken."
1244          << G4endl;
1245       }
1246     }
1247   }
1248   else PrintInvalidPointers ();
1249 }
1250 
1251 void G4VisManager::CreateViewer
1252 (const G4String& name, const G4String& XGeometry)
1253 {
1254 
1255   if (!fInitialised) Initialise ();
1256 
1257   if (!fpSceneHandler) {
1258     PrintInvalidPointers ();
1259     return;
1260   }
1261 
1262   G4VViewer* p = fpGraphicsSystem -> CreateViewer (*fpSceneHandler, name);
1263 
1264   if (!p) {
1265     if (fVerbosity >= errors) {
1266       G4warn << "ERROR in G4VisManager::CreateViewer: null pointer during "
1267        << fpGraphicsSystem -> GetName ()
1268        << " viewer creation.\n  No action taken."
1269        << G4endl;
1270     }
1271     return;
1272   }
1273 
1274   if (p -> GetViewId() < 0) {
1275     if (fVerbosity >= errors) {
1276       G4warn << "ERROR in G4VisManager::CreateViewer during "
1277        << fpGraphicsSystem -> GetName ()
1278        << " viewer instantiation.\n  No action taken."
1279        << G4endl;
1280     }
1281     return;
1282   }
1283 
1284   // Viewer is created, now we can set geometry parameters
1285   // Before 12/2008, it was done in G4VViewer.cc but it did not have to be there!
1286     
1287   G4ViewParameters initialvp = p -> GetViewParameters();
1288   initialvp.SetXGeometryString(XGeometry); //parse string and store parameters
1289   p -> SetViewParameters(initialvp);
1290   p -> Initialise ();  // (Viewer itself may change view parameters further.)
1291   if (p -> GetViewId() < 0) {
1292     if (fVerbosity >= errors) {
1293       G4warn << "ERROR in G4VisManager::CreateViewer during "
1294        << fpGraphicsSystem -> GetName ()
1295        << " viewer initialisation.\n  No action taken."
1296        << G4endl;
1297     }
1298     return;
1299   }
1300 
1301   fpViewer = p;                             // Make current.
1302   fpSceneHandler -> AddViewerToList (fpViewer);
1303   fpSceneHandler -> SetCurrentViewer (fpViewer);
1304   if (fVerbosity >= confirmations) {
1305     G4cout << "G4VisManager::CreateViewer: new viewer created."
1306      << G4endl;
1307   }
1308 
1309   const G4ViewParameters& vp = fpViewer->GetViewParameters();
1310   if (fVerbosity >= parameters) {
1311     G4cout << " view parameters are:\n  " << vp << G4endl;
1312   }
1313 
1314   if (vp.IsCulling () && vp.IsCullingInvisible ()) {
1315     static G4bool warned = false;
1316     if (fVerbosity >= confirmations) {
1317       if (!warned) {
1318   G4cout <<
1319   "NOTE: objects with visibility flag set to \"false\""
1320   " will not be drawn!"
1321   "\n  \"/vis/viewer/set/culling global false\" to Draw such objects."
1322   "\n  Also see other \"/vis/viewer/set\" commands."
1323          << G4endl;
1324   warned = true;
1325       }
1326     }
1327   }
1328   if (vp.IsCullingCovered ()) {
1329     static G4bool warned = false;
1330     if (fVerbosity >= warnings) {
1331       if (!warned) {
1332   G4warn <<
1333   "WARNING: covered objects in solid mode will not be rendered!"
1334   "\n  \"/vis/viewer/set/culling coveredDaughters false\" to reverse this."
1335   "\n  Also see other \"/vis/viewer/set\" commands."
1336          << G4endl;
1337   warned = true;
1338       }
1339     }
1340   }
1341 }
1342 
1343 void G4VisManager::GeometryHasChanged () {
1344   if (fVerbosity >= confirmations) {
1345     G4cout << "G4VisManager::GeometryHasChanged() called." << G4endl;
1346   }
1347 
1348   // Change the world...
1349   G4VPhysicalVolume* pWorld =
1350     G4TransportationManager::GetTransportationManager ()
1351     -> GetNavigatorForTracking () -> GetWorldVolume ();
1352   if (!pWorld) {
1353     if (fVerbosity >= warnings) {
1354       G4warn << "WARNING: There is no world volume!" << G4endl;
1355     }
1356   }
1357 
1358   // Check scenes.
1359   G4SceneList& sceneList = fSceneList;
1360   std::size_t iScene, nScenes = sceneList.size ();
1361   for (iScene = 0; iScene < nScenes; ++iScene) {
1362     G4Scene* pScene = sceneList [iScene];
1363     std::vector<G4Scene::Model>& modelList = pScene -> SetRunDurationModelList ();
1364     if (modelList.size ()) {
1365       G4bool modelInvalid;
1366       do {  // Remove, if required, one at a time.
1367   modelInvalid = false;
1368   std::vector<G4Scene::Model>::iterator iterModel;
1369   for (iterModel = modelList.begin();
1370        iterModel != modelList.end();
1371        ++iterModel) {
1372     modelInvalid = !(iterModel->fpModel->Validate(fVerbosity>=warnings));
1373     if (modelInvalid) {
1374       // Model invalid - remove and break.
1375       if (fVerbosity >= warnings) {
1376         G4warn << "WARNING: Model \""
1377          << iterModel->fpModel->GetGlobalDescription ()
1378          <<
1379     "\" is no longer valid - being removed\n  from scene \""
1380          << pScene -> GetName () << "\""
1381          << G4endl;
1382       }
1383       modelList.erase (iterModel);
1384       break;
1385     }
1386   }
1387       } while (modelInvalid);
1388 
1389       if (modelList.size () == 0) {
1390   if (fVerbosity >= warnings) {
1391     G4warn << "WARNING: No run-duration models left in this scene \""
1392      << pScene -> GetName ()
1393      << "\"."
1394      << G4endl;
1395   }
1396         if (pWorld) {
1397           if (fVerbosity >= warnings) {
1398             G4warn << "  Adding current world to \""
1399             << pScene -> GetName ()
1400             << "\"."
1401             << G4endl;
1402           }
1403           pScene->AddRunDurationModel(new G4PhysicalVolumeModel(pWorld),fVerbosity>=warnings);
1404           // (The above includes a re-calculation of the extent.)
1405           G4UImanager::GetUIpointer () ->
1406           ApplyCommand (G4String("/vis/scene/notifyHandlers " + pScene->GetName()));
1407         }
1408       }
1409       else {
1410   pScene->CalculateExtent();  // Recalculate extent
1411   G4UImanager::GetUIpointer () ->
1412     ApplyCommand (G4String("/vis/scene/notifyHandlers " + pScene->GetName()));
1413       }
1414     }
1415   }
1416 
1417   // Check the manager's current scene...
1418   if (fpScene && fpScene -> GetRunDurationModelList ().size () == 0) {
1419     if (fVerbosity >= warnings) {
1420       G4warn << "WARNING: The current scene \""
1421        << fpScene -> GetName ()
1422        << "\" has no run duration models."
1423              << "\n  Use \"/vis/scene/add/volume\" or create a new scene."
1424        << G4endl;
1425     }
1426     // Clean up
1427     if (fpSceneHandler) {
1428       fpSceneHandler->ClearTransientStore();
1429       fpSceneHandler->ClearStore();
1430       if (fpViewer) {
1431         fpViewer->NeedKernelVisit();
1432         fpViewer->SetView();
1433         fpViewer->ClearView();
1434         fpViewer->FinishView();
1435       }
1436     }
1437   }
1438 }
1439 
1440 void G4VisManager::NotifyHandlers () {
1441 
1442   if (fVerbosity >= confirmations) {
1443     G4cout << "G4VisManager::NotifyHandler() called." << G4endl;
1444   }
1445 
1446   if (IsValidView()) {
1447 
1448     // Check scenes.
1449     G4SceneList& sceneList = fSceneList;
1450     std::size_t iScene, nScenes = sceneList.size ();
1451     for (iScene = 0; iScene < nScenes; ++iScene) {
1452       G4Scene* pScene = sceneList [iScene];
1453       std::vector<G4Scene::Model>& modelList = pScene -> SetRunDurationModelList ();
1454 
1455       if (modelList.size ()) {
1456         pScene->CalculateExtent();
1457         G4UImanager::GetUIpointer () ->
1458         ApplyCommand (G4String("/vis/scene/notifyHandlers " + pScene->GetName()));
1459       }
1460     }
1461 
1462     // Check the manager's current scene...
1463     if (fpScene && fpScene -> GetRunDurationModelList ().size () == 0) {
1464       if (fVerbosity >= warnings) {
1465         G4warn << "WARNING: The current scene \""
1466         << fpScene -> GetName ()
1467         << "\" has no run duration models."
1468         << "\n  Use \"/vis/scene/add/volume\" or create a new scene."
1469         << G4endl;
1470       }
1471       fpSceneHandler->ClearTransientStore();
1472       fpSceneHandler->ClearStore();
1473       fpViewer->NeedKernelVisit();
1474       fpViewer->SetView();
1475       fpViewer->ClearView();
1476       fpViewer->FinishView();
1477     }
1478   }
1479 }
1480 
1481 G4bool G4VisManager::FilterTrajectory(const G4VTrajectory& trajectory)
1482 {
1483   return fpTrajFilterMgr->Accept(trajectory);
1484 }   
1485 
1486 G4bool G4VisManager::FilterHit(const G4VHit& hit)
1487 {
1488   return fpHitFilterMgr->Accept(hit);
1489 }   
1490 
1491 G4bool G4VisManager::FilterDigi(const G4VDigi& digi)
1492 {
1493   return fpDigiFilterMgr->Accept(digi);
1494 }   
1495 
1496 void G4VisManager::DispatchToModel(const G4VTrajectory& trajectory)
1497 {
1498   G4bool visible(true);
1499 
1500   // See if trajectory passes filter
1501   G4bool passed = FilterTrajectory(trajectory);
1502 
1503   if (!passed) {
1504     // Draw invisible trajectory if trajectory failed filter and
1505     // are filtering in soft mode
1506     if (fpTrajFilterMgr->GetMode() == FilterMode::Soft) visible = false;
1507     else {return;}
1508   }
1509 
1510   // Go on to draw trajectory
1511   assert (0 != fpTrajDrawModelMgr);
1512 
1513   const G4VTrajectoryModel* trajectoryModel = CurrentTrajDrawModel();
1514 
1515   assert (0 != trajectoryModel); // Should exist
1516 
1517   if (IsValidView()) {
1518       trajectoryModel->Draw(trajectory, visible);
1519   }
1520 }
1521 
1522 void G4VisManager::RegisterRunDurationUserVisAction
1523 (const G4String& name,
1524  G4VUserVisAction* pVisAction,
1525  const G4VisExtent& extent) {
1526   fRunDurationUserVisActions.push_back(UserVisAction(name,pVisAction));
1527   if (extent.GetExtentRadius() > 0.) {
1528     fUserVisActionExtents[pVisAction] = extent;
1529   } else {
1530     if (fVerbosity >= warnings) {
1531       G4warn <<
1532   "WARNING: No extent set for user vis action \"" << name << "\"."
1533        << G4endl;
1534     }
1535   }
1536   if (fVerbosity >= confirmations) {
1537     G4cout
1538     << "Run duration user vis action \"" << name << "\" registered"
1539     << G4endl;
1540   }
1541 }
1542 
1543 void G4VisManager::RegisterEndOfEventUserVisAction
1544 (const G4String& name,
1545  G4VUserVisAction* pVisAction,
1546  const G4VisExtent& extent) {
1547   fEndOfEventUserVisActions.push_back(UserVisAction(name,pVisAction));
1548   if (extent.GetExtentRadius() > 0.) {
1549     fUserVisActionExtents[pVisAction] = extent;
1550   } else {
1551     if (fVerbosity >= warnings) {
1552       G4warn <<
1553   "WARNING: No extent set for user vis action \"" << name << "\"."
1554        << G4endl;
1555     }
1556   }
1557   if (fVerbosity >= confirmations) {
1558     G4cout
1559     << "End of event user vis action \"" << name << "\" registered"
1560     << G4endl;
1561   }
1562 }
1563 
1564 void G4VisManager::RegisterEndOfRunUserVisAction
1565 (const G4String& name,
1566  G4VUserVisAction* pVisAction,
1567  const G4VisExtent& extent) {
1568   fEndOfRunUserVisActions.push_back(UserVisAction(name,pVisAction));
1569   if (extent.GetExtentRadius() > 0.) {
1570     fUserVisActionExtents[pVisAction] = extent;
1571   } else {
1572     if (fVerbosity >= warnings) {
1573       G4warn <<
1574   "WARNING: No extent set for user vis action \"" << name << "\"."
1575        << G4endl;
1576     }
1577   }
1578   if (fVerbosity >= confirmations) {
1579     G4cout
1580     << "End of run user vis action \"" << name << "\" registered"
1581     << G4endl;
1582   }
1583 }
1584 
1585 void G4VisManager::SetCurrentScene (G4Scene* pScene) {
1586   if (pScene != fpScene) {
1587     // A change of scene.  Therefore reset transients drawn flags.  All
1588     // memory of previous transient proceessing thereby erased...
1589     ResetTransientsDrawnFlags();
1590   }
1591   fpScene = pScene;
1592 }
1593 
1594 void G4VisManager::SetCurrentGraphicsSystem (G4VGraphicsSystem* pSystem) {
1595   fpGraphicsSystem = pSystem;
1596   if (fVerbosity >= confirmations) {
1597     G4cout << "G4VisManager::SetCurrentGraphicsSystem: system now "
1598      << pSystem -> GetName () << G4endl;
1599   }
1600   // If current scene handler is of same graphics system, leave unchanged.
1601   // Else find the most recent scene handler of same graphics system.
1602   // Or clear pointers.
1603   if (!(fpSceneHandler && fpSceneHandler -> GetGraphicsSystem () == pSystem)) {
1604     const G4SceneHandlerList& sceneHandlerList = fAvailableSceneHandlers;
1605     G4int nSH = (G4int)sceneHandlerList.size ();  // No. of scene handlers.
1606     G4int iSH;
1607     for (iSH = nSH - 1; iSH >= 0; iSH--) {
1608       if (sceneHandlerList [iSH] -> GetGraphicsSystem () == pSystem) break;
1609     }
1610     if (iSH >= 0) {
1611       fpSceneHandler = sceneHandlerList [iSH];
1612       if (fVerbosity >= confirmations) {
1613   G4cout << "  Scene Handler now "
1614          << fpSceneHandler -> GetName () << G4endl;
1615       }
1616       if (fpScene != fpSceneHandler -> GetScene ()) {
1617   fpScene = fpSceneHandler -> GetScene ();
1618   if (fVerbosity >= confirmations) {
1619     G4cout << "  Scene now \""
1620      << fpScene -> GetName () << "\"" << G4endl;
1621   }
1622       }
1623       const G4ViewerList& viewerList = fpSceneHandler -> GetViewerList ();
1624       if (viewerList.size ()) {
1625   fpViewer = viewerList [0];
1626   if (fVerbosity >= confirmations) {
1627     G4cout << "  Viewer now " << fpViewer -> GetName () << G4endl;
1628   }
1629       }
1630       else {
1631   fpViewer = 0;
1632       }
1633     }
1634     else {
1635       fpSceneHandler = 0;
1636       fpViewer = 0;
1637     }
1638   }
1639 }
1640 
1641 void G4VisManager::SetCurrentSceneHandler (G4VSceneHandler* pSceneHandler) {
1642   fpSceneHandler = pSceneHandler;
1643   if (fVerbosity >= confirmations) {
1644     G4cout << "G4VisManager::SetCurrentSceneHandler: scene handler now \""
1645      << pSceneHandler -> GetName () << "\"" << G4endl;
1646   }
1647   if (fpScene != fpSceneHandler -> GetScene ()) {
1648     fpScene = fpSceneHandler -> GetScene ();
1649     if (fVerbosity >= confirmations) {
1650       G4cout << "  Scene now \""
1651        << fpScene -> GetName () << "\"" << G4endl;
1652     }
1653   }
1654   if (fpGraphicsSystem != pSceneHandler -> GetGraphicsSystem ()) {
1655     fpGraphicsSystem = pSceneHandler -> GetGraphicsSystem ();
1656     if (fVerbosity >= confirmations) {
1657       G4cout << "  Graphics system now \""
1658        << fpGraphicsSystem -> GetName () << "\"" << G4endl;
1659     }
1660   }
1661   const G4ViewerList& viewerList = fpSceneHandler -> GetViewerList ();
1662   std::size_t nViewers = viewerList.size ();
1663   if (nViewers) {
1664     std::size_t iViewer;
1665     for (iViewer = 0; iViewer < nViewers; ++iViewer) {
1666       if (fpViewer == viewerList [iViewer]) break;
1667     }
1668     if (iViewer >= nViewers) {
1669       fpViewer = viewerList [0];
1670       if (fVerbosity >= confirmations) {
1671   G4cout << "  Viewer now \"" << fpViewer -> GetName () << "\""
1672          << G4endl;
1673       }
1674     }
1675     if (!IsValidView ()) {
1676       if (fVerbosity >= warnings) {
1677   G4warn <<
1678   "WARNING: Problem setting scene handler - please report circumstances."
1679          << G4endl;
1680       }
1681     }
1682   }
1683   else {
1684     fpViewer = 0;
1685     if (fVerbosity >= warnings) {
1686       G4warn <<
1687   "WARNING: No viewers for this scene handler - please create one."
1688        << G4endl;
1689     }
1690   }
1691 }
1692 
1693 void G4VisManager::SetCurrentViewer (G4VViewer* pViewer) {
1694   fpViewer  = pViewer;
1695   if (fpViewer == nullptr) {
1696     if (fVerbosity >= confirmations) {
1697       G4cout << "G4VisManager::SetCurrentViewer: current viewer pointer zeroed "
1698       << G4endl;
1699     }
1700     return;
1701   }
1702   if (fVerbosity >= confirmations) {
1703     G4cout << "G4VisManager::SetCurrentViewer: viewer now "
1704      << pViewer -> GetName ()
1705      << G4endl;
1706   }
1707   fpSceneHandler = fpViewer -> GetSceneHandler ();
1708   if (!fpSceneHandler) {
1709     if (fVerbosity >= warnings) {
1710       G4warn <<
1711       "WARNING: No scene handler for this viewer - please create one."
1712       << G4endl;
1713     }
1714     return;
1715   }
1716   // JA: I don't think we need this. Setview will be called when needed.
1717   // fpViewer->SetView();
1718   fpSceneHandler -> SetCurrentViewer (pViewer);
1719   fpScene = fpSceneHandler -> GetScene ();
1720   fpGraphicsSystem = fpSceneHandler -> GetGraphicsSystem ();
1721   if (!IsValidView ()) {
1722     if (fVerbosity >= warnings) {
1723       G4warn <<
1724   "WARNING: Problem setting viewer - please report circumstances."
1725        << G4endl;
1726     }
1727   }
1728 }
1729 
1730 void G4VisManager::PrintAvailableGraphicsSystems
1731 (Verbosity verbosity, std::ostream& out) const
1732 {
1733   out << "Registered graphics systems are:\n";
1734   if (fAvailableGraphicsSystems.size ()) {
1735     for (const auto& gs: fAvailableGraphicsSystems) {
1736       const G4String& name = gs->GetName();
1737       const std::vector<G4String>& nicknames = gs->GetNicknames();
1738       if (verbosity >= warnings) {
1739         // Brief output
1740         out << "  " << name << " (";
1741         for (std::size_t i = 0; i < nicknames.size(); ++i) {
1742           if (i != 0) {
1743             out << ", ";
1744           }
1745           out << nicknames[i];
1746         }
1747         out << ')';
1748       } else {
1749         // Full output
1750         out << *gs;
1751       }
1752       out << std::endl;
1753     }
1754     out <<
1755     "You may choose a graphics system (driver) with a parameter of"
1756     "\nthe command \"/vis/open\" or \"/vis/sceneHandler/create\","
1757     "\nor you may omit the driver parameter and choose at run time:"
1758     "\n- by argument in the construction of G4VisExecutive;"
1759     "\n- by environment variable \"G4VIS_DEFAULT_DRIVER\";"
1760     "\n- by entry in \"~/.g4session\";"
1761     "\n- by build flags."
1762     "\n- Note: This feature is not allowed in batch mode."
1763     "\nFor further information see \"examples/basic/B1/exampleB1.cc\""
1764     "\nand \"vis.mac\"."
1765     << std::endl;
1766   } else {
1767     out << "  NONE!!!  None registered - yet!  Mmmmm!" << std::endl;
1768   }
1769 }
1770 
1771 void G4VisManager::PrintAvailableModels (Verbosity verbosity) const
1772 {
1773   {
1774     //fpTrajDrawModelMgr->Print(G4cout);
1775     G4cout << "Registered model factories:" << G4endl;
1776     const std::vector<G4VModelFactory<G4VTrajectoryModel>*>& factoryList =
1777       fpTrajDrawModelMgr->FactoryList();
1778     if (factoryList.empty()) G4cout << "  None" << G4endl;
1779     else {
1780       std::vector<G4VModelFactory<G4VTrajectoryModel>*>::const_iterator i;
1781       for (i = factoryList.begin(); i != factoryList.end(); ++i) {
1782         (*i)->Print(G4cout);
1783       }
1784     }
1785     G4cout << "\nRegistered models:" << G4endl;
1786     const G4VisListManager<G4VTrajectoryModel>* listManager =
1787       fpTrajDrawModelMgr->ListManager();
1788     const std::map<G4String, G4VTrajectoryModel*>& modelMap =
1789       listManager->Map();
1790     if (modelMap.empty()) G4cout << "  None" << G4endl;
1791     else {
1792       std::map<G4String, G4VTrajectoryModel*>::const_iterator i;
1793       for (i = modelMap.begin(); i != modelMap.end(); ++i) {
1794   G4cout << "  " << i->second->Name();
1795   if (i->second == listManager->Current()) G4cout << " (Current)";
1796   G4cout << G4endl;
1797   if (verbosity >= parameters) i->second->Print(G4cout);
1798       }
1799     }
1800   }
1801 
1802   G4cout << G4endl;
1803 
1804   {
1805     //fpTrajFilterMgr->Print(G4cout);
1806     G4cout << "Registered filter factories:" << G4endl;
1807     const std::vector<G4VModelFactory<G4VFilter<G4VTrajectory> >*>&
1808       factoryList = fpTrajFilterMgr->FactoryList();
1809     if (factoryList.empty()) G4cout << "  None" << G4endl;
1810     else {
1811       std::vector<G4VModelFactory<G4VFilter<G4VTrajectory> >*>::const_iterator i;
1812       for (i = factoryList.begin(); i != factoryList.end(); ++i) {
1813         (*i)->Print(G4cout);
1814       }
1815     }
1816 
1817     G4cout << "\nRegistered filters:" << G4endl;
1818     const std::vector<G4VFilter<G4VTrajectory>*>&
1819       filterList = fpTrajFilterMgr->FilterList();
1820     if (filterList.empty()) G4cout << "  None" << G4endl;
1821     else {
1822       std::vector<G4VFilter<G4VTrajectory>*>::const_iterator i;
1823       for (i = filterList.begin(); i != filterList.end(); ++i) {
1824   G4cout << "  " << (*i)->GetName() << G4endl;
1825   if (verbosity >= parameters) (*i)->PrintAll(G4cout);
1826       }
1827     }
1828   }
1829 }
1830 
1831 void G4VisManager::PrintAvailableUserVisActions (Verbosity) const
1832 {
1833   G4cout <<
1834     "You have successfully registered the following user vis actions."
1835    << G4endl;
1836   G4cout << "Run Duration User Vis Actions:";
1837   if (fRunDurationUserVisActions.empty()) G4cout << " none" << G4endl;
1838   else {
1839     G4cout << G4endl;
1840     for (std::size_t i = 0; i < fRunDurationUserVisActions.size(); ++i) {
1841       const G4String& name = fRunDurationUserVisActions[i].fName;
1842       G4cout << "  " << name << G4endl;
1843     }
1844   }
1845 
1846   G4cout << "End of Event User Vis Actions:";
1847   if (fEndOfEventUserVisActions.empty()) G4cout << " none" << G4endl;
1848   else {
1849     G4cout << G4endl;
1850     for (std::size_t i = 0; i < fEndOfEventUserVisActions.size(); ++i) {
1851       const G4String& name = fEndOfEventUserVisActions[i].fName;
1852       G4cout << "  " << name << G4endl;
1853     }
1854   }
1855 
1856   G4cout << "End of Run User Vis Actions:";
1857   if (fEndOfRunUserVisActions.empty()) G4cout << " none" << G4endl;
1858   else {
1859     G4cout << G4endl;
1860     for (std::size_t i = 0; i < fEndOfRunUserVisActions.size(); ++i) {
1861       const G4String& name = fEndOfRunUserVisActions[i].fName;
1862       G4cout << "  " << name << G4endl;
1863     }
1864   }
1865 }
1866 
1867 void G4VisManager::PrintAvailableColours (Verbosity) const {
1868   G4cout <<
1869   "Some /vis commands (optionally) take a string to specify colour."
1870   "\nThey are also available in your C++ code, e.g:"
1871   "\n  G4Colour niceColour;  // Default - white"
1872   "\n  G4Colour::GetColour(\"pink\", niceColour);"
1873   "\n  logical->SetVisAttributes(niceColour);"
1874   "\nSee G4Colour.hh."
1875   "\nAvailable colours";
1876   for (const auto& i: G4Colour::GetMap()) {
1877     G4cout << ", " << i.first;
1878   }
1879   G4cout << G4endl;
1880 }
1881 
1882 void G4VisManager::PrintInvalidPointers () const {
1883   if (fVerbosity >= errors) {
1884     G4warn << "ERROR: G4VisManager::PrintInvalidPointers:";
1885     if (!fpGraphicsSystem) {
1886       G4warn << "\n null graphics system pointer.";
1887     }
1888     else {
1889       G4warn << "\n  Graphics system is " << fpGraphicsSystem -> GetName ()
1890        << " but:";
1891       if (!fpScene)
1892   G4warn <<
1893     "\n  Null scene pointer. Use \"/vis/drawVolume\" or"
1894     " \"/vis/scene/create\".";
1895       if (!fpSceneHandler)
1896   G4warn <<
1897     "\n  Null scene handler pointer. Use \"/vis/open\" or"
1898     " \"/vis/sceneHandler/create\".";
1899       if (!fpViewer )
1900   G4warn <<
1901     "\n  Null viewer pointer. Use \"/vis/viewer/create\".";
1902     }
1903     G4warn << G4endl;
1904   }
1905 }
1906 
1907 
1908 G4ThreadFunReturnType G4VisManager::G4VisSubThread(G4ThreadFunArgType p)
1909 {
1910 #ifdef G4MULTITHREADED
1911   G4VisManager* pVisManager = (G4VisManager*)p;
1912   G4VSceneHandler* pSceneHandler = pVisManager->GetCurrentSceneHandler();
1913   if (!pSceneHandler) return 0;
1914   G4Scene* pScene = pSceneHandler->GetScene();
1915   if (!pScene) return 0;
1916   G4VViewer* pViewer = pVisManager->GetCurrentViewer();
1917   if (!pViewer) return 0;
1918 
1919   G4UImanager::GetUIpointer()->SetUpForSpecialThread("G4VIS");
1920 
1921   // Set up geometry and navigation for a thread
1922   G4GeometryWorkspace::GetPool()->CreateAndUseWorkspace();
1923   G4SolidsWorkspace::GetPool()->CreateAndUseWorkspace();
1924   G4Navigator* navigator = G4TransportationManager::GetTransportationManager()
1925                              ->GetNavigatorForTracking();
1926   navigator->SetWorldVolume(
1927     G4RunManagerFactory::GetMasterRunManagerKernel()->GetCurrentWorld());
1928 
1929   pViewer->SwitchToVisSubThread();
1930 
1931   while (true) {
1932 
1933     G4MUTEXLOCK(&mtVisSubThreadMutex);
1934     std::size_t eventQueueSize = mtVisEventQueue.size();
1935     G4MUTEXUNLOCK(&mtVisSubThreadMutex);
1936     // G4cout << "Event queue size (A): " << eventQueueSize << G4endl;
1937 
1938     while (eventQueueSize) {
1939 
1940       G4MUTEXLOCK(&mtVisSubThreadMutex);
1941       const G4Event* event = mtVisEventQueue.front();
1942       G4MUTEXUNLOCK(&mtVisSubThreadMutex);
1943 
1944       // Here comes the event drawing
1945       pVisManager->SetTransientsDrawnThisEvent(false);
1946       pSceneHandler->SetTransientsDrawnThisEvent(false);
1947 
1948       // We are about to draw the event (trajectories, etc.), but first we
1949       // have to clear the previous event(s) if necessary.  If this event
1950       // needs to be drawn afresh, e.g., the first event or any event when
1951       // "accumulate" is not requested, the old event has to be cleared.
1952       // We have postponed this so that, for normal viewers like OGL, the
1953       // previous event(s) stay on screen until this new event comes
1954       // along.  For a file-writing viewer the geometry has to be drawn.
1955       // See, for example, G4HepRepFileSceneHandler::ClearTransientStore.
1956       pVisManager->ClearTransientStoreIfMarked();
1957 
1958       // Now draw the event...
1959       pSceneHandler->DrawEvent(event);
1960       ++pVisManager->fNoOfEventsDrawnThisRun;
1961 
1962       // Then pop and release event
1963       G4MUTEXLOCK(&mtVisSubThreadMutex);
1964       mtVisEventQueue.pop_front();
1965       pVisManager->EndOfEventCleanup(event);
1966       eventQueueSize = mtVisEventQueue.size();
1967       G4MUTEXUNLOCK(&mtVisSubThreadMutex);
1968     }
1969 
1970     G4MUTEXLOCK(&mtVisSubThreadMutex);
1971     G4bool runInProgress = mtRunInProgress;
1972     G4MUTEXUNLOCK(&mtVisSubThreadMutex);
1973     if (!runInProgress) {
1974       // EndOfRun on master thread has signalled end of run.  There is
1975       // nothing to draw so...
1976       break;
1977     }
1978 
1979     // Run still in progress but nothing to draw, so wait a while.
1980     std::this_thread::sleep_for(std::chrono::milliseconds(100));
1981   }
1982 
1983   // Inform viewer that we have finished all sub-thread drawing
1984   pViewer->DoneWithVisSubThread();
1985   pViewer->MovingToMasterThread();
1986 #else
1987   G4ConsumeParameters(p);
1988 #endif
1989   return nullptr;
1990 }
1991 
1992 void G4VisManager::BeginOfRun ()
1993 {
1994   if (fIgnoreStateChanges) return;
1995   if (G4Threading::IsWorkerThread()) return;
1996   if (!GetConcreteInstance()) return;
1997 
1998   // Check if view is valid
1999   // Protect IsValidView() with fpSceneHandler to prevent warning messages in batch
2000   isValidViewForRun = fpSceneHandler && IsValidView();
2001   if (!isValidViewForRun) return;
2002 
2003   // For a fake run...
2004   G4RunManager* runManager = G4RunManagerFactory::GetMasterRunManager();
2005   G4int nEventsToBeProcessed = runManager->GetNumberOfEventsToBeProcessed();
2006   isFakeRun = (nEventsToBeProcessed == 0);
2007   if (isFakeRun) return;
2008 
2009   // Is the run manager a sub-event type?
2010   G4RunManager::RMType rmType = runManager->GetRunManagerType();
2011   isSubEventRunManagerType =
2012   (rmType == G4RunManager::subEventMasterRM) ||
2013   (rmType == G4RunManager::subEventWorkerRM);
2014 
2015   fNKeepForPostProcessingRequests = 0;
2016   fNKeepTheEventRequests = 0;
2017   fEventKeepingSuspended = false;
2018   fTransientsDrawnThisRun = false;
2019   if (fpSceneHandler) fpSceneHandler->SetTransientsDrawnThisRun(false);
2020   fNoOfEventsDrawnThisRun = 0;
2021 
2022   // Check to see if the user has created a trajectory model. If not, create
2023   // a default one. To avoid code duplication the following function is used
2024   // and its result (a const G4VTrajectoryModel*) is thrown away at this point.
2025   // The function is called again later when needed.
2026   CurrentTrajDrawModel();
2027 
2028   // Actions only in MT mode
2029   if (G4Threading::IsMultithreadedApplication()) {
2030 
2031     // Inform viewer that we have finished all master thread drawing for now...
2032     if (fpViewer) fpViewer->DoneWithMasterThread();
2033 
2034     // Start vis sub-thread
2035     G4MUTEXLOCK(&mtVisSubThreadMutex);
2036     mtRunInProgress = true;
2037     G4MUTEXUNLOCK(&mtVisSubThreadMutex);
2038     mtVisSubThread = new G4Thread;
2039     // Launch vis thread
2040     G4THREADCREATE(mtVisSubThread,G4VisSubThread,this);
2041 
2042     // Tricky things for some viewers (e.g., Qt):
2043     // - Launch the vis thread
2044     // - Wait for the vis thread to set its QThread
2045     // - Then move current QOpenGL context (if Qt) to this Qthread
2046     // - Go ahead
2047     if (fpViewer) fpViewer->MovingToVisSubThread();
2048   }
2049 }
2050 
2051 void G4VisManager::BeginOfEvent ()
2052 {
2053   if (fIgnoreStateChanges) return;
2054   if (!GetConcreteInstance()) return;
2055   if (!isValidViewForRun) return;
2056   if (isFakeRun) return;
2057 
2058   // Some instructions that should NOT be in multithreaded version.
2059   // TODO: protect with G4Threading::IsMultithreadedApplication() instead?
2060 #ifndef G4MULTITHREADED
2061   // These instructions are in G4VisSubThread for multithreading.
2062   fTransientsDrawnThisEvent = false;
2063   if (fpSceneHandler) fpSceneHandler->SetTransientsDrawnThisEvent(false);
2064 #endif
2065 }
2066 
2067 // Here begins a sequence of functions that deal with end-of-event.
2068 // Sequential/Serial mode:
2069 //   EndOfEvent is invoked by a state change on the master thread:
2070 //     EndOfEvent pulls the event from the Event Manager and calls EndOfEventKernel.
2071 //     EndOfEventKernel draws the event and calls EndOfEventCleanup.
2072 //     (Unless the user him/herself has requested keeping, EndOfEventKernel
2073 //      also sets KeepTheEvent for a certain number (adjustable, default 100)
2074 //      of events. The run manager keeps these events until the beginning of the next
2075 //      run, so, for example, the vis manager can redraw/rebuild if required.)
2076 // Multithreading/Tasking:
2077 //   EndOfEvent is invoked by a state change on each worker thread:
2078 //     EndOfEvent pulls the event from the (worker) Event Manager and calls EndOfEventKernel.
2079 //     EndOfEventKernel pushes it to the vis queue, and returns. The events
2080 //       are marked KeepForPostProcessing so that the worker does not delete them until
2081 //       PostProcessingFinished in EndOfEventCleanup.
2082 //     EndOfEventKernel also sets KeepTheEvent as above.
2083 //     The vis sub thread pulls events from the vis queue, draws them and
2084 //       calls EndOfEventCleanup, which calls PostProcessingFinished.
2085 // Sub-event parallelism:
2086 //   EndOfEvent is invoked by a state change on each worker thread, as for
2087 //     multithreading, but actually, the event is not complete. Workers are
2088 //     employed on sub-events, so we ignore such state changes. Only
2089 //     G4SubEvtRunManager knows when the event is complete.
2090 //   Events are pushed to the Vis Manager by G4SubEvtRunManager via EventReadyForVis
2091 //      when the event is ready. EventReadyForVis calls EndOfEventKernel and
2092 //      processsing continues as for Multithreading.
2093 
2094 void G4VisManager::EndOfEvent ()
2095 {
2096   if (fIgnoreStateChanges) return;
2097   if (!GetConcreteInstance()) return;
2098   if (!isValidViewForRun) return;
2099   if (isFakeRun) return;
2100 
2101   // For sub-event parallelism, this is not the true end of event
2102   if (isSubEventRunManagerType) return;
2103 
2104   G4AutoLock al(&visEndOfEventMutex);
2105 
2106   G4RunManager* runManager = G4RunManagerFactory::GetMasterRunManager();
2107 
2108   const G4Run* currentRun = runManager->GetCurrentRun();
2109   if (!currentRun) return;
2110 
2111   // This gets the appropriate event manager
2112   G4EventManager* eventManager = G4EventManager::GetEventManager();
2113   const G4Event* currentEvent = eventManager->GetConstCurrentEvent();
2114   if (!currentEvent) return;
2115 
2116   // Call kernel
2117   EndOfEventKernel(currentEvent);
2118 }
2119 
2120 void G4VisManager::EventReadyForVis(const G4Event* event)
2121 // This is invoked by G4SubEvtRunManager.
2122 // The event is passed to EndOfEventKernel.
2123 {
2124   if (fIgnoreStateChanges) return;
2125   if (!GetConcreteInstance()) return;
2126   if (!isValidViewForRun) return;
2127   if (isFakeRun) return;
2128 
2129   G4AutoLock al(&visEndOfEventMutex);
2130   EndOfEventKernel(event);
2131 }
2132 
2133 void G4VisManager::EndOfEventKernel (const G4Event* currentEvent)
2134 {
2135   // Note: we are still subject to the mutex lock with visEndOfEventMutex
2136 
2137   // Discard event if fDrawEventOnlyIfToBeKept flag is set unless the
2138   // user has requested the event to be kept.
2139   if (fDrawEventOnlyIfToBeKept) {
2140     if (!currentEvent->KeepTheEventFlag()) return;
2141   }
2142 
2143   if (G4Threading::IsMultithreadedApplication()) {
2144 
2145     // Wait if too many events in the queue.
2146     G4MUTEXLOCK(&mtVisSubThreadMutex);
2147     std::size_t eventQueueSize = mtVisEventQueue.size();
2148     G4MUTEXUNLOCK(&mtVisSubThreadMutex);
2149 
2150     G4bool eventQueueFull = false;
2151     while (fMaxEventQueueSize > 0 && (G4int)eventQueueSize >= fMaxEventQueueSize) {
2152 
2153       if (fWaitOnEventQueueFull) {
2154 
2155         static G4bool warned = false;
2156         if (!warned) {
2157           G4warn <<
2158           "WARNING: The number of events in the visualisation queue has exceeded"
2159           "\n  the maximum, "
2160           << fMaxEventQueueSize <<
2161           ".\n  If, during a multithreaded run, the simulation gets ahead of the"
2162           "\n  visualisation by more than this maximum, the simulation is delayed"
2163           "\n  until the vis sub-thread has drawn a few more events and removed them"
2164           "\n  from the queue.  You may change this maximum number of events with"
2165           "\n  \"/vis/multithreading/maxEventQueueSize <N>\", where N is the maximum"
2166           "\n  number you wish to allow.  N <= 0 means \"unlimited\"."
2167           "\n  Alternatively you may choose to discard events for drawing by setting"
2168           "\n  \"/vis/multithreading/actionOnEventQueueFull discard\"."
2169           "\n  To avoid visualisation altogether: \"/vis/disable\"."
2170           "\n  And maybe \"/tracking/storeTrajectories 0\"."
2171           << G4endl;
2172           warned = true;
2173         }
2174 
2175         // Wait a while to give event drawing time to reduce the queue...
2176         std::this_thread::sleep_for(std::chrono::milliseconds(100));
2177 
2178       } else {
2179 
2180         static G4bool warned = false;
2181         if (!warned) {
2182           G4warn <<
2183           "WARNING: The number of events in the visualisation queue has exceeded"
2184           "\n  the maximum, "
2185           << fMaxEventQueueSize <<
2186           ".\n  Some events have been discarded for drawing.  You may change this"
2187           "\n  behaviour with \"/vis/multithreading/actionOnEventQueueFull wait\"."
2188           "\n  To avoid visualisation altogether: \"/vis/disable\"."
2189           "\n  And maybe \"/tracking/storeTrajectories 0\"."
2190           << G4endl;
2191           warned = true;
2192         }
2193 
2194         eventQueueFull = true;  // Causes event to be discarded for drawing.
2195         break;
2196       }
2197 
2198       G4MUTEXLOCK(&mtVisSubThreadMutex);
2199       eventQueueSize = mtVisEventQueue.size();
2200       G4MUTEXUNLOCK(&mtVisSubThreadMutex);
2201     }
2202 
2203     if (!eventQueueFull) {
2204 
2205       if (RequiredToBeKeptForVis(currentEvent->GetEventID())) {
2206         currentEvent->KeepTheEvent();
2207         fNKeepTheEventRequests++;  // Counts number of calls to KeepTheEvent
2208       }
2209 
2210       G4MUTEXLOCK(&mtVisSubThreadMutex);
2211 
2212       // Keep while post processing (i.e., drawing), relinquished in EndOfEventCleanup
2213       // Make sure the event is not deleted by the run manager while in the vis queue
2214       currentEvent->KeepForPostProcessing();
2215       fNKeepForPostProcessingRequests++;
2216 
2217       // Put event on vis queue
2218       mtVisEventQueue.push_back(currentEvent);
2219 
2220       G4MUTEXUNLOCK(&mtVisSubThreadMutex);
2221     }
2222 
2223   } else {
2224 
2225     // Sequential mode
2226 
2227     // We are about to draw the event (trajectories, etc.), but first we
2228     // have to clear the previous event(s) if necessary.  If this event
2229     // needs to be drawn afresh, e.g., the first event or any event when
2230     // "accumulate" is not requested, the old event has to be cleared.
2231     // We have postponed this so that, for normal viewers like OGL, the
2232     // previous event(s) stay on screen until this new event comes
2233     // along.  For a file-writing viewer the geometry has to be drawn.
2234     // See, for example, G4HepRepFileSceneHandler::ClearTransientStore.
2235     ClearTransientStoreIfMarked();
2236 
2237     // Keep while post processing (i.e., drawing), relinquished in EndOfEventCleanup
2238     currentEvent->KeepForPostProcessing();
2239     fNKeepForPostProcessingRequests++;
2240 
2241     if (RequiredToBeKeptForVis(currentEvent->GetEventID())) {
2242       currentEvent->KeepTheEvent();
2243       fNKeepTheEventRequests++;  // Counts number of calls to KeepTheEvent
2244     }
2245 
2246     // Now draw the event...
2247     fpSceneHandler->DrawEvent(currentEvent);
2248     ++fNoOfEventsDrawnThisRun;
2249 
2250     EndOfEventCleanup(currentEvent);
2251   }
2252 }
2253 
2254 void G4VisManager::EndOfEventCleanup(const G4Event* currentEvent)
2255 {
2256   if (fpScene->GetRefreshAtEndOfEvent()) {
2257 
2258     fpSceneHandler->SetMarkForClearingTransientStore(true);
2259 
2260     // ShowView guarantees the view is flushed to the screen.  It also
2261     // triggers other features such as picking (if enabled) and allows
2262     // file-writing viewers to close the file.
2263     fpViewer->ShowView();
2264   }
2265 
2266   currentEvent->PostProcessingFinished();
2267 }
2268 
2269 G4bool G4VisManager::RequiredToBeKeptForVis (G4int eventID)
2270 // i.e., kept by the run manager at least to the beginning of the next run.
2271 {
2272   G4bool requiredToBeKept = false;
2273 
2274   G4RunManager* runManager = G4RunManagerFactory::GetMasterRunManager();
2275   G4int nEventsToBeProcessed = runManager->GetNumberOfEventsToBeProcessed();
2276 
2277   if (fpScene->GetRefreshAtEndOfEvent()) {  // Refreshing every event
2278 
2279     // Keep last event only
2280     if (eventID == nEventsToBeProcessed - 1) {
2281       requiredToBeKept = true;
2282     }
2283 
2284   } else {  //  Accumulating events
2285 
2286     G4int maxNumberOfKeptEvents = fpScene->GetMaxNumberOfKeptEvents();
2287 
2288     if (maxNumberOfKeptEvents >= 0) {  // Event keeping requested
2289 
2290       if (fNKeepTheEventRequests < maxNumberOfKeptEvents) {
2291 
2292         // Not yet reached the limit
2293         requiredToBeKept = true;
2294 
2295       } else {
2296 
2297         // We have already reached the limit
2298         fEventKeepingSuspended = true;
2299         static G4bool warned = false;
2300         if (!warned) {
2301           if (fVerbosity >= warnings) {
2302             G4warn <<
2303             "WARNING: G4VisManager::EndOfEvent: Automatic event keeping suspended."
2304             << G4endl;
2305             if (maxNumberOfKeptEvents > 0) {
2306               G4warn <<
2307               "\n  The number of events exceeds the maximum, "
2308               << maxNumberOfKeptEvents <<
2309               ", that may be kept by\n  the vis manager."
2310               << G4endl;
2311             }
2312           }
2313           warned = true;
2314         }
2315       }
2316 
2317     } else {  // Indefinite event keeping (maxNumberOfKeptEvents < 0)
2318 
2319       requiredToBeKept = true;
2320 
2321     }
2322   }
2323 
2324   return requiredToBeKept;
2325 }
2326 
2327 void G4VisManager::EndOfRun ()
2328 {
2329   if (fIgnoreStateChanges) return;
2330   if (G4Threading::IsWorkerThread()) return;
2331   if (!GetConcreteInstance()) return;
2332   if (!isValidViewForRun) return;
2333   if (isFakeRun) return;
2334 
2335   G4RunManager* runManager = G4RunManagerFactory::GetMasterRunManager();
2336 
2337   // For a fake run...
2338   G4int nEventsToBeProcessed = runManager->GetNumberOfEventsToBeProcessed();
2339   if (nEventsToBeProcessed == 0) return;
2340 
2341   const G4Run* currentRun = runManager->GetCurrentRun();
2342   if (!currentRun) return;
2343 
2344   //  G4AutoLock al(&visEndOfRunMutex);  ???
2345   if (G4Threading::IsMultithreadedApplication()) {
2346     // Reset flag so that sub-thread exits when it has finished processing.
2347     G4MUTEXLOCK(&mtVisSubThreadMutex);
2348     mtRunInProgress = false;
2349     G4MUTEXUNLOCK(&mtVisSubThreadMutex);
2350     // Wait for sub-thread to finish.
2351     G4THREADJOIN(*mtVisSubThread);
2352     delete mtVisSubThread;
2353     fpViewer->SwitchToMasterThread();
2354   }
2355 
2356   if (G4Threading::IsMultithreadedApplication()) {
2357     G4int noOfEventsRequested = runManager->GetNumberOfEventsToBeProcessed();
2358     if (fNoOfEventsDrawnThisRun != noOfEventsRequested) {
2359       if (!fWaitOnEventQueueFull && fVerbosity >= warnings) {
2360         G4warn
2361         << "WARNING: Number of events drawn this run, "
2362         << fNoOfEventsDrawnThisRun << ", is different to number requested, "
2363         << noOfEventsRequested <<
2364         ".\n  (This is because you requested \"/vis/multithreading/actionOnEventQueueFull discard\".)"
2365         << G4endl;
2366       }
2367     }
2368   }
2369 
2370   G4int nKeptEvents = currentRun->GetNumberOfKeptEvents();
2371   if (fVerbosity >= warnings && nKeptEvents > 0) {
2372     G4warn << nKeptEvents;
2373     if (nKeptEvents == 1) G4warn << " event has";
2374     else G4warn << " events have";
2375     G4warn << " been kept for refreshing and/or reviewing." << G4endl;
2376     if (nKeptEvents != fNKeepTheEventRequests) {
2377       if (fNKeepTheEventRequests == 0) {
2378         G4warn << "No keep requests were";
2379       } else if (fNKeepTheEventRequests == 1) {
2380         G4warn << "1 keep request was";
2381       } else {
2382         G4warn << fNKeepTheEventRequests << " keep requests were";
2383       }
2384       G4warn << " made by the vis manager.";
2385       if (fNKeepTheEventRequests == 0) {
2386         G4warn <<
2387         "\n  The kept events are those you have asked to be kept in your user action(s).";
2388       } else {
2389         G4warn <<
2390         "\n  The same or further events may have been kept by you in your user action(s)."
2391         "\n  To turn off event keeping by the vis manager: /vis/drawOnlyToBeKeptEvents"
2392         "\n  or use /vis/scene/endOfEventAction <refresh|accumulate> 0";
2393       }
2394       G4warn << G4endl;
2395     }
2396     G4warn <<
2397   "  \"/vis/reviewKeptEvents\" to review one by one."
2398   "\n  To see accumulated, \"/vis/enable\", then \"/vis/viewer/flush\" or \"/vis/viewer/rebuild\"."
2399     << G4endl;
2400   }
2401 
2402   if (fVerbosity >= warnings) PrintListOfPlots();
2403 
2404   if (fEventKeepingSuspended && fVerbosity >= warnings) {
2405     G4warn <<
2406     "WARNING: G4VisManager::EndOfRun: Automatic event keeping was suspended."
2407     << G4endl;
2408     if (fpScene->GetMaxNumberOfKeptEvents() > 0) {
2409       G4warn <<
2410       "  The number of events in the run exceeded the maximum, "
2411       << fpScene->GetMaxNumberOfKeptEvents() <<
2412       ", that may be\n  kept by the vis manager." <<
2413       "\n  The number of events kept by the vis manager can be changed with"
2414       "\n  \"/vis/scene/endOfEventAction accumulate <N>\", where N is the"
2415       "\n  maximum number you wish to allow.  N < 0 means \"unlimited\"."
2416       << G4endl;
2417     }
2418   }
2419 
2420 //    // if (!fpSceneHandler->GetMarkForClearingTransientStore()) {
2421 //    // is here.  It prevents ShowView at end of run, which seems to be OK
2422 //    // for sequential mode, but MT mode seems to need it (I have not
2423 //    // figured out why). ???? JA ????
2424 //    if (!fpSceneHandler->GetMarkForClearingTransientStore()) {
2425   if (fpScene->GetRefreshAtEndOfRun()) {
2426     fpSceneHandler->DrawEndOfRunModels();
2427     // An extra refresh for auto-refresh viewers.
2428     // ???? I DON'T WHY THIS IS NECESSARY ???? JA ????
2429     if (fpViewer->GetViewParameters().IsAutoRefresh()) {
2430       fpViewer->RefreshView();
2431     }
2432     // ShowView guarantees the view is flushed to the screen.  It also
2433     // triggers other features such picking (if enabled) and allows
2434     // file-writing viewers to close the file.
2435     fpViewer->ShowView();
2436     fpSceneHandler->SetMarkForClearingTransientStore(true);
2437   } else {
2438     if (fpGraphicsSystem->GetFunctionality() ==
2439         G4VGraphicsSystem::fileWriter) {
2440       if (fVerbosity >= warnings) {
2441         G4warn << "\"/vis/viewer/update\" to close file." << G4endl;
2442       }
2443     }
2444   }
2445 //}
2446   fEventRefreshing = false;
2447 }
2448 
2449 void G4VisManager::ClearTransientStoreIfMarked(){
2450   // Assumes valid view.
2451   if (fpSceneHandler->GetMarkForClearingTransientStore()) {
2452     fpSceneHandler->SetMarkForClearingTransientStore(false);
2453     fpSceneHandler->ClearTransientStore();
2454   }
2455   // Record if transients drawn.  These local flags are only set
2456   // *after* ClearTransientStore.  In the code in G4VSceneHandler
2457   // triggered by ClearTransientStore, use these flags so that
2458   // event refreshing is not done too early.
2459   fTransientsDrawnThisEvent = fpSceneHandler->GetTransientsDrawnThisEvent();
2460   fTransientsDrawnThisRun = fpSceneHandler->GetTransientsDrawnThisRun();
2461 }
2462 
2463 void G4VisManager::ResetTransientsDrawnFlags()
2464 {
2465   fTransientsDrawnThisRun = false;
2466   fTransientsDrawnThisEvent = false;
2467   G4SceneHandlerListConstIterator i;
2468   for (i = fAvailableSceneHandlers.begin();
2469        i != fAvailableSceneHandlers.end(); ++i) {
2470     (*i)->SetTransientsDrawnThisEvent(false);
2471     (*i)->SetTransientsDrawnThisRun(false);
2472   }
2473 }
2474 
2475 G4String G4VisManager::ViewerShortName (const G4String& viewerName) const {
2476   const G4String& viewerShortName = viewerName.substr(0, viewerName.find (' '));
2477   return G4StrUtil::strip_copy(viewerShortName);
2478 }
2479 
2480 G4VViewer* G4VisManager::GetViewer (const G4String& viewerName) const {
2481   G4String viewerShortName = ViewerShortName (viewerName);
2482   std::size_t nHandlers = fAvailableSceneHandlers.size ();
2483   std::size_t iHandler, iViewer;
2484   G4VViewer* viewer = 0;
2485   G4bool found = false;
2486   for (iHandler = 0; iHandler < nHandlers; iHandler++) {
2487     G4VSceneHandler* sceneHandler = fAvailableSceneHandlers [iHandler];
2488     const G4ViewerList& viewerList = sceneHandler -> GetViewerList ();
2489     for (iViewer = 0; iViewer < viewerList.size (); iViewer++) {
2490       viewer = viewerList [iViewer];
2491       if (viewerShortName == viewer -> GetShortName ()) {
2492   found = true;
2493   break;
2494       }
2495     }
2496     if (found) break;
2497   }
2498   if (found) return viewer;
2499   else return 0;
2500 }
2501 
2502 std::vector<G4String> G4VisManager::VerbosityGuidanceStrings;
2503 
2504 G4String G4VisManager::VerbosityString(Verbosity verbosity) {
2505   G4String rs;
2506   switch (verbosity) {
2507   case         quiet: rs = "quiet (0)"; break;
2508   case       startup: rs = "startup (1)"; break;
2509   case        errors: rs = "errors (2)"; break;
2510   case      warnings: rs = "warnings (3)"; break;
2511   case confirmations: rs = "confirmations (4)"; break;
2512   case    parameters: rs = "parameters (5)"; break;
2513   case           all: rs = "all (6)"; break;
2514   }
2515   return rs;
2516 }
2517 
2518 void G4VisManager::PrintAvailableVerbosity(std::ostream& os)
2519 {
2520   os << "Available verbosity options:";
2521   for (std::size_t i = 0; i < VerbosityGuidanceStrings.size(); ++i) {
2522     os << '\n' << VerbosityGuidanceStrings[i];
2523   }
2524   os << "\nCurrent verbosity: " << G4VisManager::VerbosityString(fVerbosity);
2525   os << std::endl;
2526 }
2527 
2528 G4VisManager::Verbosity
2529 G4VisManager::GetVerbosityValue(const G4String& verbosityString) {
2530   G4String ss = G4StrUtil::to_lower_copy(verbosityString);
2531   Verbosity verbosity;
2532   if      (ss[0] == 'q') verbosity = quiet;
2533   else if (ss[0] == 's') verbosity = startup;
2534   else if (ss[0] == 'e') verbosity = errors;
2535   else if (ss[0] == 'w') verbosity = warnings;
2536   else if (ss[0] == 'c') verbosity = confirmations;
2537   else if (ss[0] == 'p') verbosity = parameters;
2538   else if (ss[0] == 'a') verbosity = all;
2539   else {
2540     // Could be an integer
2541     G4int intVerbosity;
2542     std::istringstream is(ss);
2543     is >> intVerbosity;
2544     if (!is) {
2545       G4warn << "ERROR: G4VisManager::GetVerbosityValue: invalid verbosity \""
2546        << verbosityString << "\"\n";
2547       PrintAvailableVerbosity(G4warn);
2548       // Return existing verbosity
2549       return fVerbosity;
2550     }
2551     else {
2552       verbosity = GetVerbosityValue(intVerbosity);
2553     }
2554   }
2555   return verbosity;
2556 }
2557 
2558 G4VisManager::Verbosity G4VisManager::GetVerbosityValue(G4int intVerbosity) {
2559   Verbosity verbosity;
2560   if      (intVerbosity < quiet) verbosity = quiet;
2561   else if (intVerbosity > all)   verbosity = all;
2562   else                           verbosity = Verbosity(intVerbosity);
2563   return verbosity;
2564 }
2565 
2566 G4VisManager::Verbosity G4VisManager::GetVerbosity () {
2567   return fVerbosity;
2568 }
2569 
2570 void G4VisManager::SetVerboseLevel (G4int intVerbosity) {
2571   fVerbosity = GetVerbosityValue(intVerbosity);
2572 }
2573 
2574 void G4VisManager::SetVerboseLevel (const G4String& verbosityString) {
2575   fVerbosity = GetVerbosityValue(verbosityString);
2576 }
2577 
2578 G4bool G4VisManager::IsValidView () {
2579 
2580   if (!fInitialised) Initialise ();
2581 
2582   static G4bool noGSPrinting = true;
2583   if (!fpGraphicsSystem) {
2584     // Limit printing - we do not want printing if the user simply does
2585     // not want to use graphics, e.g., in batch mode.
2586     if (noGSPrinting) {
2587       noGSPrinting = false;
2588       if (fVerbosity >= warnings) {
2589   G4warn <<
2590   "WARNING: G4VisManager::IsValidView(): Attempt to draw when no graphics system"
2591   "\n  has been instantiated.  Use \"/vis/open\" or \"/vis/sceneHandler/create\"."
2592   "\n  Alternatively, to avoid this message, suppress instantiation of vis"
2593   "\n  manager (G4VisExecutive) and ensure drawing code is executed only if"
2594   "\n  G4VVisManager::GetConcreteInstance() is non-zero."
2595         << G4endl;
2596       }
2597     }
2598     return false;
2599   }
2600 
2601   if ((!fpScene) || (!fpSceneHandler) || (!fpViewer)) {
2602     if (fVerbosity >= errors) {
2603       G4warn <<
2604   "ERROR: G4VisManager::IsValidView(): Current view is not valid."
2605        << G4endl;
2606       PrintInvalidPointers ();
2607     }
2608     return false;
2609   }
2610 
2611   if (fpScene != fpSceneHandler -> GetScene ()) {
2612     if (fVerbosity >= errors) {
2613       G4warn << "ERROR: G4VisManager::IsValidView ():";
2614       if (fpSceneHandler -> GetScene ()) {
2615   G4warn <<
2616     "\n  The current scene \""
2617          << fpScene -> GetName ()
2618          << "\" is not handled by"
2619     "\n  the current scene handler \""
2620          << fpSceneHandler -> GetName ()
2621          << "\""
2622     "\n  (it currently handles scene \""
2623          << fpSceneHandler -> GetScene () -> GetName ()
2624          << "\")."
2625     "\n  Either:"
2626     "\n  (a) attach it to the scene handler with"
2627     "\n      /vis/sceneHandler/attach "
2628          << fpScene -> GetName ()
2629          << ", or"
2630     "\n  (b) create a new scene handler with "
2631     "\n      /vis/sceneHandler/create <graphics-system>,"
2632     "\n      in which case it should pick up the the new scene."
2633          << G4endl;
2634       }
2635       else {
2636   G4warn << "\n  Scene handler \""
2637          << fpSceneHandler -> GetName ()
2638          << "\" has null scene pointer."
2639     "\n  Attach a scene with /vis/sceneHandler/attach [<scene-name>]"
2640          << G4endl;
2641       }
2642     }
2643     return false;
2644   }
2645 
2646   const G4ViewerList& viewerList = fpSceneHandler -> GetViewerList ();
2647   if (viewerList.size () == 0) {
2648     if (fVerbosity >= errors) {
2649       G4warn <<
2650   "ERROR: G4VisManager::IsValidView (): the current scene handler\n  \""
2651        << fpSceneHandler -> GetName ()
2652        << "\" has no viewers.  Do /vis/viewer/create."
2653        << G4endl;
2654     }
2655     return false;
2656   }
2657 
2658   G4bool isValid = true;
2659   if (fpScene -> IsEmpty ()) {  // Add world by default if possible...
2660     G4bool warn(fVerbosity >= warnings);
2661     G4bool successful = fpScene -> AddWorldIfEmpty (warn);
2662     if (!successful || fpScene -> IsEmpty ()) {        // If still empty...
2663       if (fVerbosity >= errors) {
2664   G4warn << "ERROR: G4VisManager::IsValidView ():";
2665   G4warn <<
2666     "\n  Attempt at some drawing operation when scene is empty."
2667     "\n  Maybe the geometry has not yet been defined."
2668     "  Try /run/initialize."
2669           "\n  Or use \"/vis/scene/add/extent\"."
2670          << G4endl;
2671       }
2672       isValid = false;
2673     }
2674     else {
2675       G4UImanager::GetUIpointer()->ApplyCommand ("/vis/scene/notifyHandlers");
2676       if (fVerbosity >= warnings) {
2677   G4warn <<
2678     "WARNING: G4VisManager: the scene was empty, \"world\" has been"
2679     "\n  added and the scene handlers notified.";
2680   G4warn << G4endl;
2681       }
2682     }
2683   }
2684   return isValid;
2685 }
2686 
2687 void
2688 G4VisManager::RegisterModelFactories() 
2689 {
2690   if (fVerbosity >= warnings) {
2691     G4warn<<"G4VisManager: No model factories registered with G4VisManager."<<G4endl;
2692     G4warn<<"G4VisManager::RegisterModelFactories() should be overridden in derived"<<G4endl;
2693     G4warn<<"class. See G4VisExecutive for an example."<<G4endl;
2694   }
2695 }
2696 
2697 void G4VisManager::SetUpForAThread()
2698 {
2699   // TODO: protect with G4Threading::IsMultithreadedApplication() instead?
2700 #ifdef G4MULTITHREADED
2701   new G4VisStateDependent(this); 
2702 #endif
2703 }
2704 
2705 void G4VisManager::IgnoreStateChanges(G4bool val)
2706 {
2707   fIgnoreStateChanges = val; 
2708 }
2709