Geant4  10.03
G4VSceneHandler.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4VSceneHandler.cc 101714 2016-11-22 08:53:13Z gcosmo $
28 //
29 //
30 // John Allison 19th July 1996
31 // Abstract interface class for graphics scenes.
32 
33 #include "G4VSceneHandler.hh"
34 
35 #include "G4ios.hh"
36 #include <sstream>
37 
38 #include "G4VisManager.hh"
39 #include "G4VGraphicsSystem.hh"
40 #include "G4VViewer.hh"
41 #include "G4VSolid.hh"
42 #include "G4RotationMatrix.hh"
43 #include "G4ThreeVector.hh"
44 #include "G4VPhysicalVolume.hh"
45 #include "G4Material.hh"
46 #include "G4Polyline.hh"
47 #include "G4Scale.hh"
48 #include "G4Text.hh"
49 #include "G4Circle.hh"
50 #include "G4Square.hh"
51 #include "G4Polymarker.hh"
52 #include "G4Polyhedron.hh"
53 #include "G4Visible.hh"
54 #include "G4VisAttributes.hh"
55 #include "G4VModel.hh"
56 #include "G4TrajectoriesModel.hh"
57 #include "G4Box.hh"
58 #include "G4Cons.hh"
59 #include "G4Orb.hh"
60 #include "G4Para.hh"
61 #include "G4Sphere.hh"
62 #include "G4Torus.hh"
63 #include "G4Trap.hh"
64 #include "G4Trd.hh"
65 #include "G4Tubs.hh"
66 #include "G4Ellipsoid.hh"
67 #include "G4Polycone.hh"
68 #include "G4Polyhedra.hh"
69 #include "G4DisplacedSolid.hh"
70 #include "G4LogicalVolume.hh"
71 #include "G4PhysicalVolumeModel.hh"
72 #include "G4ModelingParameters.hh"
73 #include "G4VTrajectory.hh"
74 #include "G4VTrajectoryPoint.hh"
75 #include "G4HitsModel.hh"
76 #include "G4VHit.hh"
77 #include "G4VDigi.hh"
78 #include "G4ScoringManager.hh"
80 #include "Randomize.hh"
81 #include "G4StateManager.hh"
82 #include "G4RunManager.hh"
83 #ifdef G4MULTITHREADED
84 #include "G4MTRunManager.hh"
85 #endif
86 #include "G4Run.hh"
87 #include "G4Transform3D.hh"
88 #include "G4AttHolder.hh"
89 #include "G4AttDef.hh"
90 #include "G4VVisCommand.hh"
91 #include "G4PhysicalConstants.hh"
92 #include "G4SystemOfUnits.hh"
93 
95  fSystem (system),
96  fSceneHandlerId (id),
97  fViewCount (0),
98  fpViewer (0),
99  fpScene (0),
100  fMarkForClearingTransientStore (true), // Ready for first
101  // ClearTransientStoreIfMarked(),
102  // e.g., at end of run (see
103  // G4VisManager.cc).
104  fReadyForTransients (true), // Only false while processing scene.
105  fProcessingSolid (false),
106  fProcessing2D (false),
107  fpModel (0),
108  fNestingDepth (0),
109  fpVisAttribs (0)
110 {
112  fpScene = pVMan -> GetCurrentScene ();
113  if (name == "") {
114  std::ostringstream ost;
115  ost << fSystem.GetName () << '-' << fSceneHandlerId;
116  fName = ost.str();
117  }
118  else {
119  fName = name;
120  }
123 }
124 
126  G4VViewer* last;
127  while( ! fViewerList.empty() ) {
128  last = fViewerList.back();
129  fViewerList.pop_back();
130  delete last;
131  }
132 }
133 
135 {
136  if (fpScene) {
137  return fpScene->GetExtent();
138  } else {
140  }
141 }
142 
143 void G4VSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
144  const G4VisAttributes& visAttribs) {
145  fObjectTransformation = objectTransformation;
146  fpVisAttribs = &visAttribs;
147  fProcessingSolid = true;
148 }
149 
151  fpVisAttribs = 0;
152  fProcessingSolid = false;
153  if (fReadyForTransients) {
156  }
157 }
158 
160 (const G4Transform3D& objectTransformation) {
161  //static G4int count = 0;
162  //G4cout << "G4VSceneHandler::BeginPrimitives: " << count++ << G4endl;
163  fNestingDepth++;
164  if (fNestingDepth > 1)
166  ("G4VSceneHandler::BeginPrimitives",
167  "visman0101", FatalException,
168  "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
169  fObjectTransformation = objectTransformation;
170 }
171 
173  if (fNestingDepth <= 0)
174  G4Exception("G4VSceneHandler::EndPrimitives",
175  "visman0102", FatalException, "Nesting error.");
176  fNestingDepth--;
177  if (fReadyForTransients) {
180  }
181 }
182 
184 (const G4Transform3D& objectTransformation) {
185  fNestingDepth++;
186  if (fNestingDepth > 1)
188  ("G4VSceneHandler::BeginPrimitives2D",
189  "visman0103", FatalException,
190  "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
191  fObjectTransformation = objectTransformation;
192  fProcessing2D = true;
193 }
194 
196  if (fNestingDepth <= 0)
197  G4Exception("G4VSceneHandler::EndPrimitives2D",
198  "visman0104", FatalException, "Nesting error.");
199  fNestingDepth--;
200  if (fReadyForTransients) {
203  }
204  fProcessing2D = false;
205 }
206 
208 }
209 
211 {
212  fpModel = 0;
213 }
214 
216 
218 
219 template <class T> void G4VSceneHandler::AddSolidT
220 (const T& solid)
221 {
222  // Get and check applicable vis attributes.
223  fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
224  RequestPrimitives (solid);
225 }
226 
227 template <class T> void G4VSceneHandler::AddSolidWithAuxiliaryEdges
228 (const T& solid)
229 {
230  // Get and check applicable vis attributes.
231  fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
232  // Draw with auxiliary edges unless otherwise specified.
233  if (!fpVisAttribs->IsForceAuxEdgeVisible()) {
234  // Create a vis atts object for the modified vis atts.
235  // It is static so that we may return a reliable pointer to it.
236  static G4VisAttributes visAttsWithAuxEdges;
237  // Initialise it with the current vis atts and reset the pointer.
238  visAttsWithAuxEdges = *fpVisAttribs;
239  // Force auxiliary edges visible.
240  visAttsWithAuxEdges.SetForceAuxEdgeVisible();
241  fpVisAttribs = &visAttsWithAuxEdges;
242  }
243  RequestPrimitives (solid);
244 }
245 
246 void G4VSceneHandler::AddSolid (const G4Box& box) {
247  AddSolidT (box);
248  // If your graphics system is sophisticated enough to handle a
249  // particular solid shape as a primitive, in your derived class write a
250  // function to override this.
251  // Your function might look like this...
252  // void G4MySceneHandler::AddSolid (const G4Box& box) {
253  // Get and check applicable vis attributes.
254  // fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
255  // Do not draw if not visible.
256  // if (fpVisAttribs->IsVisible()) {
257  // Get parameters of appropriate object, e.g.:
258  // G4double dx = box.GetXHalfLength ();
259  // G4double dy = box.GetYHalfLength ();
260  // G4double dz = box.GetZHalfLength ();
261  // ...
262  // and Draw or Store in your display List.
263 }
264 
265 void G4VSceneHandler::AddSolid (const G4Cons& cons) {
266  AddSolidT (cons);
267 }
268 
269 void G4VSceneHandler::AddSolid (const G4Orb& orb) {
271 }
272 
273 void G4VSceneHandler::AddSolid (const G4Para& para) {
274  AddSolidT (para);
275 }
276 
277 void G4VSceneHandler::AddSolid (const G4Sphere& sphere) {
279 }
280 
281 void G4VSceneHandler::AddSolid (const G4Torus& torus) {
283 }
284 
285 void G4VSceneHandler::AddSolid (const G4Trap& trap) {
286  AddSolidT (trap);
287 }
288 
289 void G4VSceneHandler::AddSolid (const G4Trd& trd) {
290  AddSolidT (trd);
291 }
292 
293 void G4VSceneHandler::AddSolid (const G4Tubs& tubs) {
294  AddSolidT (tubs);
295 }
296 
297 void G4VSceneHandler::AddSolid (const G4Ellipsoid& ellipsoid) {
298  AddSolidWithAuxiliaryEdges (ellipsoid);
299 }
300 
301 void G4VSceneHandler::AddSolid (const G4Polycone& polycone) {
302  AddSolidT (polycone);
303 }
304 
305 void G4VSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
306  AddSolidT (polyhedra);
307 }
308 
309 void G4VSceneHandler::AddSolid (const G4VSolid& solid) {
310  AddSolidT (solid);
311 }
312 
314  G4TrajectoriesModel* trajectoriesModel =
315  dynamic_cast<G4TrajectoriesModel*>(fpModel);
316  if (trajectoriesModel)
317  traj.DrawTrajectory();
318  else {
320  ("G4VSceneHandler::AddCompound(const G4VTrajectory&)",
321  "visman0105", FatalException, "Not a G4TrajectoriesModel.");
322  }
323 }
324 
326  // Cast away const because Draw is non-const!!!!
327  const_cast<G4VHit&>(hit).Draw();
328 }
329 
331  // Cast away const because Draw is non-const!!!!
332  const_cast<G4VDigi&>(digi).Draw();
333 }
334 
336  //G4cout << "AddCompound: hits: " << &hits << G4endl;
337  G4bool scoreMapHits = false;
339  if (scoringManager) {
340  size_t nMeshes = scoringManager->GetNumberOfMesh();
341  for (size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
342  G4VScoringMesh* mesh = scoringManager->GetMesh(iMesh);
343  if (mesh && mesh->IsActive()) {
344  MeshScoreMap scoreMap = mesh->GetScoreMap();
345  const G4String& mapNam = const_cast<G4THitsMap<G4double>&>(hits).GetName();
346  for(MeshScoreMap::const_iterator i = scoreMap.begin();
347  i != scoreMap.end(); ++i) {
348  const G4String& scoreMapName = i->first;
349  if (scoreMapName == mapNam) {
350  G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
351  scoreMapHits = true;
352  mesh->DrawMesh(scoreMapName, &colorMap);
353  }
354  }
355  }
356  }
357  }
358  if (scoreMapHits) {
359  static G4bool first = true;
360  if (first) {
361  first = false;
362  G4cout <<
363  "Scoring map drawn with default parameters."
364  "\n To get gMocren file for gMocren browser:"
365  "\n /vis/open gMocrenFile"
366  "\n /vis/viewer/flush"
367  "\n Many other options available with /score/draw... commands."
368  "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
369  << G4endl;
370  }
371  } else { // Not score map hits. Just call DrawAllHits.
372  // Cast away const because DrawAllHits is non-const!!!!
373  const_cast<G4THitsMap<G4double>&>(hits).DrawAllHits();
374  }
375 }
376 
378  //G4cout << "AddCompound: hits: " << &hits << G4endl;
379  G4bool scoreMapHits = false;
381  if (scoringManager) {
382  size_t nMeshes = scoringManager->GetNumberOfMesh();
383  for (size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
384  G4VScoringMesh* mesh = scoringManager->GetMesh(iMesh);
385  if (mesh && mesh->IsActive()) {
386  MeshScoreMap scoreMap = mesh->GetScoreMap();
387  for(MeshScoreMap::const_iterator i = scoreMap.begin();
388  i != scoreMap.end(); ++i) {
389  const G4String& scoreMapName = i->first;
390  const G4THitsMap<G4StatDouble>* foundHits = i->second;
391  if (foundHits == &hits) {
392  G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
393  scoreMapHits = true;
394  mesh->DrawMesh(scoreMapName, &colorMap);
395  }
396  }
397  }
398  }
399  }
400  if (scoreMapHits) {
401  static G4bool first = true;
402  if (first) {
403  first = false;
404  G4cout <<
405  "Scoring map drawn with default parameters."
406  "\n To get gMocren file for gMocren browser:"
407  "\n /vis/open gMocrenFile"
408  "\n /vis/viewer/flush"
409  "\n Many other options available with /score/draw... commands."
410  "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
411  << G4endl;
412  }
413  } else { // Not score map hits. Just call DrawAllHits.
414  // Cast away const because DrawAllHits is non-const!!!!
415  const_cast<G4THitsMap<G4StatDouble>&>(hits).DrawAllHits();
416  }
417 }
418 
420  fViewerList.push_back (pViewer);
421 }
422 
424 
425  const G4double margin(0.01);
426  // Fractional margin - ensures scale is comfortably inside viewing
427  // volume.
428  const G4double oneMinusMargin (1. - margin);
429 
430  const G4VisExtent& sceneExtent = fpScene->GetExtent();
431 
432  // Useful constants...
433  const G4double length(scale.GetLength());
434  const G4double halfLength(length / 2.);
435  const G4double tickLength(length / 20.);
436  const G4double piBy2(halfpi);
437 
438  // Get size of scene...
439  const G4double xmin = sceneExtent.GetXmin();
440  const G4double xmax = sceneExtent.GetXmax();
441  const G4double ymin = sceneExtent.GetYmin();
442  const G4double ymax = sceneExtent.GetYmax();
443  const G4double zmin = sceneExtent.GetZmin();
444  const G4double zmax = sceneExtent.GetZmax();
445 
446  // Create (empty) polylines having the same vis attributes...
447  G4Polyline scaleLine, tick11, tick12, tick21, tick22;
448  G4VisAttributes visAtts(*scale.GetVisAttributes()); // Long enough life.
449  scaleLine.SetVisAttributes(&visAtts);
450  tick11.SetVisAttributes(&visAtts);
451  tick12.SetVisAttributes(&visAtts);
452  tick21.SetVisAttributes(&visAtts);
453  tick22.SetVisAttributes(&visAtts);
454 
455  // Add points to the polylines to represent an scale parallel to the
456  // x-axis centred on the origin...
457  G4Point3D r1(G4Point3D(-halfLength, 0., 0.));
458  G4Point3D r2(G4Point3D( halfLength, 0., 0.));
459  scaleLine.push_back(r1);
460  scaleLine.push_back(r2);
461  G4Point3D ticky(0., tickLength, 0.);
462  G4Point3D tickz(0., 0., tickLength);
463  tick11.push_back(r1 + ticky);
464  tick11.push_back(r1 - ticky);
465  tick12.push_back(r1 + tickz);
466  tick12.push_back(r1 - tickz);
467  tick21.push_back(r2 + ticky);
468  tick21.push_back(r2 - ticky);
469  tick22.push_back(r2 + tickz);
470  tick22.push_back(r2 - tickz);
471  G4Point3D textPosition(0., tickLength, 0.);
472 
473  // Transform appropriately...
474 
475  G4Transform3D transformation;
476  if (scale.GetAutoPlacing()) {
477  G4Transform3D rotation;
478  switch (scale.GetDirection()) {
479  case G4Scale::x:
480  break;
481  case G4Scale::y:
482  rotation = G4RotateZ3D(piBy2);
483  break;
484  case G4Scale::z:
485  rotation = G4RotateY3D(piBy2);
486  break;
487  }
488  G4double sxmid;
489  G4double symid;
490  G4double szmid;
491  sxmid = xmin + oneMinusMargin * (xmax - xmin);
492  symid = ymin + margin * (ymax - ymin);
493  szmid = zmin + oneMinusMargin * (zmax - zmin);
494  switch (scale.GetDirection()) {
495  case G4Scale::x:
496  sxmid -= halfLength;
497  break;
498  case G4Scale::y:
499  symid += halfLength;
500  break;
501  case G4Scale::z:
502  szmid -= halfLength;
503  break;
504  }
505  G4Translate3D translation(sxmid, symid, szmid);
506  transformation = translation * rotation;
507  } else {
508  if (fpModel) transformation = fpModel->GetTransformation();
509  }
510 
511  // Draw...
512  // We would like to call BeginPrimitives(transformation) here but
513  // calling BeginPrimitives from within an AddPrimitive is not
514  // allowed! So we have to do our own transformation...
515  AddPrimitive(scaleLine.transform(transformation));
516  AddPrimitive(tick11.transform(transformation));
517  AddPrimitive(tick12.transform(transformation));
518  AddPrimitive(tick21.transform(transformation));
519  AddPrimitive(tick22.transform(transformation));
520  G4Text text(scale.GetAnnotation(),textPosition.transform(transformation));
522  text.SetVisAttributes(va);
523  text.SetScreenSize(scale.GetAnnotationSize());
524  AddPrimitive(text);
525 }
526 
527 void G4VSceneHandler::AddPrimitive (const G4Polymarker& polymarker) {
528  switch (polymarker.GetMarkerType()) {
529  default:
530  case G4Polymarker::dots:
531  {
532  for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
533  G4Circle dot (polymarker);
534  dot.SetPosition (polymarker[iPoint]);
535  dot.SetWorldSize (0.);
536  dot.SetScreenSize (0.1); // Very small circle.
537  AddPrimitive (dot);
538  }
539  }
540  break;
542  {
543  for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
544  G4Circle circle (polymarker);
545  circle.SetPosition (polymarker[iPoint]);
546  AddPrimitive (circle);
547  }
548  }
549  break;
551  {
552  for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
553  G4Square square (polymarker);
554  square.SetPosition (polymarker[iPoint]);
555  AddPrimitive (square);
556  }
557  }
558  break;
559  }
560 }
561 
563  fViewerList.remove(pViewer);
564 }
565 
567  fpScene = pScene;
568  // Notify all viewers that a kernel visit is required.
570  for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
571  (*i) -> SetNeedKernelVisit (true);
572  }
573 }
574 
576  G4Polyhedron::SetNumberOfRotationSteps (GetNoOfSides (fpVisAttribs));
577  G4Polyhedron* pPolyhedron = solid.GetPolyhedron ();
578  G4Polyhedron::ResetNumberOfRotationSteps ();
579  if (pPolyhedron) {
580  pPolyhedron -> SetVisAttributes (fpVisAttribs);
582  AddPrimitive (*pPolyhedron);
583  EndPrimitives ();
584  }
585  else {
587  if (verbosity >= G4VisManager::errors) {
588  G4cerr <<
589  "ERROR: G4VSceneHandler::RequestPrimitives"
590  "\n Polyhedron not available for " << solid.GetName () <<
591  "\n Touchable path: ";
592  G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
593  if (pPVModel) {
594  G4cerr << pPVModel->GetFullPVPath();
595  }
596  G4cerr <<
597  "\n This means it cannot be visualized on most systems (try RayTracer)."
598  "\n 1) The solid may not have implemented the CreatePolyhedron method."
599  "\n 2) For Boolean solids, the BooleanProcessor, which attempts to create"
600  "\n the resultant polyhedron, may have failed."
601  << G4endl;
602  }
603  }
604 }
605 
607 
608  // Assumes graphics database store has already been cleared if
609  // relevant for the particular scene handler.
610 
611  if (!fpScene) return;
612 
613  G4VisManager* visManager = G4VisManager::GetInstance();
614 
615  if (!visManager->GetConcreteInstance()) return;
616 
617  G4VisManager::Verbosity verbosity = visManager->GetVerbosity();
618 
619  fReadyForTransients = false;
620 
621  // Reset fMarkForClearingTransientStore. (Leaving
622  // fMarkForClearingTransientStore true causes problems with
623  // recomputing transients below.) Restore it again at end...
624  G4bool tmpMarkForClearingTransientStore = fMarkForClearingTransientStore;
626 
627  // Traverse geometry tree and send drawing primitives to window(s).
628 
629  const std::vector<G4Scene::Model>& runDurationModelList =
630  fpScene -> GetRunDurationModelList ();
631 
632  if (runDurationModelList.size ()) {
633  if (verbosity >= G4VisManager::confirmations) {
634  G4cout << "Traversing scene data..." << G4endl;
635  }
636 
637  BeginModeling ();
638 
639  // Create modeling parameters from view parameters...
641 
642  for (size_t i = 0; i < runDurationModelList.size (); i++) {
643  if (runDurationModelList[i].fActive) {
644  G4VModel* pModel = runDurationModelList[i].fpModel;
645  // Note: this is not the place to take action on
646  // pModel->GetTransformation(). The model must take care of
647  // this in pModel->DescribeYourselfTo(*this). See, for example,
648  // G4PhysicalVolumeModel and /vis/scene/add/logo.
649  pModel -> SetModelingParameters (pMP);
650  SetModel (pModel); // Store for use by derived class.
651  pModel -> DescribeYourselfTo (*this);
652  pModel -> SetModelingParameters (0);
653  }
654  }
655 
656  delete pMP;
657  EndModeling ();
658  }
659 
660  fReadyForTransients = true;
661 
662  // Refresh event from end-of-event model list.
663  // Allow only in Idle or GeomClosed state...
665  G4ApplicationState state = stateManager->GetCurrentState();
666  if (state == G4State_Idle || state == G4State_GeomClosed) {
667 
668  visManager->SetEventRefreshing(true);
669 
670  if (visManager->GetRequestedEvent()) {
671  DrawEvent(visManager->GetRequestedEvent());
672 
673  } else {
674 
676 #ifdef G4MULTITHREADED
678  { runManager = G4MTRunManager::GetMasterRunManager(); }
679 #endif
680  if (runManager) {
681  const G4Run* run = runManager->GetCurrentRun();
682  const std::vector<const G4Event*>* events =
683  run? run->GetEventVector(): 0;
684  size_t nKeptEvents = 0;
685  if (events) nKeptEvents = events->size();
686  if (nKeptEvents) {
687 
689 
690  if (verbosity >= G4VisManager::confirmations) {
691  G4cout << "Refreshing event..." << G4endl;
692  }
693  const G4Event* event = 0;
694  if (events && events->size()) event = events->back();
695  if (event) DrawEvent(event);
696 
697  } else { // Accumulating events.
698 
699  if (verbosity >= G4VisManager::confirmations) {
700  G4cout << "Refreshing events in run..." << G4endl;
701  }
702  for (const auto& event: *events) {
703  if (event) DrawEvent(event);
704  }
705 
706  if (!fpScene->GetRefreshAtEndOfRun()) {
707  if (verbosity >= G4VisManager::warnings) {
708  G4cout <<
709  "WARNING: Cannot refresh events accumulated over more"
710  "\n than one runs. Refreshed just the last run."
711  << G4endl;
712  }
713  }
714  }
715  }
716  }
717  }
718  visManager->SetEventRefreshing(false);
719  }
720 
721  // Refresh end-of-run model list.
722  // Allow only in Idle or GeomClosed state...
723  if (state == G4State_Idle || state == G4State_GeomClosed) {
725  }
726 
727  fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
728 }
729 
731 {
732  const std::vector<G4Scene::Model>& EOEModelList =
733  fpScene -> GetEndOfEventModelList ();
734  size_t nModels = EOEModelList.size();
735  if (nModels) {
737  pMP->SetEvent(event);
738  for (size_t i = 0; i < nModels; i++) {
739  if (EOEModelList[i].fActive) {
740  G4VModel* pModel = EOEModelList[i].fpModel;
741  pModel -> SetModelingParameters(pMP);
742  SetModel (pModel);
743  pModel -> DescribeYourselfTo (*this);
744  pModel -> SetModelingParameters(0);
745  }
746  }
747  delete pMP;
748  SetModel (0);
749  }
750 }
751 
753 {
754  const std::vector<G4Scene::Model>& EORModelList =
755  fpScene -> GetEndOfRunModelList ();
756  size_t nModels = EORModelList.size();
757  if (nModels) {
759  pMP->SetEvent(0);
760  for (size_t i = 0; i < nModels; i++) {
761  if (EORModelList[i].fActive) {
762  G4VModel* pModel = EORModelList[i].fpModel;
763  pModel -> SetModelingParameters(pMP);
764  SetModel (pModel);
765  pModel -> DescribeYourselfTo (*this);
766  pModel -> SetModelingParameters(0);
767  }
768  }
769  delete pMP;
770  SetModel (0);
771  }
772 }
773 
775 {
776  // Create modeling parameters from View Parameters...
777  if (!fpViewer) return NULL;
778 
779  const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
780 
781  // Convert drawing styles...
782  G4ModelingParameters::DrawingStyle modelDrawingStyle =
784  switch (vp.GetDrawingStyle ()) {
785  default:
787  modelDrawingStyle = G4ModelingParameters::wf;
788  break;
790  modelDrawingStyle = G4ModelingParameters::hlr;
791  break;
793  modelDrawingStyle = G4ModelingParameters::hsr;
794  break;
796  modelDrawingStyle = G4ModelingParameters::hlhsr;
797  break;
798  }
799 
800  // Decide if covered daughters are really to be culled...
801  G4bool reallyCullCovered =
802  vp.IsCullingCovered() // Culling daughters depends also on...
803  && !vp.IsSection () // Sections (DCUT) not requested.
804  && !vp.IsCutaway () // Cutaways not requested.
805  ;
806 
807  G4ModelingParameters* pModelingParams = new G4ModelingParameters
809  modelDrawingStyle,
810  vp.IsCulling (),
811  vp.IsCullingInvisible (),
812  vp.IsDensityCulling (),
813  vp.GetVisibleDensity (),
814  reallyCullCovered,
815  vp.GetNoOfSides ()
816  );
817 
818  pModelingParams->SetWarning
820 
821  pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
822  pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
823 
824  pModelingParams->SetSectionSolid(CreateSectionSolid());
825  pModelingParams->SetCutawaySolid(CreateCutawaySolid());
826  // The polyhedron objects are deleted in the modeling parameters destructor.
827 
829 
830  return pModelingParams;
831 }
832 
834 {
835  G4VSolid* sectioner = 0;
837  if (vp.IsSection () ) {
839  G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
840  G4VSolid* sectionBox =
841  new G4Box("_sectioner", safe, safe, 1.e-5 * radius); // Thin in z-plane.
842  const G4Plane3D& sp = vp.GetSectionPlane ();
843  G4double a = sp.a();
844  G4double b = sp.b();
845  G4double c = sp.c();
846  G4double d = sp.d();
847  G4Transform3D transform = G4TranslateZ3D(-d);
848  const G4Normal3D normal(a,b,c);
849  if (normal != G4Normal3D(0,0,1)) {
850  const G4double angle = std::acos(normal.dot(G4Normal3D(0,0,1)));
851  const G4Vector3D axis = G4Normal3D(0,0,1).cross(normal);
852  transform = G4Rotate3D(angle, axis) * transform;
853  }
854  sectioner = new G4DisplacedSolid
855  ("_displaced_sectioning_box", sectionBox, transform);
856  }
857  return sectioner;
858 }
859 
861 {
862  // To be reviewed.
863  return 0;
864  /*** An alternative way of getting a cutaway is to use
865  Command /vis/scene/add/volume
866  Guidance :
867  Adds a physical volume to current scene, with optional clipping volume.
868  If physical-volume-name is "world" (the default), the top of the
869  main geometry tree (material world) is added. If "worlds", the
870  top of all worlds - material world and parallel worlds, if any - are
871  added. Otherwise a search of all worlds is made, taking the first
872  matching occurence only. To see a representation of the geometry
873  hierarchy of the worlds, try "/vis/drawTree [worlds]" or one of the
874  driver/browser combinations that have the required functionality, e.g., HepRep.
875  If clip-volume-type is specified, the subsequent parameters are used to
876  to define a clipping volume. For example,
877  "/vis/scene/add/volume ! ! ! -box km 0 1 0 1 0 1" will draw the world
878  with the positive octant cut away. (If the Boolean Processor issues
879  warnings try replacing 0 by 0.000000001 or something.)
880  If clip-volume-type is prepended with '-', the clip-volume is subtracted
881  (cutaway). (This is the default if there is no prepended character.)
882  If '*' is prepended, the intersection of the physical-volume and the
883  clip-volume is made. (You can make a section/DCUT with a thin box, for
884  example).
885  For "box", the parameters are xmin,xmax,ymin,ymax,zmin,zmax.
886  Only "box" is programmed at present.
887  ***/
888 }
889 
890 void G4VSceneHandler::LoadAtts(const G4Visible& visible, G4AttHolder* holder)
891 {
892  // Load G4Atts from G4VisAttributes, if any...
893  const G4VisAttributes* va = visible.GetVisAttributes();
894  if (va) {
895  const std::map<G4String,G4AttDef>* vaDefs =
896  va->GetAttDefs();
897  if (vaDefs) {
898  holder->AddAtts(visible.GetVisAttributes()->CreateAttValues(), vaDefs);
899  }
900  }
901 
902  G4PhysicalVolumeModel* pPVModel =
903  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
904  if (pPVModel) {
905  // Load G4Atts from G4PhysicalVolumeModel...
906  const std::map<G4String,G4AttDef>* pvDefs = pPVModel->GetAttDefs();
907  if (pvDefs) {
908  holder->AddAtts(pPVModel->CreateCurrentAttValues(), pvDefs);
909  }
910  }
911 
912  G4TrajectoriesModel* trajModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
913  if (trajModel) {
914  // Load G4Atts from trajectory model...
915  const std::map<G4String,G4AttDef>* trajModelDefs = trajModel->GetAttDefs();
916  if (trajModelDefs) {
917  holder->AddAtts(trajModel->CreateCurrentAttValues(), trajModelDefs);
918  }
919  // Load G4Atts from trajectory...
920  const G4VTrajectory* traj = trajModel->GetCurrentTrajectory();
921  if (traj) {
922  const std::map<G4String,G4AttDef>* trajDefs = traj->GetAttDefs();
923  if (trajDefs) {
924  holder->AddAtts(traj->CreateAttValues(), trajDefs);
925  }
926  G4int nPoints = traj->GetPointEntries();
927  for (G4int i = 0; i < nPoints; ++i) {
928  G4VTrajectoryPoint* trajPoint = traj->GetPoint(i);
929  if (trajPoint) {
930  const std::map<G4String,G4AttDef>* pointDefs = trajPoint->GetAttDefs();
931  if (pointDefs) {
932  holder->AddAtts(trajPoint->CreateAttValues(), pointDefs);
933  }
934  }
935  }
936  }
937  }
938 
939  G4HitsModel* hitsModel = dynamic_cast<G4HitsModel*>(fpModel);
940  if (hitsModel) {
941  // Load G4Atts from hit...
942  const G4VHit* hit = hitsModel->GetCurrentHit();
943  const std::map<G4String,G4AttDef>* hitsDefs = hit->GetAttDefs();
944  if (hitsDefs) {
945  holder->AddAtts(hit->CreateAttValues(), hitsDefs);
946  }
947  }
948 }
949 
951  const G4VisAttributes* pVA = text.GetVisAttributes ();
952  if (!pVA) {
954  }
955  const G4Colour& colour = pVA -> GetColour ();
956  return colour;
957 }
958 
960 {
961  G4double lineWidth = pVisAttribs->GetLineWidth();
962  if (lineWidth < 1.) lineWidth = 1.;
963  lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
964  if (lineWidth < 1.) lineWidth = 1.;
965  return lineWidth;
966 }
967 
969 (const G4VisAttributes* pVisAttribs) {
970  // Drawing style is normally determined by the view parameters, but
971  // it can be overriddden by the ForceDrawingStyle flag in the vis
972  // attributes.
974  fpViewer->GetViewParameters().GetDrawingStyle();
975  if (pVisAttribs -> IsForceDrawingStyle ()) {
977  pVisAttribs -> GetForcedDrawingStyle ();
978  // This is complicated because if hidden line and surface removal
979  // has been requested we wish to preserve this sometimes.
980  switch (forcedStyle) {
981  case (G4VisAttributes::solid):
982  switch (style) {
983  case (G4ViewParameters::hlr):
984  style = G4ViewParameters::hlhsr;
985  break;
987  style = G4ViewParameters::hsr;
988  break;
990  case (G4ViewParameters::hsr):
991  default:
992  break;
993  }
994  break;
996  default:
997  // But if forced style is wireframe, do it, because one of its
998  // main uses is in displaying the consituent solids of a Boolean
999  // solid and their surfaces overlap with the resulting Booean
1000  // solid, making a mess if hlr is specified.
1002  break;
1003  }
1004  }
1005  return style;
1006 }
1007 
1009  G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
1010  if (pVisAttribs -> IsForceAuxEdgeVisible()) {
1011  isAuxEdgeVisible = pVisAttribs->IsForcedAuxEdgeVisible();
1012  }
1013  return isAuxEdgeVisible;
1014 }
1015 
1017 (const G4VMarker& marker,
1018  G4VSceneHandler::MarkerSizeType& markerSizeType)
1019 {
1020  G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
1021  const G4VMarker& defaultMarker =
1022  fpViewer -> GetViewParameters().GetDefaultMarker();
1023  G4double size = userSpecified ?
1024  marker.GetWorldSize() : defaultMarker.GetWorldSize();
1025  if (size) {
1026  // Draw in world coordinates.
1027  markerSizeType = world;
1028  }
1029  else {
1030  size = userSpecified ?
1031  marker.GetScreenSize() : defaultMarker.GetScreenSize();
1032  // Draw in screen coordinates.
1033  markerSizeType = screen;
1034  }
1035  size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
1036  if (markerSizeType == screen && size < 1.) size = 1.;
1037  return size;
1038 }
1039 
1041 {
1042  // No. of sides (lines segments per circle) is normally determined
1043  // by the view parameters, but it can be overriddden by the
1044  // ForceLineSegmentsPerCircle in the vis attributes.
1045  G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
1046  if (pVisAttribs) {
1047  if (pVisAttribs->IsForceLineSegmentsPerCircle())
1048  lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
1049  if (lineSegmentsPerCircle < pVisAttribs->GetMinLineSegmentsPerCircle()) {
1050  lineSegmentsPerCircle = pVisAttribs->GetMinLineSegmentsPerCircle();
1051  G4cout <<
1052  "G4VSceneHandler::GetNoOfSides: attempt to set the"
1053  "\nnumber of line segements per circle < " << lineSegmentsPerCircle
1054  << "; forced to " << pVisAttribs->GetMinLineSegmentsPerCircle() << G4endl;
1055  }
1056  }
1057  return lineSegmentsPerCircle;
1058 }
1059 
1060 std::ostream& operator << (std::ostream& os, const G4VSceneHandler& sh) {
1061 
1062  os << "Scene handler " << sh.fName << " has "
1063  << sh.fViewerList.size () << " viewer(s):";
1064  for (size_t i = 0; i < sh.fViewerList.size (); i++) {
1065  os << "\n " << *(sh.fViewerList [i]);
1066  }
1067 
1068  if (sh.fpScene) {
1069  os << "\n " << *sh.fpScene;
1070  }
1071  else {
1072  os << "\n This scene handler currently has no scene.";
1073  }
1074 
1075  return os;
1076 }
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:665
G4bool GetTransientsDrawnThisRun() const
Direction GetDirection() const
G4String GetName() const
void SetWorldSize(G4double)
virtual void ClearStore()
virtual ~G4VSceneHandler()
void AddAtts(const std::vector< G4AttValue > *values, const std::map< G4String, G4AttDef > *defs)
Definition: G4AttHolder.hh:65
Definition: G4Para.hh:77
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
virtual void AddSolid(const G4Box &)
Definition: G4Text.hh:73
static const G4VisExtent NullExtent
Definition: G4VisExtent.hh:80
virtual G4VSolid * CreateSectionSolid()
MarkerType GetMarkerType() const
G4bool GetAutoPlacing() const
HepGeom::RotateY3D G4RotateY3D
G4double GetXmin() const
Definition: G4VisExtent.hh:89
const std::map< G4String, G4AttDef > * GetAttDefs() const
G4double GetVisibleDensity() const
virtual void BeginModeling()
G4ModelingParameters * CreateModelingParameters()
virtual void BeginPrimitives(const G4Transform3D &objectTransformation)
static G4VVisManager * GetConcreteInstance()
G4double GetWorldSize() const
G4double GetLineWidth() const
const G4Point3D & GetExplodeCentre() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
G4int first(char) const
std::vector< G4AttValue > * CreateCurrentAttValues() const
void SetEventRefreshing(G4bool)
void LoadAtts(const G4Visible &, G4AttHolder *)
const std::map< G4String, G4AttDef > * GetAttDefs() const
G4bool IsCullingInvisible() const
Definition: G4Box.hh:64
G4VViewer * fpViewer
G4double GetExplodeFactor() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
Definition: G4Tubs.hh:85
void RemoveViewerFromList(G4VViewer *pView)
virtual void PostAddSolid()
const G4ViewParameters & GetViewParameters() const
G4double GetXmax() const
Definition: G4VisExtent.hh:90
const G4Transform3D & GetTransformation() const
G4bool IsDensityCulling() const
static G4double angle[DIM]
const G4String & GetName() const
const std::vector< const G4Event * > * GetEventVector() const
Definition: G4Run.hh:115
G4Transform3D fObjectTransformation
const G4Point3D & GetExtentCentre() const
Definition: G4VisExtent.cc:63
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
Definition: G4VHit.hh:48
const char * name(G4int ptype)
Definition: G4Trd.hh:72
const G4VisAttributes * GetVisAttributes() const
G4int GetNoOfSides(const G4VisAttributes *)
virtual const G4VisExtent & GetExtent() const
G4bool GetRefreshAtEndOfEvent() const
int G4int
Definition: G4Types.hh:78
virtual void AddPrimitive(const G4Polyline &)=0
HepGeom::RotateZ3D G4RotateZ3D
const G4Run * GetCurrentRun() const
virtual G4VSolid * CreateCutawaySolid()
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4double GetScreenSize() const
G4VSceneHandler(G4VGraphicsSystem &system, G4int id, const G4String &name="")
virtual void DrawTrajectory() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
void remove(G4VViewer *)
Definition: G4ViewerList.cc:31
void SetCutawaySolid(G4VSolid *pCutawaySolid)
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
static const G4Colour & GetCurrentTextColour()
static G4StateManager * GetStateManager()
G4double GetLineWidth(const G4VisAttributes *)
virtual int GetPointEntries() const =0
const G4int fSceneHandlerId
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:73
void AddViewerToList(G4VViewer *pView)
G4GLOB_DLL std::ostream G4cout
G4bool IsAuxEdgeVisible() const
static G4ScoringManager * GetScoringManagerIfExist()
const G4VisExtent & GetExtent() const
static G4int GetMinLineSegmentsPerCircle()
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4Colour & GetColour()
G4double GetYmax() const
Definition: G4VisExtent.hh:92
void SetEvent(const G4Event *pEvent)
bool G4bool
Definition: G4Types.hh:79
static G4MTRunManager * GetMasterRunManager()
Definition: G4Cons.hh:83
G4bool GetRefreshAtEndOfRun() const
static G4VisManager * GetInstance()
const G4VHit * GetCurrentHit() const
Definition: G4HitsModel.hh:58
Definition: G4Run.hh:46
G4ViewerList fViewerList
virtual void EndModeling()
virtual void EndPrimitives()
HepGeom::Transform3D G4Transform3D
G4double GetZmax() const
Definition: G4VisExtent.hh:94
G4bool IsCullingCovered() const
G4ApplicationState GetCurrentState() const
G4bool IsForcedAuxEdgeVisible() const
G4bool IsMultithreadedApplication()
Definition: G4Threading.cc:152
std::ostream & operator<<(std::ostream &os, const G4VSceneHandler &sh)
Definition: G4Orb.hh:61
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation)
G4bool GetTransientsDrawnThisEvent() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
virtual void SetScene(G4Scene *)
G4int GetForcedLineSegmentsPerCircle() const
void SetSectionSolid(G4VSolid *pSectionSolid)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
G4bool IsSection() const
void SetVisAttributesModifiers(const std::vector< VisAttributesModifier > &)
G4bool IsCutaway() const
HepGeom::Rotate3D G4Rotate3D
const G4VisAttributes * fpVisAttribs
const G4String & GetName() const
G4bool fMarkForClearingTransientStore
virtual void AddCompound(const G4VTrajectory &)
const G4VTrajectory * GetCurrentTrajectory() const
virtual std::vector< G4AttValue > * CreateAttValues() const
void SetPosition(const G4Point3D &)
std::vector< G4AttValue > * CreateCurrentAttValues() const
const std::vector< G4AttValue > * CreateAttValues() const
void DrawEvent(const G4Event *)
G4double GetAnnotationSize() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
Definition: G4VHit.hh:60
G4int GetNoOfSides() const
virtual void ProcessScene()
G4bool fTransientsDrawnThisEvent
const std::map< G4String, G4AttDef > * GetAttDefs() const
G4double GetLength() const
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
void SetWarning(G4bool)
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
G4double GetZmin() const
Definition: G4VisExtent.hh:93
static Verbosity GetVerbosity()
const G4Event * GetRequestedEvent() const
G4bool GetAuxEdgeVisible(const G4VisAttributes *)
HepGeom::Translate3D G4Translate3D
MeshScoreMap GetScoreMap() const
DrawingStyle GetDrawingStyle() const
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4double GetYmin() const
Definition: G4VisExtent.hh:91
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4VHit.hh:67
void SetExplodeFactor(G4double explodeFactor)
HepGeom::Plane3D< G4double > G4Plane3D
Definition: G4Plane3D.hh:37
#define G4endl
Definition: G4ios.hh:61
const G4Plane3D & GetSectionPlane() const
void SetForceAuxEdgeVisible(G4bool=true)
virtual void RequestPrimitives(const G4VSolid &solid)
std::vector< G4VViewer * >::iterator G4ViewerListIterator
Definition: G4ViewerList.hh:43
static constexpr double halfpi
Definition: G4SIunits.hh:77
void AddSolidWithAuxiliaryEdges(const T &solid)
void SetModel(G4VModel *)
G4bool IsActive() const
size_t GetNumberOfMesh() const
const G4VisAttributes * GetDefaultVisAttributes() const
double G4double
Definition: G4Types.hh:76
virtual void ClearTransientStore()
G4bool fTransientsDrawnThisRun
void SetExplodeCentre(const G4Point3D &explodeCentre)
G4VGraphicsSystem & fSystem
G4bool IsForceLineSegmentsPerCircle() const
const G4String & GetAnnotation() const
G4VScoringMesh * GetMesh(G4int i) const
G4bool IsCulling() const
G4Polyline & transform(const G4Transform3D &)
Definition: G4Polyline.cc:38
std::map< G4String, RunScore * > MeshScoreMap
G4ApplicationState
virtual void EndPrimitives2D()
G4GLOB_DLL std::ostream G4cerr
HepGeom::TranslateZ3D G4TranslateZ3D
void SetScreenSize(G4double)
void AddSolidT(const T &solid)
HepGeom::Normal3D< G4double > G4Normal3D
Definition: G4Normal3D.hh:35
const G4Colour & GetTextColour(const G4Text &)