Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/visualization/OpenInventor/src/G4OpenInventorSceneHandler.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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