Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4OpenInventorSceneHandler.cc
Go to the documentation of this file.
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$
28 //
29 //
30 // Jeff Kallenbach 01 Aug 1996
31 // OpenInventor stored scene - creates OpenInventor display lists.
32 #ifdef G4VIS_BUILD_OI_DRIVER
33 
34 // this :
36 
37 #include <Inventor/SoPath.h>
38 #include <Inventor/SoNodeKitPath.h>
39 #include <Inventor/nodes/SoCoordinate3.h>
40 #include <Inventor/nodes/SoCoordinate4.h>
41 #include <Inventor/nodes/SoSeparator.h>
42 #include <Inventor/nodes/SoDrawStyle.h>
43 #include <Inventor/nodes/SoLightModel.h>
44 #include <Inventor/nodes/SoMaterial.h>
45 #include <Inventor/nodes/SoLineSet.h>
46 #include <Inventor/nodes/SoCube.h>
47 #include <Inventor/nodes/SoFont.h>
48 #include <Inventor/nodes/SoText2.h>
49 #include <Inventor/nodes/SoFaceSet.h>
50 #include <Inventor/nodes/SoNormal.h>
51 #include <Inventor/nodes/SoNormalBinding.h>
52 #include <Inventor/nodes/SoComplexity.h>
53 #include <Inventor/nodes/SoNurbsSurface.h>
54 #include <Inventor/nodes/SoTranslation.h>
55 #include <Inventor/nodes/SoTransform.h>
56 #include <Inventor/nodes/SoResetTransform.h>
57 #include <Inventor/nodes/SoMatrixTransform.h>
58 
59 #define USE_SOPOLYHEDRON
60 
61 #ifndef USE_SOPOLYHEDRON
62 #include "HEPVis/nodes/SoBox.h"
63 #include "HEPVis/nodes/SoTubs.h"
64 #include "HEPVis/nodes/SoCons.h"
65 #include "HEPVis/nodes/SoTrd.h"
66 #include "HEPVis/nodes/SoTrap.h"
67 #endif
69 typedef HEPVis_SoMarkerSet SoMarkerSet;
72 
73 #include "SoG4Polyhedron.h"
74 #include "SoG4LineSet.h"
75 #include "SoG4MarkerSet.h"
76 
77 #include "G4Scene.hh"
78 #include "G4NURBS.hh"
79 #include "G4OpenInventor.hh"
81 #include "G4ThreeVector.hh"
82 #include "G4Point3D.hh"
83 #include "G4Normal3D.hh"
84 #include "G4Transform3D.hh"
85 #include "G4Polyline.hh"
86 #include "G4Text.hh"
87 #include "G4Circle.hh"
88 #include "G4Square.hh"
89 #include "G4Polymarker.hh"
90 #include "G4Polyhedron.hh"
91 #include "G4Box.hh"
92 #include "G4Tubs.hh"
93 #include "G4Cons.hh"
94 #include "G4Trap.hh"
95 #include "G4Trd.hh"
96 #include "G4ModelingParameters.hh"
97 #include "G4VPhysicalVolume.hh"
98 #include "G4LogicalVolume.hh"
99 #include "G4Material.hh"
100 #include "G4VisAttributes.hh"
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;
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 = line.size();
206  SbVec3f* pCoords = new SbVec3f[nPoints];
207 
208  for (G4int iPoint = 0; iPoint < nPoints ; iPoint++) {
209  pCoords[iPoint].setValue((float)line[iPoint].x(),
210  (float)line[iPoint].y(),
211  (float)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 long *)&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;
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 = 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((float)polymarker[iPoint].x(),
270  (float)polymarker[iPoint].y(),
271  (float)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;
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;
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;
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((float)c.GetRed(),
384  (float)c.GetGreen(),
385  (float)c.GetBlue(),
386  (float)(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;
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((float)centre.x(),
481  (float)centre.y(),
482  (float)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;
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 //
589 // Method for handling G4NURBS objects for drawing solids.
590 // Knots and Ctrl Pnts MUST be arrays of GLfloats.
591 //
592 void G4OpenInventorSceneHandler::AddPrimitive (const G4NURBS& nurb) {
593 
594  if (fProcessing2D) {
595  static G4bool warned = false;
596  if (!warned) {
597  warned = true;
599  ("G4OpenInventorSceneHandler::AddPrimitive (const G4NURBS&)",
600  "OpenInventor-0006", JustWarning,
601  "2D NURBS not implemented. Ignored.");
602  }
603  return;
604  }
605 
606  // Get vis attributes - pick up defaults if none.
607  const G4VisAttributes* pVA =
608  fpViewer -> GetApplicableVisAttributes (nurb.GetVisAttributes ());
609 
610  AddProperties(pVA); // Colour, etc.
611  AddTransform(); // Transformation
612 
613  G4float *u_knot_array, *u_knot_array_ptr;
614  u_knot_array = u_knot_array_ptr = new G4float [nurb.GetnbrKnots(G4NURBS::U)];
615  G4NURBS::KnotsIterator u_iterator (nurb, G4NURBS::U);
616  while (u_iterator.pick (u_knot_array_ptr++)){}
617 
618  G4float *v_knot_array, *v_knot_array_ptr;
619  v_knot_array = v_knot_array_ptr = new G4float [nurb.GetnbrKnots(G4NURBS::V)];
620  G4NURBS::KnotsIterator v_iterator (nurb, G4NURBS::V);
621  while (v_iterator.pick (v_knot_array_ptr++)){}
622 
623  G4float *ctrl_pnt_array, *ctrl_pnt_array_ptr;
624  ctrl_pnt_array = ctrl_pnt_array_ptr =
625  new G4float [nurb.GettotalnbrCtrlPts () * G4NURBS::NofC*sizeof(G4float)];
626  G4NURBS::CtrlPtsCoordsIterator c_p_iterator (nurb);
627  while (c_p_iterator.pick (ctrl_pnt_array_ptr++)){}
628 
629  SoSeparator *surfSep = new SoSeparator();
630 
631  //
632  // Set up NURBS
633  //
634  SoComplexity *complexity = new SoComplexity;
635  SoCoordinate4 *ctlPts = new SoCoordinate4;
636  SoNurbsSurface *oi_nurb = new SoNurbsSurface;
637 
638  complexity->value = (float)0.6;
639  G4int nPoints = nurb.GettotalnbrCtrlPts ();
640  SbVec4f* points = new SbVec4f[nPoints];
641  for (G4int iPoint = 0; iPoint < nPoints ; iPoint++) {
642  points[iPoint].setValue(
643  ctrl_pnt_array[iPoint*4 + 0],
644  ctrl_pnt_array[iPoint*4 + 1],
645  ctrl_pnt_array[iPoint*4 + 2],
646  ctrl_pnt_array[iPoint*4 + 3]);
647  }
648  ctlPts->point.setValues (0,nPoints,points);
649  oi_nurb->numUControlPoints = nurb.GetnbrCtrlPts(G4NURBS::U);
650  oi_nurb->numVControlPoints = nurb.GetnbrCtrlPts(G4NURBS::V);
651  oi_nurb->uKnotVector.setValues(0,nurb.GetnbrKnots(G4NURBS::U),u_knot_array);
652  oi_nurb->vKnotVector.setValues(0,nurb.GetnbrKnots(G4NURBS::V),v_knot_array);
653 
654  surfSep->addChild(complexity);
655  surfSep->addChild(ctlPts);
656  surfSep->addChild(oi_nurb);
657 
658  fCurrentSeparator->addChild(surfSep);
659 
660  //
661  // Clean-up
662  //
663  delete [] u_knot_array;
664  delete [] v_knot_array;
665  delete [] ctrl_pnt_array;
666  delete [] points;
667 }
668 
669 void G4OpenInventorSceneHandler::GeneratePrerequisites()
670 {
671  // Utility for PreAddSolid and BeginPrimitives.
672 
673  // This routines prepares for adding items to the scene database. We
674  // are expecting two kinds of solids: leaf parts and non-leaf parts.
675  // For non-leaf parts, we create a detector tree kit. This has two
676  // separators. The solid itself goes in the preview separator, the
677  // full separator is forseen for daughters. This separator is not
678  // only created--it is also put in a dictionary for future use by
679  // the leaf part.
680 
681  // For leaf parts, we first locate the mother volume and find its
682  // separator through the dictionary.
683 
684  // The private member fCurrentSeparator is set to the proper
685  // location on in the scene database so that when the solid is
686  // actually added (in addthis), it is put in the right place.
687 
688  G4PhysicalVolumeModel* pPVModel =
689  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
690 
691  if (pPVModel) {
692 
693  // This call comes from a G4PhysicalVolumeModel. drawnPVPath is
694  // the path of the current drawn (non-culled) volume in terms of
695  // drawn (non-culled) ancesters. Each node is identified by a
696  // PVNodeID object, which is a physical volume and copy number. It
697  // is a vector of PVNodeIDs corresponding to the geometry hierarchy
698  // actually selected, i.e., not culled.
700  typedef std::vector<PVNodeID> PVPath;
701  const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
702  //G4int currentDepth = pPVModel->GetCurrentDepth();
703  G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
704  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
705  //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
706  // Note: pCurrentMaterial may be zero (parallel world).
707 
708  // The simplest algorithm, used by the Open Inventor Driver
709  // developers, is to rely on the fact the G4PhysicalVolumeModel
710  // traverses the geometry hierarchy in an orderly manner. The last
711  // mother, if any, will be the node to which the volume should be
712  // added. So it is enough to keep a map of scene graph nodes keyed
713  // on the volume path ID. Actually, it is enough to use the logical
714  // volume as the key. (An alternative would be to keep the PVNodeID
715  // in the tree and match the PVPath from the root down.)
716 
717  // Find mother. ri points to mother, if any...
718  PVPath::const_reverse_iterator ri;
719  G4LogicalVolume* MotherVolume = 0;
720  ri = ++drawnPVPath.rbegin();
721  if (ri != drawnPVPath.rend()) {
722  // This volume has a mother.
723  MotherVolume = ri->GetPhysicalVolume()->GetLogicalVolume();
724  }
725 
726  if (pCurrentLV->GetNoDaughters()!=0 ||
727  pCurrentPV->IsReplicated()) { //????Don't understand this???? JA
728  // This block of code is executed for non-leaf parts:
729 
730  // Make the detector tree kit:
731  SoDetectorTreeKit* detectorTreeKit = new SoDetectorTreeKit();
732 
733  SoSeparator* previewSeparator =
734  (SoSeparator*) detectorTreeKit->getPart("previewSeparator",TRUE);
735  previewSeparator->renderCaching = SoSeparator::OFF;
736 
737  SoSeparator* fullSeparator =
738  (SoSeparator*) detectorTreeKit->getPart("fullSeparator", TRUE);
739  fullSeparator->renderCaching = SoSeparator::OFF;
740 
741  if(fPreviewAndFull) detectorTreeKit->setPreviewAndFull();
742  else detectorTreeKit->setPreview(TRUE);
743 
744  // Colour, etc., for SoDetectorTreeKit. Treated differently to
745  // othere SoNodes(?). Use fpVisAttribs stored away in
746  // PreAddSolid...
747  const G4VisAttributes* pApplicableVisAttribs =
748  fpViewer->GetApplicableVisAttributes (fpVisAttribs);
749 
750  // First find the color attributes...
751  const G4Colour& g4Col = pApplicableVisAttribs->GetColour ();
752  const double red = g4Col.GetRed ();
753  const double green = g4Col.GetGreen ();
754  const double blue = g4Col.GetBlue ();
755  double transparency = 1 - g4Col.GetAlpha();
756 
757  // Drawing style...
758  G4ViewParameters::DrawingStyle drawing_style =
759  GetDrawingStyle(pApplicableVisAttribs);
760  switch (drawing_style) {
762  fModelingSolid = false;
763  break;
764  case (G4ViewParameters::hlr):
765  case (G4ViewParameters::hsr):
767  fModelingSolid = true;
768  break;
769  }
770 
771  SoMaterial* material =
772  fStyleCache->getMaterial((float)red,
773  (float)green,
774  (float)blue,
775  (float)transparency);
776  detectorTreeKit->setPart("appearance.material",material);
777 
778  SoLightModel* lightModel =
779  fModelingSolid ? fStyleCache->getLightModelPhong() :
780  fStyleCache->getLightModelBaseColor();
781  detectorTreeKit->setPart("appearance.lightModel",lightModel);
782 
783  // Add the full separator to the dictionary; it is indexed by the
784  // address of the logical volume!
785  fSeparatorMap[pCurrentLV] = fullSeparator;
786 
787  // Find out where to add this volume.
788  // If no mother can be found, it goes under root.
789  if (MotherVolume) {
790  if (fSeparatorMap.find(MotherVolume) != fSeparatorMap.end()) {
791  //printf("debug : PreAddSolid : mother %s found in map\n",
792  // MotherVolume->GetName().c_str());
793  fSeparatorMap[MotherVolume]->addChild(detectorTreeKit);
794  } else {
795  // Mother not previously encountered. Shouldn't happen, since
796  // G4PhysicalVolumeModel sends volumes as it encounters them,
797  // i.e., mothers before daughters, in its descent of the
798  // geometry tree. Error!
799  G4cout <<
800  "ERROR: G4OpenInventorSceneHandler::GeneratePrerequisites: Mother "
801  << ri->GetPhysicalVolume()->GetName()
802  << ':' << ri->GetCopyNo()
803  << " not previously encountered."
804  "\nShouldn't happen! Please report to visualization coordinator."
805  << G4endl;
806  // Continue anyway. Add to root of scene graph tree...
807  //printf("debug : PreAddSolid : mother %s not found in map !!!\n",
808  // MotherVolume->GetName().c_str());
809  fDetectorRoot->addChild(detectorTreeKit);
810  }
811  } else {
812  //printf("debug : PreAddSolid : has no mother\n");
813  fDetectorRoot->addChild(detectorTreeKit);
814  }
815 
816  fCurrentSeparator = previewSeparator;
817 
818  } else {
819  // This block of code is executed for leaf parts.
820 
821  if (MotherVolume) {
822  if (fSeparatorMap.find(MotherVolume) != fSeparatorMap.end()) {
823  fCurrentSeparator = fSeparatorMap[MotherVolume];
824  } else {
825  // Mother not previously encountered. Shouldn't happen, since
826  // G4PhysicalVolumeModel sends volumes as it encounters them,
827  // i.e., mothers before daughters, in its descent of the
828  // geometry tree. Error!
829  G4cout << "ERROR: G4OpenInventorSceneHandler::PreAddSolid: Mother "
830  << ri->GetPhysicalVolume()->GetName()
831  << ':' << ri->GetCopyNo()
832  << " not previously encountered."
833  "\nShouldn't happen! Please report to visualization coordinator."
834  << G4endl;
835  // Continue anyway. Add to root of scene graph tree...
836  fCurrentSeparator = fDetectorRoot;
837  }
838  } else {
839  fCurrentSeparator = fDetectorRoot;
840  }
841  }
842 
843  } else {
844  // Not from G4PhysicalVolumeModel, so add to root as leaf part...
845 
846  if (fReadyForTransients) {
847  fCurrentSeparator = fTransientRoot;
848  } else {
849  fCurrentSeparator = fDetectorRoot;
850  }
851  }
852 }
853 
854 void G4OpenInventorSceneHandler::AddProperties(const G4VisAttributes* visAtts)
855 {
856  // Use the applicable vis attributes...
857  const G4VisAttributes* pApplicableVisAttribs =
858  fpViewer->GetApplicableVisAttributes (visAtts);
859 
860  // First find the color attributes...
861  const G4Colour& g4Col = pApplicableVisAttribs->GetColour ();
862  const double red = g4Col.GetRed ();
863  const double green = g4Col.GetGreen ();
864  const double blue = g4Col.GetBlue ();
865  double transparency = 1 - g4Col.GetAlpha();
866 
867  // Drawing style...
868  G4ViewParameters::DrawingStyle drawing_style =
869  GetDrawingStyle(pApplicableVisAttribs);
870  switch (drawing_style) {
872  fModelingSolid = false;
873  break;
874  case (G4ViewParameters::hlr):
875  case (G4ViewParameters::hsr):
877  fModelingSolid = true;
878  break;
879  }
880 
881  // Edge visibility...
882  G4bool isAuxEdgeVisible = GetAuxEdgeVisible (pApplicableVisAttribs);
883  fReducedWireFrame = !isAuxEdgeVisible;
884 
885  SoMaterial* material =
886  fStyleCache->getMaterial((float)red,
887  (float)green,
888  (float)blue,
889  (float)transparency);
890  fCurrentSeparator->addChild(material);
891 
892  SoLightModel* lightModel =
893  fModelingSolid ? fStyleCache->getLightModelPhong() :
894  fStyleCache->getLightModelBaseColor();
895  fCurrentSeparator->addChild(lightModel);
896 }
897 
898 void G4OpenInventorSceneHandler::AddTransform(const G4Point3D& translation)
899 {
900  // AddTransform takes fObjectTransformation and "adds" a translation.
901  // Set up the geometrical transformation for the coming
902  fCurrentSeparator->addChild(fStyleCache->getResetTransform());
903 
904  SoMatrixTransform* matrixTransform = new SoMatrixTransform;
905  G4OpenInventorTransform3D oiTran
906  (fObjectTransformation * G4Translate3D(translation));
907  SbMatrix* sbMatrix = oiTran.GetSbMatrix();
908 
909  const G4Vector3D scale = fpViewer->GetViewParameters().GetScaleFactor();
910  SbMatrix sbScale;
911  sbScale.setScale
912  (SbVec3f((float)scale.x(),(float)scale.y(),(float)scale.z()));
913  sbMatrix->multRight(sbScale);
914 
915  matrixTransform->matrix.setValue(*sbMatrix);
916  delete sbMatrix;
917  fCurrentSeparator->addChild(matrixTransform);
918 }
919 #endif