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

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

cvs update

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