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

Last change on this file since 959 was 959, checked in by garnier, 15 years ago

after maj on CVS good-HEAD

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