source: trunk/source/visualization/XXX/src/G4XXXSGSceneHandler.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 13.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: G4XXXSGSceneHandler.cc,v 1.6 2006/11/05 20:41:36 allison Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
29//
30//
31// John Allison  10th March 2006
32// A template for a sophisticated graphics driver with a scene graph.
33//?? Lines beginning like this require specialisation for your driver.
34
35#ifdef G4VIS_BUILD_XXXSG_DRIVER
36
37#include "G4XXXSGSceneHandler.hh"
38
39#include "G4XXXSGViewer.hh"
40#include "G4PhysicalVolumeModel.hh"
41#include "G4VPhysicalVolume.hh"
42#include "G4LogicalVolume.hh"
43#include "G4Box.hh"
44#include "G4Polyline.hh"
45#include "G4Text.hh"
46#include "G4Circle.hh"
47#include "G4Square.hh"
48#include "G4Polyhedron.hh"
49#include "G4UnitsTable.hh"
50
51#include <sstream>
52
53G4int G4XXXSGSceneHandler::fSceneIdCount = 0;
54// Counter for XXX scene handlers.
55
56G4XXXSGSceneHandler::G4XXXSGSceneHandler(G4VGraphicsSystem& system,
57                                             const G4String& name):
58  G4VSceneHandler(system, fSceneIdCount++, name)
59{
60  fRoot = fSceneGraph.setRoot("\nroot\n");
61  fPermanentsRoot = fRoot.push_back("\npermanentsRoot\n");
62  fTransientsRoot = fRoot.push_back("\ntransientsRoot\n");
63}
64
65G4XXXSGSceneHandler::~G4XXXSGSceneHandler() {}
66
67#ifdef G4XXXSGDEBUG
68// Useful function...
69void G4XXXSGSceneHandler::PrintThings() {
70  G4cout <<
71    "  with transformation "
72         << (void*)fpObjectTransformation;
73  if (fpModel) {
74    G4cout << " from " << fpModel->GetCurrentDescription()
75           << " (tag " << fpModel->GetCurrentTag()
76           << ')';
77  } else {
78    G4cout << "(not from a model)";
79  }
80  G4PhysicalVolumeModel* pPVModel =
81    dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
82  if (pPVModel) {
83    G4cout <<
84      "\n  current physical volume: "
85           << pPVModel->GetCurrentPV()->GetName() <<
86      "\n  current logical volume: "
87           << pPVModel->GetCurrentLV()->GetName() <<
88      "\n  current depth of geometry tree: "
89           << pPVModel->GetCurrentDepth();
90  }
91  G4cout << G4endl;
92}
93#endif
94
95void G4XXXSGSceneHandler::CreateCurrentItem(const G4String& header) {
96  // Utility for PreAddSolid and BeginPrimitives.
97 
98  G4PhysicalVolumeModel* pPVModel =
99    dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
100 
101  if (pPVModel) {
102
103    // This call comes from a G4PhysicalVolumeModel.  drawnPVPath is
104    // the path of the current drawn (non-culled) volume in terms of
105    // drawn (non-culled) ancesters.  Each node is identified by a
106    // PVNodeID object, which is a physical volume and copy number.  It
107    // is a vector of PVNodeIDs corresponding to the geometry hierarchy
108    // actually selected, i.e., not culled.
109    typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
110    typedef std::vector<PVNodeID> PVPath;
111    const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
112    //G4int currentDepth = pPVModel->GetCurrentDepth();
113    //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
114    G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
115    //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
116    // Note: pCurrentMaterial may be zero (parallel world).
117
118    // The simplest algorithm, used by the Open Inventor Driver
119    // developers, is to rely on the fact the G4PhysicalVolumeModel
120    // traverses the geometry hierarchy in an orderly manner.  The last
121    // mother, if any, will be the node to which the volume should be
122    // added.  So it is enough to keep a map of scene graph nodes keyed
123    // on the volume path ID.  Actually, it is enough to use the logical
124    // volume as the key.  (An alternative would be to keep the PVNodeID
125    // in the tree and match the PVPath from the root down.)
126
127    // Find mother.  ri points to drawn mother, if any.
128    PVPath::const_reverse_iterator ri = ++drawnPVPath.rbegin();
129    if (ri != drawnPVPath.rend()) {
130      // This volume has a mother.
131      G4LogicalVolume* drawnMotherLV =
132        ri->GetPhysicalVolume()->GetLogicalVolume();
133      LVMapIterator mother = fLVMap.find(drawnMotherLV);
134      if (mother != fLVMap.end()) {
135        // This adds a child in Troy's tree...
136        fCurrentItem = mother->second.push_back(header);
137      } else {
138        // Mother not previously encountered.  Shouldn't happen, since
139        // G4PhysicalVolumeModel sends volumes as it encounters them,
140        // i.e., mothers before daughters, in its descent of the
141        // geometry tree.  Error!
142        G4cout << "ERROR: G4XXXSGSceneHandler::PreAddSolid: Mother "
143               << ri->GetPhysicalVolume()->GetName()
144               << ':' << ri->GetCopyNo()
145               << " not previously encountered."
146          "\nShouldn't happen!  Please report to visualization coordinator."
147               << G4endl;
148        // Continue anyway.  Add to root of scene graph tree...
149        fCurrentItem = fPermanentsRoot.push_back(header);
150      }
151    } else {
152      // This volume has no mother.  Must be a top level un-culled
153      // volume.  Add to root of scene graph tree...
154      fCurrentItem = fPermanentsRoot.push_back(header);
155    }
156
157    std::ostringstream oss;
158    oss << "Path of drawn PVs: ";
159    for (PVPath::const_iterator i = drawnPVPath.begin();
160         i != drawnPVPath.end(); ++i) {
161      oss << '/' << i->GetPhysicalVolume()->GetName()
162          << ':' << i->GetCopyNo();
163    }
164    oss << std::endl;
165    *fCurrentItem += oss.str();
166
167    // Store for future searches.  Overwrites previous entries for this
168    // LV, so entry is always the *last* LV.
169    fLVMap[pCurrentLV] = fCurrentItem;
170
171  } else {  // Not from a G4PhysicalVolumeModel.
172
173    // Create a place for current solid in root...
174    if (fReadyForTransients) {
175      fCurrentItem = fTransientsRoot.push_back(header);
176    } else {
177      fCurrentItem = fPermanentsRoot.push_back(header);
178    }
179  }
180}
181
182void G4XXXSGSceneHandler::PreAddSolid
183(const G4Transform3D& objectTransformation,
184 const G4VisAttributes& visAttribs)
185{ 
186  G4VSceneHandler::PreAddSolid(objectTransformation, visAttribs);
187  CreateCurrentItem(G4String("\nPreAddSolid:\n"));
188}
189
190void G4XXXSGSceneHandler::PostAddSolid()
191{
192  *fCurrentItem += "EndSolid\n";
193  G4VSceneHandler::PostAddSolid();
194}
195
196void G4XXXSGSceneHandler::BeginPrimitives
197(const G4Transform3D& objectTransformation)
198{
199  G4VSceneHandler::BeginPrimitives(objectTransformation);
200
201  // If thread of control has already passed through PreAddSolid,
202  // avoid opening a graphical data base component again.
203  if (!fProcessingSolid) {
204    CreateCurrentItem(G4String("\nBeginPrimitives:\n"));
205  }
206}
207
208void G4XXXSGSceneHandler::EndPrimitives ()
209{
210  if (!fProcessingSolid) {  // Already done if so.
211    *fCurrentItem += "EndPrimitives\n";
212  }
213  G4VSceneHandler::EndPrimitives ();
214}
215
216// Note: This function overrides G4VSceneHandler::AddSolid(const
217// G4Box&).  You may not want to do this, but this is how it's done if
218// you do.  Certain other specific solids may be treated this way -
219// see G4VSceneHandler.hh.  The simplest possible driver would *not*
220// implement these polymorphic functions, with the effect that the
221// default versions in G4VSceneHandler are used, which simply call
222// G4VSceneHandler::RequestPrimitives to turn the solid into a
223// G4Polyhedron usually.
224// Don't forget, solids can be transients too (e.g., representing a hit).
225void G4XXXSGSceneHandler::AddSolid(const G4Box& box) {
226#ifdef G4XXXSGDEBUG
227  G4cout <<
228    "G4XXXSGSceneHandler::AddSolid(const G4Box& box) called for "
229         << box.GetName()
230         << G4endl;
231#endif
232  //?? Process your box...
233  std::ostringstream oss;
234  oss << "G4Box(" <<
235    G4String
236    (G4BestUnit
237     (G4ThreeVector
238      (box.GetXHalfLength(), box.GetYHalfLength(), box.GetZHalfLength()),
239      "Length")).strip() << ')' << std::endl;
240  *fCurrentItem += oss.str();
241}
242
243void G4XXXSGSceneHandler::AddPrimitive(const G4Polyline& polyline) {
244#ifdef G4XXXSGDEBUG
245  G4cout <<
246    "G4XXXSGSceneHandler::AddPrimitive(const G4Polyline& polyline) called.\n"
247         << polyline
248         << G4endl;
249#endif
250  // Get vis attributes - pick up defaults if none.
251  //const G4VisAttributes* pVA =
252  //  fpViewer -> GetApplicableVisAttributes (polyline.GetVisAttributes ());
253  //?? Process polyline.
254  std::ostringstream oss;
255  oss << polyline << std::endl;
256  *fCurrentItem += oss.str();
257}
258
259void G4XXXSGSceneHandler::AddPrimitive(const G4Text& text) {
260#ifdef G4XXXSGDEBUG
261  G4cout <<
262    "G4XXXSGSceneHandler::AddPrimitive(const G4Text& text) called.\n"
263         << text
264         << G4endl;
265#endif
266  // Get text colour - special method since default text colour is
267  // determined by the default text vis attributes, which may be
268  // specified independent of default vis attributes of other types of
269  // visible objects.
270  //const G4Colour& c = GetTextColour (text);  // Picks up default if none.
271  //?? Process text.
272  std::ostringstream oss;
273  oss << text << std::endl;
274  *fCurrentItem += oss.str();
275}
276
277void G4XXXSGSceneHandler::AddPrimitive(const G4Circle& circle) {
278#ifdef G4XXXSGDEBUG
279  G4cout <<
280    "G4XXXSGSceneHandler::AddPrimitive(const G4Circle& circle) called.\n"
281         << circle
282         << G4endl; 
283  MarkerSizeType sizeType;
284  G4double size = GetMarkerSize (circle, sizeType);
285  switch (sizeType) {
286  default:
287  case screen:
288    // Draw in screen coordinates.
289    G4cout << "screen";
290    break;
291  case world:
292    // Draw in world coordinates.
293    G4cout << "world";
294    break;
295  }
296  G4cout << " size: " << size << G4endl;
297#endif
298  // Get vis attributes - pick up defaults if none.
299  //const G4VisAttributes* pVA =
300  //  fpViewer -> GetApplicableVisAttributes (circle.GetVisAttributes ());
301  //?? Process circle.
302  std::ostringstream oss;
303  oss << circle << std::endl;
304  *fCurrentItem += oss.str();
305}
306
307void G4XXXSGSceneHandler::AddPrimitive(const G4Square& square) {
308#ifdef G4XXXSGDEBUG
309  G4cout <<
310    "G4XXXSGSceneHandler::AddPrimitive(const G4Square& square) called.\n"
311         << square
312         << G4endl;
313  MarkerSizeType sizeType;
314  G4double size = GetMarkerSize (square, sizeType);
315  switch (sizeType) {
316  default:
317  case screen:
318    // Draw in screen coordinates.
319    G4cout << "screen";
320    break;
321  case world:
322    // Draw in world coordinates.
323    G4cout << "world";
324    break;
325  }
326  G4cout << " size: " << size << G4endl;
327#endif
328  // Get vis attributes - pick up defaults if none.
329  //const G4VisAttributes* pVA =
330  //  fpViewer -> GetApplicableVisAttributes (square.GetVisAttributes ());
331  //?? Process square.
332  std::ostringstream oss;
333  oss << square << std::endl;
334  *fCurrentItem += oss.str();
335}
336
337void G4XXXSGSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) {
338#ifdef G4XXXSGDEBUG
339  G4cout <<
340 "G4XXXSGSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron) called.\n"
341         << polyhedron
342         << G4endl;
343#endif
344  //?? Process polyhedron.
345  std::ostringstream oss;
346  oss << polyhedron;
347  *fCurrentItem += oss.str();
348
349  //?? Or... here are some ideas for decomposing into polygons...
350  //Assume all facets are convex quadrilaterals.
351  //Draw each G4Facet individually
352 
353  //Get colour, etc..
354  if (polyhedron.GetNoFacets() == 0) return;
355
356  // Get vis attributes - pick up defaults if none.
357  const G4VisAttributes* pVA =
358    fpViewer -> GetApplicableVisAttributes (polyhedron.GetVisAttributes ());
359
360  // Get view parameters that the user can force through the vis
361  // attributes, thereby over-riding the current view parameter.
362  G4ViewParameters::DrawingStyle drawing_style = GetDrawingStyle (pVA);
363  //G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pVA);
364 
365  //Get colour, etc..
366  //const G4Colour& c = pVA -> GetColour ();
367 
368  // Initial action depending on drawing style.
369  switch (drawing_style) {
370  case (G4ViewParameters::hsr):
371    {
372      break;
373    }
374  case (G4ViewParameters::hlr):
375    {
376      break;
377    }
378  case (G4ViewParameters::wireframe):
379    {
380      break;
381    }
382  default:
383    {
384      break;
385    }     
386  }
387
388  // Loop through all the facets...
389
390  // Look at G4OpenGLSceneHandler::AddPrimitive(const G4Polyhedron&)
391  // for an example of how to get facets out of a G4Polyhedron,
392  // including how to cope with triangles if that's a problem.
393}
394
395void G4XXXSGSceneHandler::AddPrimitive(const G4NURBS& nurbs) {
396#ifdef G4XXXSGDEBUG
397  G4cout <<
398    "G4XXXSGSceneHandler::AddPrimitive(const G4NURBS& nurbs) called."
399         << G4endl;
400#endif
401  //?? Don't bother implementing this.  NURBS are not functional.
402}
403
404void G4XXXSGSceneHandler::ClearStore () {
405
406  G4VSceneHandler::ClearStore ();  // Sets need kernel visit, etc.
407
408  fPermanentsRoot.clear();
409  fTransientsRoot.clear();
410  fLVMap.clear();
411}
412
413void G4XXXSGSceneHandler::ClearTransientStore () {
414
415  G4VSceneHandler::ClearTransientStore ();
416
417  fTransientsRoot.clear();
418}
419
420#endif
Note: See TracBrowser for help on using the repository browser.