source: trunk/geant4/visualization/OpenInventor/src/G4OpenInventorSceneHandler.cc @ 552

Last change on this file since 552 was 529, checked in by garnier, 17 years ago

r658@mac-90108: laurentgarnier | 2007-06-25 12:02:16 +0200
import de visualisation

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