source: trunk/source/visualization/OpenInventor/src/G4OpenInventorSceneHandler.cc@ 930

Last change on this file since 930 was 929, checked in by garnier, 17 years ago

juste pour sauver le boulot ....ne compile pas....

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