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