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

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

update

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