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