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