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

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

otpimisations

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