source: trunk/source/visualization/management/src/G4VSceneHandler.cc

Last change on this file was 1346, checked in by garnier, 14 years ago

before tag

  • Property svn:mime-type set to text/cpp
File size: 30.4 KB
Line 
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,v 1.96 2010/11/05 16:00:11 allison Exp $
28// GEANT4 tag $Name:  $
29//
30//
31// John Allison  19th July 1996
32// Abstract interface class for graphics scenes.
33
34#include "G4VSceneHandler.hh"
35
36#include "G4ios.hh"
37#include <sstream>
38
39#include "G4VisManager.hh"
40#include "G4VGraphicsSystem.hh"
41#include "G4VViewer.hh"
42#include "G4VSolid.hh"
43#include "G4RotationMatrix.hh"
44#include "G4ThreeVector.hh"
45#include "G4VPhysicalVolume.hh"
46#include "G4Material.hh"
47#include "G4Polyline.hh"
48#include "G4Scale.hh"
49#include "G4Text.hh"
50#include "G4Circle.hh"
51#include "G4Square.hh"
52#include "G4Polymarker.hh"
53#include "G4Polyhedron.hh"
54#include "G4NURBS.hh"
55#include "G4Visible.hh"
56#include "G4VisAttributes.hh"
57#include "G4VModel.hh"
58#include "G4TrajectoriesModel.hh"
59#include "G4Box.hh"
60#include "G4Cons.hh"
61#include "G4Tubs.hh"
62#include "G4Trd.hh"
63#include "G4Trap.hh"
64#include "G4Sphere.hh"
65#include "G4Para.hh"
66#include "G4Torus.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"
79#include "G4DefaultLinearColorMap.hh"
80#include "Randomize.hh"
81#include "G4StateManager.hh"
82#include "G4RunManager.hh"
83#include "G4Run.hh"
84#include "G4Transform3D.hh"
85#include "G4AttHolder.hh"
86#include "G4AttDef.hh"
87
88G4VSceneHandler::G4VSceneHandler (G4VGraphicsSystem& system, G4int id, const G4String& name):
89  fSystem                (system),
90  fSceneHandlerId        (id),
91  fViewCount             (0),
92  fpViewer               (0),
93  fpScene                (0),
94  fMarkForClearingTransientStore (true),  // Ready for first
95                                          // ClearTransientStoreIfMarked(),
96                                          // e.g., at end of run (see
97                                          // G4VisManager.cc).
98  fReadyForTransients    (true),  // Only false while processing scene.
99  fProcessingSolid       (false),
100  fSecondPassRequested   (false),
101  fSecondPass            (false),
102  fpModel                (0),
103  fpObjectTransformation (0),
104  fNestingDepth          (0),
105  fpVisAttribs           (0)
106{
107  G4VisManager* pVMan = G4VisManager::GetInstance ();
108  fpScene = pVMan -> GetCurrentScene ();
109  if (name == "") {
110    std::ostringstream ost;
111    ost << fSystem.GetName () << '-' << fSceneHandlerId;
112    fName = ost.str();
113  }
114  else {
115    fName = name;
116  }
117  fTransientsDrawnThisEvent = pVMan->GetTransientsDrawnThisEvent();
118  fTransientsDrawnThisRun = pVMan->GetTransientsDrawnThisRun();
119}
120
121G4VSceneHandler::~G4VSceneHandler () {
122  G4VViewer* last;
123  while( ! fViewerList.empty() ) {
124    last = fViewerList.back();
125    fViewerList.pop_back();
126    delete last;
127  }
128}
129
130void G4VSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
131                                  const G4VisAttributes& visAttribs) {
132  fpObjectTransformation = &objectTransformation;
133  fpVisAttribs = &visAttribs;
134  fProcessingSolid = true;
135}
136
137void G4VSceneHandler::PostAddSolid () {
138  fpObjectTransformation = 0;
139  fpVisAttribs = 0;
140  fProcessingSolid = false;
141  if (fReadyForTransients) {
142    fTransientsDrawnThisEvent = true;
143    fTransientsDrawnThisRun = true;
144  }
145}
146
147void G4VSceneHandler::BeginPrimitives
148(const G4Transform3D& objectTransformation) {
149  fNestingDepth++;
150  if (fNestingDepth > 1)
151    G4Exception("G4VSceneHandler::BeginPrimitives: Nesting detected."
152                "\n  It is illegal to nest Begin/EndPrimitives.");
153  fpObjectTransformation = &objectTransformation;
154}
155
156void G4VSceneHandler::EndPrimitives () {
157  if (fNestingDepth <= 0)
158    G4Exception("G4VSceneHandler::EndPrimitives: Nesting error");
159  fNestingDepth--;
160  fpObjectTransformation = 0;
161  if (fReadyForTransients) {
162    fTransientsDrawnThisEvent = true;
163    fTransientsDrawnThisRun = true;
164  }
165}
166
167void G4VSceneHandler::BeginPrimitives2D
168(const G4Transform3D& objectTransformation) {
169  fNestingDepth++;
170  if (fNestingDepth > 1)
171    G4Exception("G4VSceneHandler::BeginPrimitives2D: Nesting detected."
172                "\n  It is illegal to nest Begin/EndPrimitives.");
173  fpObjectTransformation = &objectTransformation;
174}
175
176void G4VSceneHandler::EndPrimitives2D () {
177  if (fNestingDepth <= 0)
178    G4Exception("G4VSceneHandler::EndPrimitives2D: Nesting error");
179  fNestingDepth--;
180  fpObjectTransformation = 0;
181  if (fReadyForTransients) {
182    fTransientsDrawnThisEvent = true;
183    fTransientsDrawnThisRun = true;
184  }
185}
186
187void G4VSceneHandler::BeginModeling () {
188}
189
190void G4VSceneHandler::EndModeling ()
191{
192  fpModel = 0;
193}
194
195void G4VSceneHandler::ClearStore () {
196  // if (fpViewer) fpViewer -> NeedKernelVisit (true);
197  // ?? Viewer is supposed to be smart enough to know when to visit
198  // kernel, but a problem in OpenGL Stored seems to require a forced
199  // kernel visit triggered by the above code.  John Allison Aug 2001
200  // Feb 2005 - commented out.  Let's fix OpenGL if necessary.
201}
202
203void G4VSceneHandler::ClearTransientStore () {
204}
205
206void G4VSceneHandler::AddSolid (const G4Box& box) {
207  RequestPrimitives (box);
208// If your graphics system is sophisticated enough to handle a
209//  particular solid shape as a primitive, in your derived class write a
210//  function to override this.  (Note: some compilers warn that your
211//  function "hides" this one.  That's OK.)
212// Your function might look like this...
213// void G4MyScene::AddSolid (const G4Box& box) {
214// Get parameters of appropriate object, e.g.:
215//   G4double dx = box.GetXHalfLength ();
216//   G4double dy = box.GetYHalfLength ();
217//   G4double dz = box.GetZHalfLength ();
218// and Draw or Store in your display List.
219}
220
221void G4VSceneHandler::AddSolid (const G4Tubs& tubs) {
222  RequestPrimitives (tubs);
223}
224
225void G4VSceneHandler::AddSolid (const G4Cons& cons) {
226  RequestPrimitives (cons);
227}
228
229void G4VSceneHandler::AddSolid (const G4Trd& trd) {
230  RequestPrimitives (trd);
231}
232
233void G4VSceneHandler::AddSolid (const G4Trap& trap) {
234  RequestPrimitives (trap);
235}
236
237void G4VSceneHandler::AddSolid (const G4Sphere& sphere) {
238  RequestPrimitives (sphere );
239}
240
241void G4VSceneHandler::AddSolid (const G4Para& para) {
242  RequestPrimitives (para);
243}
244
245void G4VSceneHandler::AddSolid (const G4Torus& torus) {
246  RequestPrimitives (torus);
247}
248
249void G4VSceneHandler::AddSolid (const G4Polycone& polycone) {
250  RequestPrimitives (polycone);
251}
252
253void G4VSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
254  RequestPrimitives (polyhedra);
255}
256
257void G4VSceneHandler::AddSolid (const G4VSolid& solid) {
258  RequestPrimitives (solid);
259}
260
261void G4VSceneHandler::AddCompound (const G4VTrajectory& traj) {
262  G4TrajectoriesModel* trajectoriesModel =
263    dynamic_cast<G4TrajectoriesModel*>(fpModel);
264  if (!trajectoriesModel) G4Exception
265    ("G4VSceneHandler::AddCompound(const G4VTrajectory&): Not a G4TrajectoriesModel.");
266  else {
267    if (trajectoriesModel->IsDrawingModeSet()) {
268      traj.DrawTrajectory(trajectoriesModel->GetDrawingMode());
269    } else {
270      traj.DrawTrajectory();
271    }
272  }
273}
274
275void G4VSceneHandler::AddCompound (const G4VHit& hit) {
276  // Cast away const because Draw is non-const!!!!
277  const_cast<G4VHit&>(hit).Draw();
278}
279
280void G4VSceneHandler::AddCompound (const G4VDigi& digi) {
281  // Cast away const because Draw is non-const!!!!
282  const_cast<G4VDigi&>(digi).Draw();
283}
284
285void G4VSceneHandler::AddCompound (const G4THitsMap<G4double>& hits) {
286  //G4cout << "AddCompound: hits: " << &hits << G4endl;
287  G4bool scoreMapHits = false;
288  G4ScoringManager* scoringManager = G4ScoringManager::GetScoringManagerIfExist();
289  if (scoringManager) {
290    size_t nMeshes = scoringManager->GetNumberOfMesh();
291    for (size_t i = 0; i < nMeshes; ++i) {
292      G4VScoringMesh* mesh = scoringManager->GetMesh(i);
293      if (mesh && mesh->IsActive()) {
294        MeshScoreMap scoreMap = mesh->GetScoreMap();
295        for(MeshScoreMap::const_iterator i = scoreMap.begin();
296            i != scoreMap.end(); ++i) {
297          const G4String& scoreMapName = i->first;
298          const G4THitsMap<G4double>* foundHits = i->second;
299          if (foundHits == &hits) {
300            G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
301            scoreMapHits = true;
302            mesh->DrawMesh(scoreMapName, &colorMap);
303          }
304        }
305      }
306    }
307  }
308  if (scoreMapHits) {
309    static G4bool first = true;
310    if (first) {
311      first = false;
312      G4cout <<
313        "Scoring map drawn with default parameters."
314        "\n  To get gMocren file for gMocren browser:"
315        "\n    /vis/open gMocrenFile"
316        "\n    /vis/viewer/flush"
317        "\n  Many other options available with /score/draw... commands."
318        "\n  You might want to \"/vis/viewer/set/autoRefresh false\"."
319             << G4endl;
320    }
321  } else {  // Not score map hits.  Just call DrawAllHits.
322    // Cast away const because DrawAllHits is non-const!!!!
323    const_cast<G4THitsMap<G4double>&>(hits).DrawAllHits();
324  }
325}
326
327void G4VSceneHandler::AddViewerToList (G4VViewer* pViewer) {
328  fViewerList.push_back (pViewer);
329}
330
331void G4VSceneHandler::AddPrimitive (const G4Scale& scale) {
332
333  const G4double margin(0.01);
334  // Fractional margin - ensures scale is comfortably inside viewing
335  // volume.
336  const G4double oneMinusMargin (1. - margin);
337
338  const G4VisExtent& sceneExtent = fpScene->GetExtent();
339
340  // Useful constants...
341  const G4double length(scale.GetLength());
342  const G4double halfLength(length / 2.);
343  const G4double tickLength(length / 20.);
344  const G4double piBy2(halfpi);
345
346  // Get size of scene...
347  const G4double xmin = sceneExtent.GetXmin();
348  const G4double xmax = sceneExtent.GetXmax();
349  const G4double ymin = sceneExtent.GetYmin();
350  const G4double ymax = sceneExtent.GetYmax();
351  const G4double zmin = sceneExtent.GetZmin();
352  const G4double zmax = sceneExtent.GetZmax();
353
354  // Create (empty) polylines having the same vis attributes...
355  G4Polyline scaleLine, tick11, tick12, tick21, tick22;
356  G4VisAttributes visAtts(*scale.GetVisAttributes());  // Long enough life.
357  scaleLine.SetVisAttributes(&visAtts);
358  tick11.SetVisAttributes(&visAtts);
359  tick12.SetVisAttributes(&visAtts);
360  tick21.SetVisAttributes(&visAtts);
361  tick22.SetVisAttributes(&visAtts);
362
363  // Add points to the polylines to represent an scale parallel to the
364  // x-axis centred on the origin...
365  G4Point3D r1(G4Point3D(-halfLength, 0., 0.));
366  G4Point3D r2(G4Point3D( halfLength, 0., 0.));
367  scaleLine.push_back(r1);
368  scaleLine.push_back(r2);
369  G4Point3D ticky(0., tickLength, 0.);
370  G4Point3D tickz(0., 0., tickLength);
371  tick11.push_back(r1 + ticky);
372  tick11.push_back(r1 - ticky);
373  tick12.push_back(r1 + tickz);
374  tick12.push_back(r1 - tickz);
375  tick21.push_back(r2 + ticky);
376  tick21.push_back(r2 - ticky);
377  tick22.push_back(r2 + tickz);
378  tick22.push_back(r2 - tickz);
379  G4Point3D textPosition(0., tickLength, 0.);
380
381  // Transform appropriately...
382
383  G4Transform3D transformation;
384  if (scale.GetAutoPlacing()) {
385    G4Transform3D rotation;
386    switch (scale.GetDirection()) {
387    case G4Scale::x:
388      break;
389    case G4Scale::y:
390      rotation = G4RotateZ3D(piBy2);
391      break;
392    case G4Scale::z:
393      rotation = G4RotateY3D(piBy2);
394      break;
395    }
396    G4double sxmid(scale.GetXmid());
397    G4double symid(scale.GetYmid());
398    G4double szmid(scale.GetZmid());
399    sxmid = xmin + oneMinusMargin * (xmax - xmin);
400    symid = ymin + margin * (ymax - ymin);
401    szmid = zmin + oneMinusMargin * (zmax - zmin);
402    switch (scale.GetDirection()) {
403    case G4Scale::x:
404      sxmid -= halfLength;
405      break;
406    case G4Scale::y:
407      symid += halfLength;
408      break;
409    case G4Scale::z:
410      szmid -= halfLength;
411      break;
412    }
413    G4Translate3D translation(sxmid, symid, szmid);
414    transformation = translation * rotation;
415  } else {
416    if (fpModel) transformation = fpModel->GetTransformation();
417  }
418
419  // Draw...
420  // We would like to call BeginPrimitives(transformation) here but
421  // calling BeginPrimitives from within an AddPrimitive is not
422  // allowed!  So we have to do our own transformation...
423  AddPrimitive(scaleLine.transform(transformation));
424  AddPrimitive(tick11.transform(transformation));
425  AddPrimitive(tick12.transform(transformation));
426  AddPrimitive(tick21.transform(transformation));
427  AddPrimitive(tick22.transform(transformation));
428  G4Text text(scale.GetAnnotation(),textPosition.transform(transformation));
429  text.SetScreenSize(12.);
430  AddPrimitive(text);
431}
432
433void G4VSceneHandler::AddPrimitive (const G4Polymarker& polymarker) {
434  switch (polymarker.GetMarkerType()) {
435  default:
436  case G4Polymarker::dots:
437    {
438      for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
439        G4Circle dot (polymarker);
440        dot.SetPosition (polymarker[iPoint]);
441        dot.SetWorldSize  (0.);
442        dot.SetScreenSize (0.1);  // Very small circle.
443        AddPrimitive (dot);
444      }
445    }
446    break;
447  case G4Polymarker::circles:
448    {
449      for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
450        G4Circle circle (polymarker);
451        circle.SetPosition (polymarker[iPoint]);
452        AddPrimitive (circle);
453      }
454    }
455    break;
456  case G4Polymarker::squares:
457    {
458      for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
459        G4Square square (polymarker);
460        square.SetPosition (polymarker[iPoint]);
461        AddPrimitive (square);
462      }
463    }
464    break;
465  }
466}
467
468void G4VSceneHandler::RemoveViewerFromList (G4VViewer* pViewer) {
469  fViewerList.remove(pViewer);
470}
471
472void G4VSceneHandler::SetScene (G4Scene* pScene) {
473  fpScene = pScene;
474  // Notify all viewers that a kernel visit is required.
475  G4ViewerListIterator i;
476  for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
477    (*i) -> SetNeedKernelVisit (true);
478  }
479}
480
481void G4VSceneHandler::RequestPrimitives (const G4VSolid& solid) {
482  BeginPrimitives (*fpObjectTransformation);
483  G4NURBS* pNURBS = 0;
484  G4Polyhedron* pPolyhedron = 0;
485  switch (fpViewer -> GetViewParameters () . GetRepStyle ()) {
486  case G4ViewParameters::nurbs:
487    pNURBS = solid.CreateNURBS ();
488    if (pNURBS) {
489      pNURBS -> SetVisAttributes (fpVisAttribs);
490      AddPrimitive (*pNURBS);
491      delete pNURBS;
492      break;
493    }
494    else {
495      G4VisManager::Verbosity verbosity = G4VisManager::GetVerbosity();
496      if (verbosity >= G4VisManager::errors) {
497        G4cout <<
498          "ERROR: G4VSceneHandler::RequestPrimitives"
499          "\n  NURBS not available for "
500               << solid.GetName () << G4endl;
501        G4cout << "Trying polyhedron." << G4endl;
502      }
503    }
504    // Dropping through to polyhedron...
505  case G4ViewParameters::polyhedron:
506  default:
507    G4Polyhedron::SetNumberOfRotationSteps (GetNoOfSides (fpVisAttribs));
508    pPolyhedron = solid.GetPolyhedron ();
509    G4Polyhedron::ResetNumberOfRotationSteps ();
510    if (pPolyhedron) {
511      pPolyhedron -> SetVisAttributes (fpVisAttribs);
512      AddPrimitive (*pPolyhedron);
513    }
514    else {
515      G4VisManager::Verbosity verbosity = G4VisManager::GetVerbosity();
516      if (verbosity >= G4VisManager::errors) {
517        G4cout <<
518          "ERROR: G4VSceneHandler::RequestPrimitives"
519          "\n  Polyhedron not available for " << solid.GetName () <<
520          ".\n  This means it cannot be visualized on most systems."
521          "\n  Contact the Visualization Coordinator." << G4endl;
522      }
523    }
524    break;
525  }
526  EndPrimitives ();
527}
528
529void G4VSceneHandler::ProcessScene (G4VViewer&) {
530
531  if (!fpScene) return;
532
533  G4VisManager* visManager = G4VisManager::GetInstance();
534
535  if (!visManager->GetConcreteInstance()) return;
536
537  G4VisManager::Verbosity verbosity = visManager->GetVerbosity();
538
539  fReadyForTransients = false;
540
541  // Clear stored scene, if any, i.e., display lists, scene graphs.
542  ClearStore ();
543
544  // Reset fMarkForClearingTransientStore.  No need to clear transient
545  // store since it has just been cleared above.  (Leaving
546  // fMarkForClearingTransientStore true causes problems with
547  // recomputing transients below.)  Restore it again at end...
548  G4bool tmpMarkForClearingTransientStore = fMarkForClearingTransientStore;
549  fMarkForClearingTransientStore = false;
550
551  // Traverse geometry tree and send drawing primitives to window(s).
552
553  const std::vector<G4VModel*>& runDurationModelList =
554    fpScene -> GetRunDurationModelList ();
555
556  if (runDurationModelList.size ()) {
557    if (verbosity >= G4VisManager::confirmations) {
558      G4cout << "Traversing scene data..." << G4endl;
559    }
560
561    BeginModeling ();
562
563    // Create modeling parameters from view parameters...
564    G4ModelingParameters* pMP = CreateModelingParameters ();
565
566    for (size_t i = 0; i < runDurationModelList.size (); i++) {
567      G4VModel* pModel = runDurationModelList[i];
568      // Note: this is not the place to take action on
569      // pModel->GetTransformation().  The model must take care of
570      // this in pModel->DescribeYourselfTo(*this).  See, for example,
571      // G4PhysicalVolumeModel and /vis/scene/add/logo.
572      pModel -> SetModelingParameters (pMP);
573      SetModel (pModel);  // Store for use by derived class.
574      pModel -> DescribeYourselfTo (*this);
575      pModel -> SetModelingParameters (0);
576    }
577
578    // Repeat if required...
579    if (fSecondPassRequested) {
580      fSecondPass = true;
581      for (size_t i = 0; i < runDurationModelList.size (); i++) {
582        G4VModel* pModel = runDurationModelList[i];
583        pModel -> SetModelingParameters (pMP);
584        SetModel (pModel);  // Store for use by derived class.
585        pModel -> DescribeYourselfTo (*this);
586        pModel -> SetModelingParameters (0);
587      }
588      fSecondPass = false;
589      fSecondPassRequested = false;
590    }
591
592    delete pMP;
593    EndModeling ();
594
595  }
596
597  fpViewer->FinishView();  // Flush streams and/or swap buffers.
598
599  fReadyForTransients = true;
600
601  // Refresh event from end-of-event model list.
602  // Allow only in Idle or GeomClosed state...
603  G4StateManager* stateManager = G4StateManager::GetStateManager();
604  G4ApplicationState state = stateManager->GetCurrentState();
605  if (state == G4State_Idle || state == G4State_GeomClosed) {
606
607    visManager->SetEventRefreshing(true);
608
609    if (visManager->GetRequestedEvent()) {
610      DrawEvent(visManager->GetRequestedEvent());
611
612    } else {
613
614      G4RunManager* runManager = G4RunManager::GetRunManager();
615      if (runManager) {
616        const G4Run* run = runManager->GetCurrentRun();
617        const std::vector<const G4Event*>* events =
618          run? run->GetEventVector(): 0;
619        size_t nKeptEvents = 0;
620        if (events) nKeptEvents = events->size();
621        if (nKeptEvents) {
622
623          if (fpScene->GetRefreshAtEndOfEvent()) {
624
625            if (verbosity >= G4VisManager::confirmations) {
626              G4cout << "Refreshing event..." << G4endl;
627            }
628            const G4Event* event = 0;
629            if (events && events->size()) event = events->back();
630            if (event) DrawEvent(event);
631
632          } else {  // Accumulating events.
633
634            if (verbosity >= G4VisManager::confirmations) {
635              G4cout << "Refreshing events in run..." << G4endl;
636            }
637            for (size_t i = 0; i < nKeptEvents; ++i) {
638              const G4Event* event = (*events)[i];
639              if (event) DrawEvent(event);
640            }
641
642            if (!fpScene->GetRefreshAtEndOfRun()) {
643              if (verbosity >= G4VisManager::warnings) {
644                G4cout <<
645                  "WARNING: Cannot refresh events accumulated over more"
646                  "\n  than one runs.  Refreshed just the last run."
647                       << G4endl;
648              }
649            }
650          }
651        }
652      }
653    }
654    visManager->SetEventRefreshing(false);
655  }
656
657  // Refresh end-of-run model list.
658  // Allow only in Idle or GeomClosed state...
659  if (state == G4State_Idle || state == G4State_GeomClosed) {
660    DrawEndOfRunModels();
661  }
662
663  fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
664}
665
666void G4VSceneHandler::DrawEvent(const G4Event* event)
667{
668  const std::vector<G4VModel*>& EOEModelList =
669    fpScene -> GetEndOfEventModelList ();
670  size_t nModels = EOEModelList.size();
671  if (nModels) {
672    G4ModelingParameters* pMP = CreateModelingParameters();
673    pMP->SetEvent(event);
674    for (size_t i = 0; i < nModels; i++) {
675      G4VModel* pModel = EOEModelList [i];
676      pModel -> SetModelingParameters(pMP);
677      SetModel (pModel);
678      pModel -> DescribeYourselfTo (*this);
679      pModel -> SetModelingParameters(0);
680    }
681    delete pMP;
682    SetModel (0);
683  }
684}
685
686void G4VSceneHandler::DrawEndOfRunModels()
687{
688  const std::vector<G4VModel*>& EORModelList =
689    fpScene -> GetEndOfRunModelList ();
690  size_t nModels = EORModelList.size();
691  if (nModels) {
692    G4ModelingParameters* pMP = CreateModelingParameters();
693    pMP->SetEvent(0);
694    for (size_t i = 0; i < nModels; i++) {
695      G4VModel* pModel = EORModelList [i];
696      pModel -> SetModelingParameters(pMP);
697      SetModel (pModel);
698      pModel -> DescribeYourselfTo (*this);
699      pModel -> SetModelingParameters(0);
700    }
701    delete pMP;
702    SetModel (0);
703  }
704}
705
706G4ModelingParameters* G4VSceneHandler::CreateModelingParameters ()
707{
708  // Create modeling parameters from View Parameters...
709  const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
710
711  // Convert drawing styles...
712  G4ModelingParameters::DrawingStyle modelDrawingStyle =
713    G4ModelingParameters::wf;
714  switch (vp.GetDrawingStyle ()) {
715  default:
716  case G4ViewParameters::wireframe:
717    modelDrawingStyle = G4ModelingParameters::wf;
718    break;
719  case G4ViewParameters::hlr:
720    modelDrawingStyle = G4ModelingParameters::hlr;
721    break;
722  case G4ViewParameters::hsr:
723    modelDrawingStyle = G4ModelingParameters::hsr;
724    break;
725  case G4ViewParameters::hlhsr:
726    modelDrawingStyle = G4ModelingParameters::hlhsr;
727    break;
728  }
729
730  // Decide if covered daughters are really to be culled...
731  G4bool reallyCullCovered =
732    vp.IsCullingCovered()   // Culling daughters depends also on...
733    && !vp.IsSection ()     // Sections (DCUT) not requested.
734    && !vp.IsCutaway ()     // Cutaways not requested.
735    ;
736
737  G4ModelingParameters* pModelingParams = new G4ModelingParameters
738    (vp.GetDefaultVisAttributes (),
739     modelDrawingStyle,
740     vp.IsCulling (),
741     vp.IsCullingInvisible (),
742     vp.IsDensityCulling (),
743     vp.GetVisibleDensity (),
744     reallyCullCovered,
745     vp.GetNoOfSides ()
746     );
747
748  pModelingParams->SetWarning
749    (G4VisManager::GetVerbosity() >= G4VisManager::warnings);
750
751  pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
752  pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
753
754  pModelingParams->SetSectionSolid(CreateSectionSolid());
755  pModelingParams->SetCutawaySolid(CreateCutawaySolid());
756  // The polyhedron objects are deleted in the modeling parameters destructor.
757
758  return pModelingParams;
759}
760
761G4VSolid* G4VSceneHandler::CreateSectionSolid()
762{
763  G4VSolid* sectioner = 0;
764  const G4ViewParameters& vp = fpViewer->GetViewParameters();
765  if (vp.IsSection () ) {
766    G4double radius = fpScene->GetExtent().GetExtentRadius();
767    G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
768    G4VSolid* sectionBox =
769      new G4Box("_sectioner", safe, safe, 1.e-5 * radius);  // Thin in z-plane.
770    const G4Plane3D& s = vp.GetSectionPlane ();
771    G4double a = s.a();
772    G4double b = s.b();
773    G4double c = s.c();
774    G4double d = s.d();
775    G4Transform3D transform = G4TranslateZ3D(-d);
776    const G4Normal3D normal(a,b,c);
777    if (normal != G4Normal3D(0,0,1)) {
778      const G4double angle = std::acos(normal.dot(G4Normal3D(0,0,1)));
779      const G4Vector3D axis = G4Normal3D(0,0,1).cross(normal);
780      transform = G4Rotate3D(angle, axis) * transform;
781    }
782    sectioner = new G4DisplacedSolid
783      ("_displaced_sectioning_box", sectionBox, transform);
784  }
785  return sectioner;
786}
787
788G4VSolid* G4VSceneHandler::CreateCutawaySolid()
789{
790  return 0;
791}
792
793void G4VSceneHandler::LoadAtts(const G4Visible& visible, G4AttHolder* holder)
794{
795  // Load G4Atts from G4VisAttributes, if any...
796  const G4VisAttributes* va = visible.GetVisAttributes();
797  if (va) {
798    const std::map<G4String,G4AttDef>* vaDefs =
799      va->GetAttDefs();
800    if (vaDefs) {
801      holder->AddAtts(visible.GetVisAttributes()->CreateAttValues(), vaDefs);
802    }
803  }
804
805  G4PhysicalVolumeModel* pPVModel =
806    dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
807  if (pPVModel) {
808    // Load G4Atts from G4PhysicalVolumeModel...
809    const std::map<G4String,G4AttDef>* defs = pPVModel->GetAttDefs();
810    if (defs) {
811      holder->AddAtts(pPVModel->CreateCurrentAttValues(), defs);
812    }
813  }
814
815  G4TrajectoriesModel* trajModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
816  if (trajModel) {
817    // Load G4Atts from trajectory...
818    const G4VTrajectory* traj = trajModel->GetCurrentTrajectory();
819    const std::map<G4String,G4AttDef>* defs = traj->GetAttDefs();
820    if (defs) {
821      holder->AddAtts(traj->CreateAttValues(), defs);
822    }
823    G4int nPoints = traj->GetPointEntries();
824    for (G4int i = 0; i < nPoints; ++i) {
825      G4VTrajectoryPoint* trajPoint = traj->GetPoint(i);
826      const std::map<G4String,G4AttDef>* defs = trajPoint->GetAttDefs();
827      if (defs) {
828        holder->AddAtts(trajPoint->CreateAttValues(), defs);
829      }
830    }
831  }
832
833  G4HitsModel* hitsModel = dynamic_cast<G4HitsModel*>(fpModel);
834  if (hitsModel) {
835    // Load G4Atts from hit...
836    const G4VHit* hit = hitsModel->GetCurrentHit();
837    const std::map<G4String,G4AttDef>* defs = hit->GetAttDefs();
838    if (defs) {
839      holder->AddAtts(hit->CreateAttValues(), defs);
840    }
841  }
842}
843
844const G4Colour& G4VSceneHandler::GetColour (const G4Visible& visible) {
845  // Colour is determined by the applicable vis attributes.
846  const G4Colour& colour = fpViewer ->
847    GetApplicableVisAttributes (visible.GetVisAttributes ()) -> GetColour ();
848  return colour;
849}
850
851const G4Colour& G4VSceneHandler::GetTextColour (const G4Text& text) {
852  const G4VisAttributes* pVA = text.GetVisAttributes ();
853  if (!pVA) {
854    pVA = fpViewer -> GetViewParameters (). GetDefaultTextVisAttributes ();
855  }
856  const G4Colour& colour = pVA -> GetColour ();
857  return colour;
858}
859
860G4double G4VSceneHandler::GetLineWidth(const G4VisAttributes* pVisAttribs)
861{
862  G4double lineWidth = pVisAttribs->GetLineWidth();
863  if (lineWidth < 1.) lineWidth = 1.;
864  lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
865  if (lineWidth < 1.) lineWidth = 1.;
866  return lineWidth;
867}
868
869G4ViewParameters::DrawingStyle G4VSceneHandler::GetDrawingStyle
870(const G4VisAttributes* pVisAttribs) {
871  // Drawing style is normally determined by the view parameters, but
872  // it can be overriddden by the ForceDrawingStyle flag in the vis
873  // attributes.
874  G4ViewParameters::DrawingStyle style =
875    fpViewer->GetViewParameters().GetDrawingStyle();
876  if (pVisAttribs -> IsForceDrawingStyle ()) {
877    G4VisAttributes::ForcedDrawingStyle forcedStyle =
878      pVisAttribs -> GetForcedDrawingStyle ();
879    // This is complicated because if hidden line and surface removal
880    // has been requested we wish to preserve this sometimes.
881    switch (forcedStyle) {
882    case (G4VisAttributes::solid):
883      switch (style) {
884      case (G4ViewParameters::hlr):
885        style = G4ViewParameters::hlhsr;
886        break;
887      case (G4ViewParameters::wireframe):
888        style = G4ViewParameters::hsr;
889        break;
890      case (G4ViewParameters::hlhsr):
891      case (G4ViewParameters::hsr):
892      default:
893        break;
894      }
895      break;
896    case (G4VisAttributes::wireframe):
897    default:
898      // But if forced style is wireframe, do it, because one of its
899      // main uses is in displaying the consituent solids of a Boolean
900      // solid and their surfaces overlap with the resulting Booean
901      // solid, making a mess if hlr is specified.
902      style = G4ViewParameters::wireframe;
903      break;
904    }
905  }
906  return style;
907}
908
909G4bool G4VSceneHandler::GetAuxEdgeVisible (const G4VisAttributes* pVisAttribs) {
910  G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
911  if (pVisAttribs -> IsForceAuxEdgeVisible()) isAuxEdgeVisible = true;
912  return isAuxEdgeVisible;
913}
914
915G4double G4VSceneHandler::GetMarkerSize
916(const G4VMarker& marker,
917 G4VSceneHandler::MarkerSizeType& markerSizeType)
918{
919  G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
920  const G4VMarker& defaultMarker =
921    fpViewer -> GetViewParameters().GetDefaultMarker();
922  G4double size = userSpecified ?
923    marker.GetWorldSize() : defaultMarker.GetWorldSize();
924  if (size) {
925    // Draw in world coordinates.
926    markerSizeType = world;
927  }
928  else {
929    size = userSpecified ?
930      marker.GetScreenSize() : defaultMarker.GetScreenSize();
931    // Draw in screen coordinates.
932    markerSizeType = screen;
933  }
934  if (size <= 1.) size = 1.;
935  size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
936  if (size <= 1.) size = 1.;
937  return size;
938}
939
940G4int G4VSceneHandler::GetNoOfSides(const G4VisAttributes* pVisAttribs)
941{
942  // No. of sides (lines segments per circle) is normally determined
943  // by the view parameters, but it can be overriddden by the
944  // ForceLineSegmentsPerCircle in the vis attributes.
945  G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
946  if (pVisAttribs) {
947    if (pVisAttribs->IsForceLineSegmentsPerCircle())
948      lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
949    const G4int nSegmentsMin = 12;
950    if (lineSegmentsPerCircle < nSegmentsMin) {
951      lineSegmentsPerCircle = nSegmentsMin;
952      G4cout <<
953        "G4VSceneHandler::GetNoOfSides: attempt to set the"
954        "\nnumber of line segements per circle < " << nSegmentsMin
955             << "; forced to " << lineSegmentsPerCircle << G4endl;
956    }
957  }
958  return lineSegmentsPerCircle;
959}
960
961std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s) {
962
963  os << "Scene handler " << s.fName << " has "
964     << s.fViewerList.size () << " viewer(s):";
965  for (size_t i = 0; i < s.fViewerList.size (); i++) {
966    os << "\n  " << *(s.fViewerList [i]);
967  }
968
969  if (s.fpScene) {
970    os << "\n  " << *s.fpScene;
971  }
972  else {
973    os << "\n  This scene handler currently has no scene.";
974  }
975
976  return os;
977}
Note: See TracBrowser for help on using the repository browser.