source: tags/Visualization_after-vis09-02-01-tag/OpenInventor/src/G4OpenInventorSceneHandler.cc @ 958

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

test pour GL_POINTS au lieu de glBitmap

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