source: trunk/source/visualization/OpenInventor/src/G4OpenInventorSceneHandler.cc @ 924

Last change on this file since 924 was 850, checked in by garnier, 16 years ago

geant4.8.2 beta

  • Property svn:mime-type set to text/cpp
File size: 26.5 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: G4OpenInventorSceneHandler.cc,v 1.54 2008/04/04 13:40:04 allison Exp $
28// GEANT4 tag $Name: HEAD $
29//
30//
31// Jeff Kallenbach 01 Aug 1996
32// OpenInventor stored scene - creates OpenInventor display lists.
33#ifdef G4VIS_BUILD_OI_DRIVER
34
35// this :
36#include "G4OpenInventorSceneHandler.hh"
37
38#include <Inventor/SoPath.h>
39#include <Inventor/SoNodeKitPath.h>
40#include <Inventor/nodes/SoCoordinate3.h>
41#include <Inventor/nodes/SoCoordinate4.h>
42#include <Inventor/nodes/SoSeparator.h>
43#include <Inventor/nodes/SoDrawStyle.h>
44#include <Inventor/nodes/SoLightModel.h>
45#include <Inventor/nodes/SoMaterial.h>
46#include <Inventor/nodes/SoLineSet.h>
47#include <Inventor/nodes/SoCube.h>
48#include <Inventor/nodes/SoSphere.h>
49#include <Inventor/nodes/SoFont.h>
50#include <Inventor/nodes/SoText2.h>
51#include <Inventor/nodes/SoFaceSet.h>
52#include <Inventor/nodes/SoNormal.h>
53#include <Inventor/nodes/SoNormalBinding.h>
54#include <Inventor/nodes/SoComplexity.h>
55#include <Inventor/nodes/SoNurbsSurface.h>
56#include <Inventor/nodes/SoTranslation.h>
57#include <Inventor/nodes/SoTransform.h>
58#include <Inventor/nodes/SoResetTransform.h>
59#include <Inventor/nodes/SoMatrixTransform.h>
60
61#define USE_SOPOLYHEDRON
62
63#ifndef USE_SOPOLYHEDRON
64#include "HEPVis/nodes/SoBox.h"
65#include "HEPVis/nodes/SoTubs.h"
66#include "HEPVis/nodes/SoCons.h"
67#include "HEPVis/nodes/SoTrd.h"
68#include "HEPVis/nodes/SoTrap.h"
69#endif
70#include "HEPVis/nodes/SoMarkerSet.h"
71typedef HEPVis_SoMarkerSet SoMarkerSet;
72#include "HEPVis/nodekits/SoDetectorTreeKit.h"
73#include "HEPVis/misc/SoStyleCache.h"
74
75#include "SoG4Polyhedron.h"
76#include "SoG4LineSet.h"
77#include "SoG4MarkerSet.h"
78
79#include "G4Scene.hh"
80#include "G4NURBS.hh"
81#include "G4OpenInventor.hh"
82#include "G4OpenInventorTransform3D.hh"
83#include "G4ThreeVector.hh"
84#include "G4Point3D.hh"
85#include "G4Normal3D.hh"
86#include "G4Transform3D.hh"
87#include "G4Polyline.hh"
88#include "G4Text.hh"
89#include "G4Circle.hh"
90#include "G4Square.hh"
91#include "G4Polymarker.hh"
92#include "G4Polyhedron.hh"
93#include "G4Box.hh"
94#include "G4Tubs.hh"
95#include "G4Cons.hh"
96#include "G4Trap.hh"
97#include "G4Trd.hh"
98#include "G4ModelingParameters.hh"
99#include "G4VPhysicalVolume.hh"
100#include "G4LogicalVolume.hh"
101#include "G4Material.hh"
102#include "G4VisAttributes.hh"
103
104G4int G4OpenInventorSceneHandler::fSceneIdCount = 0;
105
106G4OpenInventorSceneHandler::G4OpenInventorSceneHandler (G4OpenInventor& system,
107                                          const G4String& name)
108:G4VSceneHandler (system, fSceneIdCount++, name)
109,fRoot(0)
110,fDetectorRoot(0)
111,fTransientRoot(0)
112,fCurrentSeparator(0)
113,fModelingSolid(false)
114,fReducedWireFrame(true)
115,fStyleCache(0)
116,fPreviewAndFull(false)
117{
118  fStyleCache = new SoStyleCache;
119  fStyleCache->ref();
120
121  fRoot = new SoSeparator;
122  fRoot->ref();
123  fRoot->setName("Root");
124 
125  fDetectorRoot = new SoSeparator;
126  fDetectorRoot->setName("StaticRoot");
127  fRoot->addChild(fDetectorRoot);
128 
129  fTransientRoot = new SoSeparator;
130  fTransientRoot->setName("TransientRoot");
131  fRoot->addChild(fTransientRoot);
132 
133  fCurrentSeparator = fTransientRoot;
134}
135
136G4OpenInventorSceneHandler::~G4OpenInventorSceneHandler ()
137{
138  fRoot->unref();
139  fStyleCache->unref();
140}
141
142void G4OpenInventorSceneHandler::ClearStore () {
143  G4VSceneHandler::ClearStore();
144
145  fDetectorRoot->removeAllChildren();
146  fSeparatorMap.clear();
147
148  fTransientRoot->removeAllChildren();
149}
150
151void G4OpenInventorSceneHandler::ClearTransientStore () {
152  G4VSceneHandler::ClearTransientStore ();
153
154  fTransientRoot->removeAllChildren();
155}
156
157//
158// Generates prerequisites for solids
159// 
160void G4OpenInventorSceneHandler::PreAddSolid
161(const G4Transform3D& objectTransformation,
162 const G4VisAttributes& visAttribs)
163{
164  G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
165  // Stores arguments away for future use, e.g., AddPrimitives.
166
167  GeneratePrerequisites();
168}
169
170//
171// Generates prerequisites for primitives
172// 
173void G4OpenInventorSceneHandler::BeginPrimitives
174(const G4Transform3D& objectTransformation) {
175
176  G4VSceneHandler::BeginPrimitives (objectTransformation);
177
178  // If thread of control has already passed through PreAddSolid,
179  // avoid opening a graphical data base component again.
180  if (!fProcessingSolid) {
181    GeneratePrerequisites();
182  }
183}
184
185//
186// Method for handling G4Polyline objects (from tracking).
187//
188void G4OpenInventorSceneHandler::AddPrimitive (const G4Polyline& line)
189{
190  AddProperties(line.GetVisAttributes());  // Transformation, colour, etc.
191
192  G4int nPoints = line.size();
193  SbVec3f* pCoords = new SbVec3f[nPoints];
194
195  for (G4int iPoint = 0; iPoint < nPoints ; iPoint++) {
196    pCoords[iPoint].setValue((float)line[iPoint].x(),
197                             (float)line[iPoint].y(),
198                             (float)line[iPoint].z());
199  }
200
201  //
202  // Point Set
203  //
204  SoCoordinate3 *polyCoords = new SoCoordinate3;
205  polyCoords->point.setValues(0,nPoints,pCoords);
206  fCurrentSeparator->addChild(polyCoords);
207 
208  //
209  // Wireframe
210  //
211  SoDrawStyle* drawStyle = fStyleCache->getLineStyle();
212  fCurrentSeparator->addChild(drawStyle);
213
214  SoG4LineSet *pLine = new SoG4LineSet;
215
216  // Loads G4Atts for picking...
217  if (fpViewer->GetViewParameters().IsPicking()) LoadAtts(line, pLine);
218
219#ifdef INVENTOR2_0
220  pLine->numVertices.setValues(0,1,(const long *)&nPoints);
221#else
222  pLine->numVertices.setValues(0,1,&nPoints);
223#endif
224
225  fCurrentSeparator->addChild(pLine);
226
227  delete [] pCoords;
228}
229
230void G4OpenInventorSceneHandler::AddPrimitive (const G4Polymarker& polymarker)
231{
232  AddProperties(polymarker.GetVisAttributes()); // Transformation, colour, etc.
233
234  G4int pointn = polymarker.size();
235  if(pointn<=0) return;
236
237  SbVec3f* points = new SbVec3f[pointn];
238  for (G4int iPoint = 0; iPoint < pointn ; iPoint++) {
239    points[iPoint].setValue((float)polymarker[iPoint].x(),
240                            (float)polymarker[iPoint].y(),
241                            (float)polymarker[iPoint].z());
242  }
243
244  SoCoordinate3* coordinate3 = new SoCoordinate3;
245  coordinate3->point.setValues(0,pointn,points);
246  fCurrentSeparator->addChild(coordinate3);
247
248  MarkerSizeType sizeType;
249  G4double screenSize = GetMarkerSize (polymarker, sizeType);
250  switch (sizeType) {
251  default:
252  case screen:
253    // Draw in screen coordinates.  OK.
254    break;
255  case world:
256    // Draw in world coordinates.   Not implemented.  Use screenSize = 10.
257    screenSize = 10.;
258    break;
259  }
260 
261  SoG4MarkerSet* markerSet = new SoG4MarkerSet;
262  markerSet->numPoints = pointn;
263
264  // Loads G4Atts for picking...
265  if (fpViewer->GetViewParameters().IsPicking())
266    LoadAtts(polymarker, markerSet);
267
268  G4VMarker::FillStyle style = polymarker.GetFillStyle();
269  switch (polymarker.GetMarkerType()) {
270  default:
271    // Are available 5_5, 7_7 and 9_9
272  case G4Polymarker::dots:
273    if (screenSize <= 5.) {
274      markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_5_5;
275    } else if (screenSize <= 7.) {
276      markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_7_7;
277    } else {
278      markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_9_9;
279    }
280    break;
281  case G4Polymarker::circles:
282    if (screenSize <= 5.) {
283      if (style == G4VMarker::filled) {
284        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_5_5;
285      } else {
286        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_5_5;
287      }
288    } else if (screenSize <= 7.) {
289      if (style == G4VMarker::filled) {
290        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_7_7;
291      } else {
292        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_7_7;
293      }
294    } else {
295      if (style == G4VMarker::filled) {
296        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_9_9;
297      } else {
298        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_9_9;
299      }
300    }
301    break;
302  case G4Polymarker::squares:
303    if (screenSize <= 5.) {
304      if (style == G4VMarker::filled) {
305        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_5_5;
306      } else {
307        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_5_5;
308      }
309    } else if (screenSize <= 7.) {
310      if (style == G4VMarker::filled) {
311        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_7_7;
312      } else {
313        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_7_7;
314      }
315    } else {
316      if (style == G4VMarker::filled) {
317        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_9_9;
318      } else {
319        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_9_9;
320      }
321    }
322  }
323  fCurrentSeparator->addChild(markerSet);
324
325  delete [] points;
326}
327
328// ********* NOTE ********* NOTE ********* NOTE ********* NOTE *********
329//
330//  This method (Text) has not been tested, as it is
331//  innaccessible from the menu in the current configuration
332//
333//  Currently draws at the origin!  How do I get it to draw at
334//  text.GetPosition()?  JA
335//
336// ********* NOTE ********* NOTE ********* NOTE ********* NOTE *********
337//
338// Method for handling G4Text objects
339//
340void G4OpenInventorSceneHandler::AddPrimitive (const G4Text& text)
341{
342  AddProperties(text.GetVisAttributes());  // Transformation, colour, etc.
343
344  //
345  // Color.  Note: text colour is worked out differently.  This
346  // over-rides the colour added in AddProperties...
347  //
348  const G4Colour& c = GetTextColour (text);
349  SoMaterial* material =
350    fStyleCache->getMaterial((float)c.GetRed(),
351                             (float)c.GetGreen(),
352                             (float)c.GetBlue(),
353                             (float)(1-c.GetAlpha()));
354  fCurrentSeparator->addChild(material);
355
356  MarkerSizeType sizeType;
357  G4double size = GetMarkerSize (text, sizeType);
358  switch (sizeType) {
359  default:
360  case screen:
361    // Draw in screen coordinates.  OK.
362    break;
363  case world:
364    // Draw in world coordinates.   Not implemented.  Use size = 20.
365    size = 20.;
366    break;
367  }
368
369  //
370  // Font
371  //
372  SoFont *g4Font = new SoFont();
373  g4Font->size = size;
374  fCurrentSeparator->addChild(g4Font);
375
376  //
377  // Text
378  //
379  SoText2 *g4String = new SoText2();
380  g4String->string.setValue(text.GetText());
381  g4String->spacing = 2.0;
382  switch (text.GetLayout()) {
383  default:
384  case G4Text::left:
385    g4String->justification = SoText2::LEFT; break;
386  case G4Text::centre:
387    g4String->justification = SoText2::CENTER; break;
388  case G4Text::right:
389    g4String->justification = SoText2::RIGHT; break;
390  }
391  fCurrentSeparator->addChild(g4String);
392}
393
394//
395// Method for handling G4Circle objects
396//
397void G4OpenInventorSceneHandler::AddPrimitive (const G4Circle& circle) {
398  AddCircleSquare(G4OICircle, circle);
399}
400
401//
402// Method for handling G4Square objects - defaults to wireframe
403//
404void G4OpenInventorSceneHandler::AddPrimitive (const G4Square& square) {
405  AddCircleSquare(G4OISquare, square);
406}
407
408void G4OpenInventorSceneHandler::AddCircleSquare
409(G4OIMarker markerType, const G4VMarker& marker)
410{
411  AddProperties(marker.GetVisAttributes());  // Transformation, colour, etc.
412
413  MarkerSizeType sizeType;
414  G4double screenSize = GetMarkerSize (marker, sizeType);
415  switch (sizeType) {
416  default:
417  case screen:
418    // Draw in screen coordinates.  OK.
419    break;
420  case world:
421    // Draw in world coordinates.   Not implemented.  Use size = 10.
422    screenSize = 10.;
423    break;
424  }
425
426  G4Point3D centre = marker.GetPosition();
427
428  // Borrowed from AddPrimitive(G4Polymarker) - inefficient? JA
429  SbVec3f* points = new SbVec3f[1];
430  points[0].setValue((float)centre.x(),
431                     (float)centre.y(),
432                     (float)centre.z());
433  SoCoordinate3* coordinate3 = new SoCoordinate3;
434  coordinate3->point.setValues(0,1,points);
435  fCurrentSeparator->addChild(coordinate3);
436
437  SoG4MarkerSet* markerSet = new SoG4MarkerSet;
438  markerSet->numPoints = 1;
439
440  // Loads G4Atts for picking...
441  if (fpViewer->GetViewParameters().IsPicking()) LoadAtts(marker, markerSet);
442
443  G4VMarker::FillStyle style = marker.GetFillStyle();
444  switch (markerType) {
445  case G4OICircle:
446    if (screenSize <= 5.) {
447      if (style == G4VMarker::filled) {
448        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_5_5;
449      } else {
450        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_5_5;
451      }
452    } else if (screenSize <= 7.) {
453      if (style == G4VMarker::filled) {
454        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_7_7;
455      } else {
456        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_7_7;
457      }
458    } else {
459      if (style == G4VMarker::filled) {
460        markerSet->markerIndex = SoMarkerSet::CIRCLE_FILLED_9_9;
461      } else {
462        markerSet->markerIndex = SoMarkerSet::CIRCLE_LINE_9_9;
463      }
464    }
465    break;
466  case G4OISquare:
467    if (screenSize <= 5.) {
468      if (style == G4VMarker::filled) {
469        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_5_5;
470      } else {
471        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_5_5;
472      }
473    } else if (screenSize <= 7.) {
474      if (style == G4VMarker::filled) {
475        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_7_7;
476      } else {
477        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_7_7;
478      }
479    } else {
480      if (style == G4VMarker::filled) {
481        markerSet->markerIndex = SoMarkerSet::SQUARE_FILLED_9_9;
482      } else {
483        markerSet->markerIndex = SoMarkerSet::SQUARE_LINE_9_9;
484      }
485    }
486  break;
487  }
488  fCurrentSeparator->addChild(markerSet);
489
490  delete [] points;
491}
492
493//
494// Method for handling G4Polyhedron objects for drawing solids.
495//
496void G4OpenInventorSceneHandler::AddPrimitive (const G4Polyhedron& polyhedron)
497{
498  if (polyhedron.GetNoFacets() == 0) return;
499
500  AddProperties(polyhedron.GetVisAttributes()); // Transformation, colour, etc.
501
502  SoG4Polyhedron* soPolyhedron = new SoG4Polyhedron(polyhedron);
503
504  // Loads G4Atts for picking...
505  if (fpViewer->GetViewParameters().IsPicking())
506    LoadAtts(polyhedron, soPolyhedron);
507
508  SbString name = "Non-geometry";
509  G4PhysicalVolumeModel* pPVModel =
510    dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
511  if (pPVModel) {
512    name = pPVModel->GetCurrentLV()->GetName().c_str();
513  }
514  SbName sbName(name);
515  soPolyhedron->setName(sbName);
516  soPolyhedron->solid.setValue(fModelingSolid);
517  soPolyhedron->reducedWireFrame.setValue(fReducedWireFrame?TRUE:FALSE);
518  fCurrentSeparator->addChild(soPolyhedron); 
519}
520
521//
522// Method for handling G4NURBS objects for drawing solids.
523// Knots and Ctrl Pnts MUST be arrays of GLfloats.
524//
525void G4OpenInventorSceneHandler::AddPrimitive (const G4NURBS& nurb) {
526
527  AddProperties(nurb.GetVisAttributes()); // Transformation, colour, etc.
528
529  G4float *u_knot_array, *u_knot_array_ptr;
530  u_knot_array = u_knot_array_ptr = new G4float [nurb.GetnbrKnots(G4NURBS::U)];
531  G4NURBS::KnotsIterator u_iterator (nurb, G4NURBS::U);
532  while (u_iterator.pick (u_knot_array_ptr++)){}
533
534  G4float *v_knot_array, *v_knot_array_ptr;
535  v_knot_array = v_knot_array_ptr = new G4float [nurb.GetnbrKnots(G4NURBS::V)];
536  G4NURBS::KnotsIterator v_iterator (nurb, G4NURBS::V);
537  while (v_iterator.pick (v_knot_array_ptr++)){}
538
539  G4float *ctrl_pnt_array, *ctrl_pnt_array_ptr;
540  ctrl_pnt_array = ctrl_pnt_array_ptr =
541    new G4float [nurb.GettotalnbrCtrlPts () * G4NURBS::NofC*sizeof(G4float)];
542  G4NURBS::CtrlPtsCoordsIterator c_p_iterator (nurb);
543  while (c_p_iterator.pick (ctrl_pnt_array_ptr++)){}
544 
545  SoSeparator *surfSep = new SoSeparator();
546
547  //
548  // Set up NURBS
549  //
550  SoComplexity *complexity = new SoComplexity;
551  SoCoordinate4 *ctlPts = new SoCoordinate4;
552  SoNurbsSurface *oi_nurb = new SoNurbsSurface;
553 
554  complexity->value = (float)0.6;
555  G4int    nPoints = nurb.GettotalnbrCtrlPts ();
556  SbVec4f* points  = new SbVec4f[nPoints];
557  for (G4int iPoint = 0; iPoint < nPoints ; iPoint++) {
558    points[iPoint].setValue(
559                            ctrl_pnt_array[iPoint*4 + 0],
560                            ctrl_pnt_array[iPoint*4 + 1],
561                            ctrl_pnt_array[iPoint*4 + 2],
562                            ctrl_pnt_array[iPoint*4 + 3]);
563  }
564  ctlPts->point.setValues (0,nPoints,points);
565  oi_nurb->numUControlPoints = nurb.GetnbrCtrlPts(G4NURBS::U);
566  oi_nurb->numVControlPoints = nurb.GetnbrCtrlPts(G4NURBS::V);
567  oi_nurb->uKnotVector.setValues(0,nurb.GetnbrKnots(G4NURBS::U),u_knot_array);
568  oi_nurb->vKnotVector.setValues(0,nurb.GetnbrKnots(G4NURBS::V),v_knot_array);
569
570  surfSep->addChild(complexity);
571  surfSep->addChild(ctlPts);
572  surfSep->addChild(oi_nurb);
573 
574  fCurrentSeparator->addChild(surfSep);
575
576  //
577  // Clean-up
578  //
579  delete [] u_knot_array;
580  delete [] v_knot_array;
581  delete [] ctrl_pnt_array;
582  delete [] points;
583}
584
585void G4OpenInventorSceneHandler::GeneratePrerequisites()
586{
587  // Utility for PreAddSolid and BeginPrimitives.
588
589  // This routines prepares for adding items to the scene database.  We
590  // are expecting two kinds of solids: leaf parts and non-leaf parts.
591  // For non-leaf parts, we create a detector tree kit.  This has two
592  // separators.  The solid itself goes in the preview separator, the
593  // full separator is forseen for daughters.  This separator is not
594  // only created--it is also put in a dictionary for future use by
595  // the leaf part.
596
597  // For leaf parts, we first locate the mother volume and find its
598  // separator through the dictionary.
599
600  // The private member fCurrentSeparator is set to the proper
601  // location on in the scene database so that when the solid is
602  // actually added (in addthis), it is put in the right place.
603
604  G4PhysicalVolumeModel* pPVModel =
605    dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
606 
607  if (pPVModel) {
608
609    // This call comes from a G4PhysicalVolumeModel.  drawnPVPath is
610    // the path of the current drawn (non-culled) volume in terms of
611    // drawn (non-culled) ancesters.  Each node is identified by a
612    // PVNodeID object, which is a physical volume and copy number.  It
613    // is a vector of PVNodeIDs corresponding to the geometry hierarchy
614    // actually selected, i.e., not culled.
615    typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
616    typedef std::vector<PVNodeID> PVPath;
617    const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
618    //G4int currentDepth = pPVModel->GetCurrentDepth();
619    G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
620    G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
621    //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
622    // Note: pCurrentMaterial may be zero (parallel world).
623
624    // The simplest algorithm, used by the Open Inventor Driver
625    // developers, is to rely on the fact the G4PhysicalVolumeModel
626    // traverses the geometry hierarchy in an orderly manner.  The last
627    // mother, if any, will be the node to which the volume should be
628    // added.  So it is enough to keep a map of scene graph nodes keyed
629    // on the volume path ID.  Actually, it is enough to use the logical
630    // volume as the key.  (An alternative would be to keep the PVNodeID
631    // in the tree and match the PVPath from the root down.)
632
633    // Find mother.  ri points to mother, if any...
634    PVPath::const_reverse_iterator ri;
635    G4LogicalVolume* MotherVolume = 0;
636    ri = ++drawnPVPath.rbegin();
637    if (ri != drawnPVPath.rend()) {
638      // This volume has a mother.
639      MotherVolume = ri->GetPhysicalVolume()->GetLogicalVolume();
640    }
641
642    if (pCurrentLV->GetNoDaughters()!=0 ||
643        pCurrentPV->IsReplicated()) {  //????Don't understand this???? JA
644      // This block of code is executed for non-leaf parts:
645
646      // Make the detector tree kit:
647      SoDetectorTreeKit* detectorTreeKit = new SoDetectorTreeKit(); 
648
649      SoSeparator* previewSeparator   = 
650        (SoSeparator*) detectorTreeKit->getPart("previewSeparator",TRUE);
651      previewSeparator->renderCaching = SoSeparator::OFF;
652
653      SoSeparator* fullSeparator = 
654        (SoSeparator*) detectorTreeKit->getPart("fullSeparator",   TRUE);
655      fullSeparator->renderCaching = SoSeparator::OFF;
656
657      if(fPreviewAndFull) detectorTreeKit->setPreviewAndFull();
658      else detectorTreeKit->setPreview(TRUE);
659
660      // Colour, etc., for SoDetectorTreeKit.  Treated differently to
661      // othere SoNodes(?).  Use fpVisAttribs stored away in
662      // PreAddSolid...
663      const G4VisAttributes* pApplicableVisAttribs =
664        fpViewer->GetApplicableVisAttributes (fpVisAttribs);
665
666      // First find the color attributes...
667      const G4Colour& g4Col =  pApplicableVisAttribs->GetColour ();
668      const double red = g4Col.GetRed ();
669      const double green = g4Col.GetGreen ();
670      const double blue = g4Col.GetBlue ();
671      double transparency = 1 - g4Col.GetAlpha();
672
673      // Drawing style...
674      G4ViewParameters::DrawingStyle drawing_style =
675        GetDrawingStyle(pApplicableVisAttribs);
676      switch (drawing_style) {
677      case (G4ViewParameters::wireframe):   
678        fModelingSolid = false;
679        break;
680      case (G4ViewParameters::hlr):
681      case (G4ViewParameters::hsr):
682      case (G4ViewParameters::hlhsr):
683        fModelingSolid = true;
684        break;
685      }
686
687      SoMaterial* material =
688        fStyleCache->getMaterial((float)red,
689                                 (float)green,
690                                 (float)blue,
691                                 (float)transparency);
692      detectorTreeKit->setPart("appearance.material",material);
693
694      SoLightModel* lightModel =
695        fModelingSolid ? fStyleCache->getLightModelPhong() :
696        fStyleCache->getLightModelBaseColor();
697      detectorTreeKit->setPart("appearance.lightModel",lightModel);
698
699      // Add the full separator to the dictionary; it is indexed by the
700      // address of the logical volume!
701      fSeparatorMap[pCurrentLV] = fullSeparator;
702
703      // Find out where to add this volume.
704      // If no mother can be found, it goes under root.
705      if (MotherVolume) {
706        if (fSeparatorMap.find(MotherVolume) != fSeparatorMap.end()) {
707          //printf("debug : PreAddSolid : mother %s found in map\n",
708          //     MotherVolume->GetName().c_str());
709          fSeparatorMap[MotherVolume]->addChild(detectorTreeKit);
710        } else {
711          // Mother not previously encountered.  Shouldn't happen, since
712          // G4PhysicalVolumeModel sends volumes as it encounters them,
713          // i.e., mothers before daughters, in its descent of the
714          // geometry tree.  Error!
715          G4cout <<
716            "ERROR: G4OpenInventorSceneHandler::GeneratePrerequisites: Mother "
717                 << ri->GetPhysicalVolume()->GetName()
718                 << ':' << ri->GetCopyNo()
719                 << " not previously encountered."
720            "\nShouldn't happen!  Please report to visualization coordinator."
721                 << G4endl;
722          // Continue anyway.  Add to root of scene graph tree...
723          //printf("debug : PreAddSolid : mother %s not found in map !!!\n",
724          //     MotherVolume->GetName().c_str());
725          fDetectorRoot->addChild(detectorTreeKit);
726        }
727      } else {
728        //printf("debug : PreAddSolid : has no mother\n");
729        fDetectorRoot->addChild(detectorTreeKit);
730      }
731
732      fCurrentSeparator = previewSeparator;
733
734    } else {
735      // This block of code is executed for leaf parts.
736
737      if (MotherVolume) {
738        if (fSeparatorMap.find(MotherVolume) != fSeparatorMap.end()) {
739          fCurrentSeparator = fSeparatorMap[MotherVolume];
740        } else {
741          // Mother not previously encountered.  Shouldn't happen, since
742          // G4PhysicalVolumeModel sends volumes as it encounters them,
743          // i.e., mothers before daughters, in its descent of the
744          // geometry tree.  Error!
745          G4cout << "ERROR: G4OpenInventorSceneHandler::PreAddSolid: Mother "
746                 << ri->GetPhysicalVolume()->GetName()
747                 << ':' << ri->GetCopyNo()
748                 << " not previously encountered."
749            "\nShouldn't happen!  Please report to visualization coordinator."
750                 << G4endl;
751          // Continue anyway.  Add to root of scene graph tree...
752          fCurrentSeparator = fDetectorRoot;
753        }
754      } else {
755        fCurrentSeparator = fDetectorRoot;
756      }
757    }
758
759  } else {
760    // Not from G4PhysicalVolumeModel, so add to root as leaf part...
761
762    if (fReadyForTransients) {
763      fCurrentSeparator = fTransientRoot;
764    } else {
765      fCurrentSeparator = fDetectorRoot;
766    }
767  }
768}
769
770void G4OpenInventorSceneHandler::AddProperties(const G4VisAttributes* visAtts)
771{
772  // Use the applicable vis attributes...
773  const G4VisAttributes* pApplicableVisAttribs =
774    fpViewer->GetApplicableVisAttributes (visAtts);
775
776  // First find the color attributes...
777  const G4Colour& g4Col =  pApplicableVisAttribs->GetColour ();
778  const double red = g4Col.GetRed ();
779  const double green = g4Col.GetGreen ();
780  const double blue = g4Col.GetBlue ();
781  double transparency = 1 - g4Col.GetAlpha();
782
783  // Drawing style...
784  G4ViewParameters::DrawingStyle drawing_style =
785    GetDrawingStyle(pApplicableVisAttribs);
786  switch (drawing_style) {
787  case (G4ViewParameters::wireframe):   
788    fModelingSolid = false;
789    break;
790  case (G4ViewParameters::hlr):
791  case (G4ViewParameters::hsr):
792  case (G4ViewParameters::hlhsr):
793    fModelingSolid = true;
794    break;
795  }
796
797  // Edge visibility...
798  G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pApplicableVisAttribs);
799  fReducedWireFrame = !isAuxEdgeVisible;
800
801  SoMaterial* material =
802    fStyleCache->getMaterial((float)red,
803                             (float)green,
804                             (float)blue,
805                             (float)transparency);
806  fCurrentSeparator->addChild(material);
807
808  SoLightModel* lightModel =
809    fModelingSolid ? fStyleCache->getLightModelPhong() :
810    fStyleCache->getLightModelBaseColor();
811  fCurrentSeparator->addChild(lightModel);
812
813  // Set up the geometrical transformation for the coming
814  fCurrentSeparator->addChild(fStyleCache->getResetTransform());
815
816  SoMatrixTransform* matrixTransform = new SoMatrixTransform;
817  G4OpenInventorTransform3D oiTran(*fpObjectTransformation);
818  SbMatrix* sbMatrix = oiTran.GetSbMatrix();
819
820  const G4Vector3D scale = fpViewer->GetViewParameters().GetScaleFactor();
821  SbMatrix sbScale;
822  sbScale.setScale
823    (SbVec3f((float)scale.x(),(float)scale.y(),(float)scale.z()));
824  sbMatrix->multRight(sbScale);
825
826  matrixTransform->matrix.setValue(*sbMatrix);
827  delete sbMatrix;
828  fCurrentSeparator->addChild(matrixTransform);
829}
830#endif
Note: See TracBrowser for help on using the repository browser.