Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhysicalVolumeModel Class Reference

#include <G4PhysicalVolumeModel.hh>

+ Inheritance diagram for G4PhysicalVolumeModel:

Classes

class  G4PhysicalVolumeModelTouchable
 
class  G4PhysicalVolumeNodeID
 

Public Types

enum  { UNLIMITED = -1 }
 
enum  ClippingMode { subtraction , intersection }
 

Public Member Functions

 G4PhysicalVolumeModel (G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false)
 
virtual ~G4PhysicalVolumeModel ()
 
void DescribeYourselfTo (G4VGraphicsScene &)
 
G4String GetCurrentDescription () const
 
G4String GetCurrentTag () const
 
G4VPhysicalVolumeGetTopPhysicalVolume () const
 
G4int GetRequestedDepth () const
 
const G4VSolidGetClippingSolid () const
 
G4int GetCurrentDepth () const
 
G4VPhysicalVolumeGetCurrentPV () const
 
G4LogicalVolumeGetCurrentLV () const
 
G4MaterialGetCurrentMaterial () const
 
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath () const
 
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath () const
 
const std::map< G4String, G4AttDef > * GetAttDefs () const
 
std::vector< G4AttValue > * CreateCurrentAttValues () const
 
void SetBaseFullPVPath (const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath)
 
void SetRequestedDepth (G4int requestedDepth)
 
void SetClippingSolid (G4VSolid *pClippingSolid)
 
void SetClippingMode (ClippingMode mode)
 
G4bool Validate (G4bool warn)
 
void CurtailDescent ()
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
virtual void DescribeYourselfTo (G4VGraphicsScene &)=0
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () const
 
virtual G4String GetCurrentDescription () const
 
virtual G4String GetCurrentTag () const
 
const G4VisExtentGetExtent () const
 
const G4StringGetGlobalDescription () const
 
const G4StringGetGlobalTag () const
 
const G4Transform3DGetTransformation () const
 
void SetModelingParameters (const G4ModelingParameters *)
 
void SetExtent (const G4VisExtent &)
 
void SetType (const G4String &)
 
void SetGlobalDescription (const G4String &)
 
void SetGlobalTag (const G4String &)
 
void SetTransformation (const G4Transform3D &)
 
virtual G4bool Validate (G4bool warn=true)
 

Protected Member Functions

void VisitGeometryAndGetVisReps (G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
 
void DescribeAndDescend (G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
 
virtual void DescribeSolid (const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
 
void CalculateExtent ()
 

Protected Attributes

G4VPhysicalVolumefpTopPV
 
G4String fTopPVName
 
G4int fTopPVCopyNo
 
G4int fRequestedDepth
 
G4bool fUseFullExtent
 
G4int fCurrentDepth
 
G4VPhysicalVolumefpCurrentPV
 
G4LogicalVolumefpCurrentLV
 
G4MaterialfpCurrentMaterial
 
G4Transform3DfpCurrentTransform
 
std::vector< G4PhysicalVolumeNodeIDfBaseFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfFullPVPath
 
std::vector< G4PhysicalVolumeNodeIDfDrawnPVPath
 
G4bool fCurtailDescent
 
G4VSolidfpClippingSolid
 
ClippingMode fClippingMode
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
G4Transform3D fTransform
 
const G4ModelingParametersfpMP
 

Detailed Description

Definition at line 82 of file G4PhysicalVolumeModel.hh.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
UNLIMITED 

Definition at line 86 of file G4PhysicalVolumeModel.hh.

◆ ClippingMode

Constructor & Destructor Documentation

◆ G4PhysicalVolumeModel()

G4PhysicalVolumeModel::G4PhysicalVolumeModel ( G4VPhysicalVolume pVPV = 0,
G4int  requestedDepth = UNLIMITED,
const G4Transform3D modelTransformation = G4Transform3D(),
const G4ModelingParameters pMP = 0,
G4bool  useFullExtent = false 
)

Definition at line 58 of file G4PhysicalVolumeModel.cc.

63 :
64 G4VModel (modelTransformation, pMP),
65 fpTopPV (pVPV),
66 fRequestedDepth (requestedDepth),
67 fUseFullExtent (useFullExtent),
68 fCurrentDepth (0),
69 fpCurrentPV (0),
70 fpCurrentLV (0),
73 fCurtailDescent (false),
76{
77 if (fpTopPV) {
78
79 fType = "G4PhysicalVolumeModel";
80 fTopPVName = fpTopPV -> GetName ();
81 fTopPVCopyNo = fpTopPV -> GetCopyNo ();
82 std::ostringstream o;
83 o << fpTopPV -> GetCopyNo ();
84 fGlobalTag = fpTopPV -> GetName () + "." + o.str();
85 fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
86
88 }
89}
G4VPhysicalVolume * fpTopPV
G4Transform3D * fpCurrentTransform
G4VPhysicalVolume * fpCurrentPV
G4String fGlobalDescription
Definition: G4VModel.hh:110
G4String fType
Definition: G4VModel.hh:108
G4String fGlobalTag
Definition: G4VModel.hh:109

◆ ~G4PhysicalVolumeModel()

G4PhysicalVolumeModel::~G4PhysicalVolumeModel ( )
virtual

Definition at line 91 of file G4PhysicalVolumeModel.cc.

92{
93 delete fpClippingSolid;
94}

Member Function Documentation

◆ CalculateExtent()

void G4PhysicalVolumeModel::CalculateExtent ( )
protected

Definition at line 96 of file G4PhysicalVolumeModel.cc.

97{
98 if (fUseFullExtent) {
99 fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
100 }
101 else {
102 G4BoundingSphereScene bsScene(this);
103 const G4int tempRequestedDepth = fRequestedDepth;
104 fRequestedDepth = -1; // Always search to all depths to define extent.
105 const G4ModelingParameters* tempMP = fpMP;
107 (0, // No default vis attributes needed.
108 G4ModelingParameters::wf, // wireframe (not relevant for this).
109 true, // Global culling.
110 true, // Cull invisible volumes.
111 false, // Density culling.
112 0., // Density (not relevant if density culling false).
113 true, // Cull daughters of opaque mothers.
114 24); // No of sides (not relevant for this operation).
115 fpMP = &mParams;
116 DescribeYourselfTo (bsScene);
117 G4double radius = bsScene.GetRadius();
118 if (radius < 0.) { // Nothing in the scene.
119 fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
120 } else {
121 // Transform back to coordinates relative to the top
122 // transformation, which is in G4VModel::fTransform. This makes
123 // it conform to all models, which are defined by a
124 // transformation and an extent relative to that
125 // transformation...
126 G4Point3D centre = bsScene.GetCentre();
127 centre.transform(fTransform.inverse());
128 fExtent = G4VisExtent(centre, radius);
129 }
130 fpMP = tempMP;
131 fRequestedDepth = tempRequestedDepth;
132 }
133}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void DescribeYourselfTo(G4VGraphicsScene &)
G4VisExtent fExtent
Definition: G4VModel.hh:111
const G4VisExtent & GetExtent() const
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
G4Transform3D fTransform
Definition: G4VModel.hh:112
Transform3D inverse() const
Definition: Transform3D.cc:142

Referenced by G4PhysicalVolumeModel(), and Validate().

◆ CreateCurrentAttValues()

std::vector< G4AttValue > * G4PhysicalVolumeModel::CreateCurrentAttValues ( ) const

Definition at line 829 of file G4PhysicalVolumeModel.cc.

830{
831 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
832 std::ostringstream oss;
833 for (size_t i = 0; i < fFullPVPath.size(); ++i) {
834 oss << fFullPVPath[i].GetPhysicalVolume()->GetName()
835 << ':' << fFullPVPath[i].GetCopyNo();
836 if (i != fFullPVPath.size() - 1) oss << '/';
837 }
838
839 if (!fpCurrentLV) {
841 ("G4PhysicalVolumeModel::CreateCurrentAttValues",
842 "modeling0004",
844 "Current logical volume not defined.");
845 return values;
846 }
847
848 values->push_back(G4AttValue("PVPath", oss.str(),""));
849 values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
850 G4VSolid* pSol = fpCurrentLV->GetSolid();
851 values->push_back(G4AttValue("Solid", pSol->GetName(),""));
852 values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
853 oss.str(""); oss << '\n' << *pSol;
854 values->push_back(G4AttValue("DmpSol", oss.str(),""));
855 oss.str(""); oss << '\n' << *fpCurrentTransform;
856 values->push_back(G4AttValue("Trans", oss.str(),""));
857 G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
858 values->push_back(G4AttValue("Material", matName,""));
860 values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
862 oss.str(""); oss << matState;
863 values->push_back(G4AttValue("State", oss.str(),""));
865 values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
866 G4Region* region = fpCurrentLV->GetRegion();
867 G4String regionName = region? region->GetName(): G4String("No region");
868 values->push_back(G4AttValue("Region", regionName,""));
869 oss.str(""); oss << fpCurrentLV->IsRootRegion();
870 values->push_back(G4AttValue("RootRegion", oss.str(),""));
871 return values;
872}
@ JustWarning
G4State
Definition: G4Material.hh:114
@ kStateUndefined
Definition: G4Material.hh:114
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4VSolid * GetSolid() const
G4String GetName() const
G4bool IsRootRegion() const
G4Region * GetRegion() const
G4double GetDensity() const
Definition: G4Material.hh:179
G4State GetState() const
Definition: G4Material.hh:180
G4double GetRadlen() const
Definition: G4Material.hh:219
const G4String & GetName() const
Definition: G4Material.hh:177
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
const G4String & GetName() const
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G4VSceneHandler::LoadAtts(), and G4XXXStoredSceneHandler::PreAddSolid().

◆ CurtailDescent()

void G4PhysicalVolumeModel::CurtailDescent ( )
inline

Definition at line 222 of file G4PhysicalVolumeModel.hh.

222{fCurtailDescent = true;}

Referenced by G4ASCIITreeSceneHandler::RequestPrimitives().

◆ DescribeAndDescend()

void G4PhysicalVolumeModel::DescribeAndDescend ( G4VPhysicalVolume pVPV,
G4int  requestedDepth,
G4LogicalVolume pLV,
G4VSolid pSol,
G4Material pMaterial,
const G4Transform3D theAT,
G4VGraphicsScene sceneHandler 
)
protected

Definition at line 334 of file G4PhysicalVolumeModel.cc.

342{
343 // Maintain useful data members...
344 fpCurrentPV = pVPV;
345 fpCurrentLV = pLV;
346 fpCurrentMaterial = pMaterial;
347
348 const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue ();
349 const G4ThreeVector& translation = pVPV -> GetTranslation ();
350 G4Transform3D theLT (G4Transform3D (objectRotation, translation));
351
352 // Compute the accumulated transformation...
353 // Note that top volume's transformation relative to the world
354 // coordinate system is specified in theAT == startingTransformation
355 // = fTransform (see DescribeYourselfTo), so first time through the
356 // volume's own transformation, which is only relative to its
357 // mother, i.e., not relative to the world coordinate system, should
358 // not be accumulated.
359 G4Transform3D theNewAT (theAT);
360 if (fCurrentDepth != 0) theNewAT = theAT * theLT;
361 fpCurrentTransform = &theNewAT;
362
363 const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes();
364 if (!pVisAttribs) pVisAttribs = fpMP->GetDefaultVisAttributes();
365 // Beware - pVisAttribs might still be zero - create a temporary default one...
366 G4bool visAttsCreated = false;
367 if (!pVisAttribs) {
368 pVisAttribs = new G4VisAttributes;
369 visAttsCreated = true;
370 }
371 // From here, can assume pVisAttribs is a valid pointer.
372
373 // Make decision to draw...
374 G4bool thisToBeDrawn = true;
375
376 // Update full path of physical volumes...
377 G4int copyNo = fpCurrentPV->GetCopyNo();
378 fFullPVPath.push_back
379 (G4PhysicalVolumeNodeID
381
382 // In case we need to copy the vis atts for modification...
383 G4bool copyForVAM = false;
384 const G4VisAttributes* pUnmodifiedVisAtts = pVisAttribs;
385 G4VisAttributes* pModifiedVisAtts = 0;
386
387 // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
388 const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
390 std::vector<G4ModelingParameters::VisAttributesModifier>::const_iterator
391 iModifier;
392 for (iModifier = vams.begin();
393 iModifier != vams.end();
394 ++iModifier) {
396 iModifier->GetPVNameCopyNoPath();
397 if (vamPath.size() == fFullPVPath.size()) {
398 // OK - there's a size match. Check it out.
399// G4cout << "Size match" << G4endl;
401 std::vector<G4PhysicalVolumeNodeID>::const_iterator iPVNodeId;
402 for (iVAMNameCopyNo = vamPath.begin(), iPVNodeId = fFullPVPath.begin();
403 iVAMNameCopyNo != vamPath.end();
404 ++iVAMNameCopyNo, ++iPVNodeId) {
405// G4cout
406// << iVAMNameCopyNo->fName
407// << ',' << iVAMNameCopyNo->fCopyNo
408// << "; " << iPVNodeId->GetPhysicalVolume()->GetName()
409// << ',' << iPVNodeId->GetPhysicalVolume()->GetCopyNo()
410// << G4endl;
411 if (!(
412 iVAMNameCopyNo->GetName() ==
413 iPVNodeId->GetPhysicalVolume()->GetName() &&
414 iVAMNameCopyNo->GetCopyNo() ==
415 iPVNodeId->GetPhysicalVolume()->GetCopyNo()
416 )) {
417 break;
418 }
419 }
420 if (iVAMNameCopyNo == vamPath.end()) {
421// G4cout << "Match found" << G4endl;
422 if (!copyForVAM) {
423 pModifiedVisAtts = new G4VisAttributes(*pUnmodifiedVisAtts);
424 pVisAttribs = pModifiedVisAtts;
425 copyForVAM = true;
426 }
427 const G4VisAttributes& transVisAtts = iModifier->GetVisAttributes();
428 switch (iModifier->GetVisAttributesSignifier()) {
430 pModifiedVisAtts->SetVisibility(transVisAtts.IsVisible());
431 break;
433 pModifiedVisAtts->SetDaughtersInvisible
434 (transVisAtts.IsDaughtersInvisible());
435 break;
437 pModifiedVisAtts->SetColour(transVisAtts.GetColour());
438 break;
440 pModifiedVisAtts->SetLineStyle(transVisAtts.GetLineStyle());
441 break;
443 pModifiedVisAtts->SetLineWidth(transVisAtts.GetLineWidth());
444 break;
446 if (transVisAtts.GetForcedDrawingStyle() ==
448 pModifiedVisAtts->SetForceWireframe
449 (transVisAtts.IsForceDrawingStyle());
450 }
451 break;
453 if (transVisAtts.GetForcedDrawingStyle() ==
455 pModifiedVisAtts->SetForceSolid
456 (transVisAtts.IsForceDrawingStyle());
457 }
458 break;
460 pModifiedVisAtts->SetForceAuxEdgeVisible
461 (transVisAtts.IsForceAuxEdgeVisible());
462 break;
464 pModifiedVisAtts->SetForceLineSegmentsPerCircle
465 (transVisAtts.GetForcedLineSegmentsPerCircle());
466 break;
467 }
468 }
469 }
470 }
471
472 // There are various reasons why this volume
473 // might not be drawn...
474 G4bool culling = fpMP->IsCulling();
475 G4bool cullingInvisible = fpMP->IsCullingInvisible();
476 G4bool markedVisible = pVisAttribs->IsVisible();
477 G4bool cullingLowDensity = fpMP->IsDensityCulling();
478 G4double density = pMaterial? pMaterial->GetDensity(): 0;
479 G4double densityCut = fpMP -> GetVisibleDensity ();
480
481 // 1) Global culling is on....
482 if (culling) {
483 // 2) Culling of invisible volumes is on...
484 if (cullingInvisible) {
485 // 3) ...and the volume is marked not visible...
486 if (!markedVisible) thisToBeDrawn = false;
487 }
488 // 4) Or culling of low density volumes is on...
489 if (cullingLowDensity) {
490 // 5) ...and density is less than cut value...
491 if (density < densityCut) thisToBeDrawn = false;
492 }
493 }
494
495 // Record thisToBeDrawn in path...
496 fFullPVPath.back().SetDrawn(thisToBeDrawn);
497
498 if (thisToBeDrawn) {
499
500 // Update path of drawn physical volumes...
501 fDrawnPVPath.push_back
502 (G4PhysicalVolumeNodeID
503 (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform,thisToBeDrawn));
504
505 if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
506 // For top-level drawn volumes, explode along radius...
508 G4Transform3D centred = centering.inverse() * theNewAT;
509 G4Scale3D oldScale;
510 G4Rotate3D oldRotation;
511 G4Translate3D oldTranslation;
512 centred.getDecomposition(oldScale, oldRotation, oldTranslation);
513 G4double explodeFactor = fpMP->GetExplodeFactor();
514 G4Translate3D newTranslation =
515 G4Translate3D(explodeFactor * oldTranslation.dx(),
516 explodeFactor * oldTranslation.dy(),
517 explodeFactor * oldTranslation.dz());
518 theNewAT = centering * newTranslation * oldRotation * oldScale;
519 }
520
521 DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
522
523 }
524
525 // Make decision to draw daughters, if any. There are various
526 // reasons why daughters might not be drawn...
527
528 // First, reasons that do not depend on culling policy...
529 G4int nDaughters = pLV->GetNoDaughters();
530 G4bool daughtersToBeDrawn = true;
531 // 1) There are no daughters...
532 if (!nDaughters) daughtersToBeDrawn = false;
533 // 2) We are at the limit if requested depth...
534 else if (requestedDepth == 0) daughtersToBeDrawn = false;
535 // 3) The user has asked that the descent be curtailed...
536 else if (fCurtailDescent) daughtersToBeDrawn = false;
537
538 // Now, reasons that depend on culling policy...
539 else {
540 G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
541 // Culling of covered daughters request. This is computed in
542 // G4VSceneHandler::CreateModelingParameters() depending on view
543 // parameters...
544 G4bool cullingCovered = fpMP->IsCullingCovered();
545 G4bool surfaceDrawing =
548 if (pVisAttribs->IsForceDrawingStyle()) {
549 switch (pVisAttribs->GetForcedDrawingStyle()) {
550 default:
551 case G4VisAttributes::wireframe: surfaceDrawing = false; break;
552 case G4VisAttributes::solid: surfaceDrawing = true; break;
553 }
554 }
555 G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
556 // 4) Global culling is on....
557 if (culling) {
558 // 5) ..and culling of invisible volumes is on...
559 if (cullingInvisible) {
560 // 6) ...and the mother requests daughters invisible
561 if (daughtersInvisible) daughtersToBeDrawn = false;
562 }
563 // 7) Or culling of covered daughters is requested...
564 if (cullingCovered) {
565 // 8) ...and surface drawing is operating...
566 if (surfaceDrawing) {
567 // 9) ...but only if mother is visible...
568 if (thisToBeDrawn) {
569 // 10) ...and opaque...
570 if (opaque) daughtersToBeDrawn = false;
571 }
572 }
573 }
574 }
575 }
576
577 // Delete modified vis atts if created...
578 if (copyForVAM) {
579 delete pModifiedVisAtts;
580 pVisAttribs = pUnmodifiedVisAtts;
581 }
582
583 // Vis atts for this volume no longer needed if created...
584 if (visAttsCreated) delete pVisAttribs;
585
586 if (daughtersToBeDrawn) {
587 for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
588 // Store daughter pVPV in local variable ready for recursion...
589 G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
590 // Descend the geometry structure recursively...
593 (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
595 }
596 }
597
598 // Reset for normal descending of next volume at this level...
599 fCurtailDescent = false;
600
601 // Pop item from paths physical volumes...
602 fFullPVPath.pop_back();
603 if (thisToBeDrawn) {
604 fDrawnPVPath.pop_back();
605 }
606}
HepGeom::Translate3D G4Translate3D
bool G4bool
Definition: G4Types.hh:67
G4double GetAlpha() const
Definition: G4Colour.hh:141
const G4VisAttributes * GetVisAttributes() const
G4int GetNoDaughters() const
const G4VisAttributes * GetDefaultVisAttributes() const
const std::vector< VisAttributesModifier > & GetVisAttributesModifiers() const
const G4Point3D & GetExplodeCentre() const
PVNameCopyNoPath::const_iterator PVNameCopyNoPathConstIterator
std::vector< PVNameCopyNo > PVNameCopyNoPath
G4bool IsCullingInvisible() const
G4bool IsExplode() const
G4bool IsCulling() const
G4bool IsDensityCulling() const
DrawingStyle GetDrawingStyle() const
G4double GetExplodeFactor() const
G4bool IsCullingCovered() const
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
virtual void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
virtual G4int GetCopyNo() const =0
void SetForceAuxEdgeVisible(G4bool)
G4double GetLineWidth() const
void SetVisibility(G4bool)
G4bool IsDaughtersInvisible() const
void SetColour(const G4Colour &)
void SetDaughtersInvisible(G4bool)
void SetForceWireframe(G4bool)
G4int GetForcedLineSegmentsPerCircle() const
void SetLineWidth(G4double)
LineStyle GetLineStyle() const
const G4Colour & GetColour() const
void SetForceSolid(G4bool)
G4bool IsVisible() const
G4bool IsForceAuxEdgeVisible() const
ForcedDrawingStyle GetForcedDrawingStyle() const
void SetLineStyle(LineStyle)
void SetForceLineSegmentsPerCircle(G4int nSegments)
G4bool IsForceDrawingStyle() const
double dy() const
Definition: Transform3D.h:282
double dz() const
Definition: Transform3D.h:285
double dx() const
Definition: Transform3D.h:279
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:174

Referenced by VisitGeometryAndGetVisReps().

◆ DescribeSolid()

void G4PhysicalVolumeModel::DescribeSolid ( const G4Transform3D theAT,
G4VSolid pSol,
const G4VisAttributes pVisAttribs,
G4VGraphicsScene sceneHandler 
)
protectedvirtual

Reimplemented in G4LogicalVolumeModel.

Definition at line 608 of file G4PhysicalVolumeModel.cc.

613{
614 sceneHandler.PreAddSolid (theAT, *pVisAttribs);
615
616 G4VSolid* pSectionSolid = fpMP->GetSectionSolid();
617 G4VSolid* pCutawaySolid = fpMP->GetCutawaySolid();
618
619 if (!fpClippingSolid && !pSectionSolid && !pCutawaySolid) {
620
621 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
622
623 } else {
624
625 // Clipping, etc., performed by Boolean operations.
626
627 // First, get polyhedron for current solid...
628 if (pVisAttribs->IsForceLineSegmentsPerCircle())
630 (pVisAttribs->GetForcedLineSegmentsPerCircle());
631 else
633 const G4Polyhedron* pOriginal = pSol->GetPolyhedron();
635
636 if (!pOriginal) {
637
638 if (fpMP->IsWarning())
639 G4cout <<
640 "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \""
641 << pSol->GetName() <<
642 "\" has no polyhedron. Cannot by clipped."
643 << G4endl;
644 pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
645
646 } else {
647
648 G4Polyhedron resultant(*pOriginal);
649 G4VisAttributes resultantVisAttribs(*pVisAttribs);
650 G4VSolid* resultantSolid = 0;
651
652 if (fpClippingSolid) {
653 switch (fClippingMode) {
654 default:
655 case subtraction:
656 resultantSolid = new G4SubtractionSolid
657 ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
658 break;
659 case intersection:
660 resultantSolid = new G4IntersectionSolid
661 ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
662 break;
663 }
664 }
665
666 if (pSectionSolid) {
667 resultantSolid = new G4IntersectionSolid
668 ("sectioned_solid", pSol, pSectionSolid, theAT.inverse());
669 }
670
671 if (pCutawaySolid) {
672 resultantSolid = new G4SubtractionSolid
673 ("cutaway_solid", pSol, pCutawaySolid, theAT.inverse());
674 }
675
676 G4Polyhedron* tmpResultant = resultantSolid->GetPolyhedron();
677 if (tmpResultant) resultant = *tmpResultant;
678 else {
679 if (fpMP->IsWarning())
680 G4cout <<
681 "WARNING: G4PhysicalVolumeModel::DescribeSolid: resultant polyhedron for"
682 "\n solid \"" << pSol->GetName() <<
683 "\" not defined due to error during Boolean processing."
684 "\n Original will be drawn in red."
685 << G4endl;
686 resultantVisAttribs.SetColour(G4Colour::Red());
687 }
688
689 delete resultantSolid;
690
691 // Finally, force polyhedron drawing...
692 resultant.SetVisAttributes(resultantVisAttribs);
693 sceneHandler.BeginPrimitives(theAT);
694 sceneHandler.AddPrimitive(resultant);
695 sceneHandler.EndPrimitives();
696 }
697 }
698 sceneHandler.PostAddSolid ();
699}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4Colour Red()
Definition: G4Colour.hh:147
G4bool IsWarning() const
G4VSolid * GetCutawaySolid() const
G4int GetNoOfSides() const
G4VSolid * GetSectionSolid() const
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void PostAddSolid()=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void EndPrimitives()=0
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:647
G4bool IsForceLineSegmentsPerCircle() const
static void SetNumberOfRotationSteps(G4int n)
static void ResetNumberOfRotationSteps()

Referenced by DescribeAndDescend().

◆ DescribeYourselfTo()

void G4PhysicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene sceneHandler)
virtual

Implements G4VModel.

Definition at line 135 of file G4PhysicalVolumeModel.cc.

137{
138 if (!fpTopPV) G4Exception
139 ("G4PhysicalVolumeModel::DescribeYourselfTo",
140 "modeling0012", FatalException, "No model.");
141
142 if (!fpMP) G4Exception
143 ("G4PhysicalVolumeModel::DescribeYourselfTo",
144 "modeling0003", FatalException, "No modeling parameters.");
145
146 // For safety...
147 fCurrentDepth = 0;
148
149 G4Transform3D startingTransformation = fTransform;
150
152 (fpTopPV,
154 startingTransformation,
155 sceneHandler);
156
157 // Clear data...
158 fCurrentDepth = 0;
159 fpCurrentPV = 0;
160 fpCurrentLV = 0;
162 if (fFullPVPath.size() != fBaseFullPVPath.size()) {
163 // They should be equal if pushing and popping is happening properly.
165 ("G4PhysicalVolumeModel::DescribeYourselfTo",
166 "modeling0013",
168 "Path at start of modeling not equal to base path. Something badly"
169 "\nwrong. Please contact visualisation coordinator.");
170 }
171 fDrawnPVPath.clear();
172}
@ FatalException
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath

Referenced by CalculateExtent(), DescribeSolid(), G4LogicalVolumeModel::DescribeYourselfTo(), G4ASCIITreeSceneHandler::EndModeling(), G4VisCommandSceneAddVolume::SetNewValue(), and Validate().

◆ GetAttDefs()

const std::map< G4String, G4AttDef > * G4PhysicalVolumeModel::GetAttDefs ( ) const

Definition at line 754 of file G4PhysicalVolumeModel.cc.

755{
756 G4bool isNew;
757 std::map<G4String,G4AttDef>* store
758 = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
759 if (isNew) {
760 (*store)["PVPath"] =
761 G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
762 (*store)["LVol"] =
763 G4AttDef("LVol","Logical Volume","Physics","","G4String");
764 (*store)["Solid"] =
765 G4AttDef("Solid","Solid Name","Physics","","G4String");
766 (*store)["EType"] =
767 G4AttDef("EType","Entity Type","Physics","","G4String");
768 (*store)["DmpSol"] =
769 G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
770 (*store)["Trans"] =
771 G4AttDef("Trans","Transformation of volume","Physics","","G4String");
772 (*store)["Material"] =
773 G4AttDef("Material","Material Name","Physics","","G4String");
774 (*store)["Density"] =
775 G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
776 (*store)["State"] =
777 G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
778 (*store)["Radlen"] =
779 G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
780 (*store)["Region"] =
781 G4AttDef("Region","Cuts Region","Physics","","G4String");
782 (*store)["RootRegion"] =
783 G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
784 }
785 return store;
786}
std::map< G4String, G4AttDef > * GetInstance(G4String storeKey, G4bool &isNew)

Referenced by G4VSceneHandler::LoadAtts(), G4XXXStoredSceneHandler::PreAddSolid(), and G4VisCommandList::SetNewValue().

◆ GetClippingSolid()

const G4VSolid * G4PhysicalVolumeModel::GetClippingSolid ( ) const
inline

Definition at line 158 of file G4PhysicalVolumeModel.hh.

159 {return fpClippingSolid;}

◆ GetCurrentDepth()

G4int G4PhysicalVolumeModel::GetCurrentDepth ( ) const
inline

◆ GetCurrentDescription()

G4String G4PhysicalVolumeModel::GetCurrentDescription ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 186 of file G4PhysicalVolumeModel.cc.

187{
188 return "G4PhysicalVolumeModel " + GetCurrentTag ();
189}

◆ GetCurrentLV()

G4LogicalVolume * G4PhysicalVolumeModel::GetCurrentLV ( ) const
inline

◆ GetCurrentMaterial()

G4Material * G4PhysicalVolumeModel::GetCurrentMaterial ( ) const
inline

◆ GetCurrentPV()

G4VPhysicalVolume * G4PhysicalVolumeModel::GetCurrentPV ( ) const
inline

◆ GetCurrentTag()

G4String G4PhysicalVolumeModel::GetCurrentTag ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 174 of file G4PhysicalVolumeModel.cc.

175{
176 if (fpCurrentPV) {
177 std::ostringstream o;
178 o << fpCurrentPV -> GetCopyNo ();
179 return fpCurrentPV -> GetName () + "." + o.str();
180 }
181 else {
182 return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
183 }
184}

Referenced by GetCurrentDescription().

◆ GetDrawnPVPath()

◆ GetFullPVPath()

const std::vector< G4PhysicalVolumeNodeID > & G4PhysicalVolumeModel::GetFullPVPath ( ) const
inline

Definition at line 173 of file G4PhysicalVolumeModel.hh.

174 {return fFullPVPath;}

◆ GetRequestedDepth()

G4int G4PhysicalVolumeModel::GetRequestedDepth ( ) const
inline

Definition at line 156 of file G4PhysicalVolumeModel.hh.

156{return fRequestedDepth;}

Referenced by G4ASCIITreeSceneHandler::EndModeling().

◆ GetTopPhysicalVolume()

G4VPhysicalVolume * G4PhysicalVolumeModel::GetTopPhysicalVolume ( ) const
inline

◆ SetBaseFullPVPath()

void G4PhysicalVolumeModel::SetBaseFullPVPath ( const std::vector< G4PhysicalVolumeNodeID > &  baseFullPVPath)
inline

Definition at line 200 of file G4PhysicalVolumeModel.hh.

202 {
203 fBaseFullPVPath = baseFullPVPath;
204 fFullPVPath = baseFullPVPath;
205 }

Referenced by G4VisCommandSceneAddVolume::SetNewValue().

◆ SetClippingMode()

void G4PhysicalVolumeModel::SetClippingMode ( ClippingMode  mode)
inline

Definition at line 215 of file G4PhysicalVolumeModel.hh.

215 {
216 fClippingMode = mode;
217 }

◆ SetClippingSolid()

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid pClippingSolid)
inline

Definition at line 211 of file G4PhysicalVolumeModel.hh.

211 {
212 fpClippingSolid = pClippingSolid;
213 }

◆ SetRequestedDepth()

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int  requestedDepth)
inline

Definition at line 207 of file G4PhysicalVolumeModel.hh.

207 {
208 fRequestedDepth = requestedDepth;
209 }

◆ Validate()

G4bool G4PhysicalVolumeModel::Validate ( G4bool  warn)
virtual

Reimplemented from G4VModel.

Definition at line 701 of file G4PhysicalVolumeModel.cc.

702{
703 G4TransportationManager* transportationManager =
705
706 size_t nWorlds = transportationManager->GetNoWorlds();
707
708 G4bool found = false;
709
710 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
711 transportationManager->GetWorldsIterator();
712 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
713 G4VPhysicalVolume* world = (*iterWorld);
714 // The idea now is to seek a PV with the same name and copy no
715 // in the hope it's the same one!!
716 G4PhysicalVolumeModel searchModel (world);
718 (&searchModel, fTopPVName, fTopPVCopyNo);
719 G4ModelingParameters mp; // Default modeling parameters for this search.
721 searchModel.SetModelingParameters (&mp);
722 searchModel.DescribeYourselfTo (searchScene);
723 G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
724 if (foundVolume) {
725 if (foundVolume != fpTopPV && warn) {
726 G4cout <<
727 "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
728 "\n copy number (\""
729 << fTopPVName << "\", copy " << fTopPVCopyNo
730 << ") still exists and is being used."
731 "\n But it is not the same volume you originally specified"
732 "\n in /vis/scene/add/."
733 << G4endl;
734 }
735 fpTopPV = foundVolume;
737 found = true;
738 }
739 }
740 if (found) return true;
741 else {
742 if (warn) {
743 G4cout <<
744 "G4PhysicalVolumeModel::Validate(): No volume of name and"
745 "\n copy number (\""
746 << fTopPVName << "\", copy " << fTopPVCopyNo
747 << ") exists."
748 << G4endl;
749 }
750 return false;
751 }
752}
void SetDefaultVisAttributes(const G4VisAttributes *pDefaultVisAttributes)
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const

◆ VisitGeometryAndGetVisReps()

void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps ( G4VPhysicalVolume pVPV,
G4int  requestedDepth,
const G4Transform3D theAT,
G4VGraphicsScene sceneHandler 
)
protected

Definition at line 191 of file G4PhysicalVolumeModel.cc.

196{
197 // Visits geometry structure to a given depth (requestedDepth), starting
198 // at given physical volume with given starting transformation and
199 // describes volumes to the scene handler.
200 // requestedDepth < 0 (default) implies full visit.
201 // theAT is the Accumulated Transformation.
202
203 // Find corresponding logical volume and (later) solid, storing in
204 // local variables to preserve re-entrancy.
205 G4LogicalVolume* pLV = pVPV -> GetLogicalVolume ();
206
207 G4VSolid* pSol;
208 G4Material* pMaterial;
209
210 if (!(pVPV -> IsReplicated ())) {
211 // Non-replicated physical volume.
212 pSol = pLV -> GetSolid ();
213 pMaterial = pLV -> GetMaterial ();
214 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
215 theAT, sceneHandler);
216 }
217 else {
218 // Replicated or parametrised physical volume.
219 EAxis axis;
220 G4int nReplicas;
221 G4double width;
222 G4double offset;
223 G4bool consuming;
224 pVPV -> GetReplicationData (axis, nReplicas, width, offset, consuming);
225 if (fCurrentDepth == 0) nReplicas = 1; // Just draw first
226 G4VPVParameterisation* pP = pVPV -> GetParameterisation ();
227 if (pP) { // Parametrised volume.
228 for (int n = 0; n < nReplicas; n++) {
229 pSol = pP -> ComputeSolid (n, pVPV);
230 pP -> ComputeTransformation (n, pVPV);
231 pSol -> ComputeDimensions (pP, n, pVPV);
232 pVPV -> SetCopyNo (n);
233 // Create a touchable of current parent for ComputeMaterial.
234 // fFullPVPath has not been updated yet so at this point it
235 // corresponds to the parent.
236 G4PhysicalVolumeModelTouchable parentTouchable(fFullPVPath);
237 pMaterial = pP -> ComputeMaterial (n, pVPV, &parentTouchable);
238 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
239 theAT, sceneHandler);
240 }
241 }
242 else { // Plain replicated volume. From geometry_guide.txt...
243 // The replica's positions are claculated by means of a linear formula.
244 // Replication may occur along:
245 //
246 // o Cartesian axes (kXAxis,kYAxis,kZAxis)
247 //
248 // The replications, of specified width have coordinates of
249 // form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1
250 // for the case of kXAxis, and are unrotated.
251 //
252 // o Radial axis (cylindrical polar) (kRho)
253 //
254 // The replications are cons/tubs sections, centred on the origin
255 // and are unrotated.
256 // They have radii of width*n+offset to width*(n+1)+offset
257 // where n=0..nReplicas-1
258 //
259 // o Phi axis (cylindrical polar) (kPhi)
260 // The replications are `phi sections' or wedges, and of cons/tubs form
261 // They have phi of offset+n*width to offset+(n+1)*width where
262 // n=0..nReplicas-1
263 //
264 pSol = pLV -> GetSolid ();
265 pMaterial = pLV -> GetMaterial ();
266 G4ThreeVector originalTranslation = pVPV -> GetTranslation ();
267 G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation ();
268 G4double originalRMin = 0., originalRMax = 0.;
269 if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
270 originalRMin = ((G4Tubs*)pSol)->GetInnerRadius();
271 originalRMax = ((G4Tubs*)pSol)->GetOuterRadius();
272 }
273 G4bool visualisable = true;
274 for (int n = 0; n < nReplicas; n++) {
275 G4ThreeVector translation; // Null.
276 G4RotationMatrix rotation; // Null - life long enough for visualizing.
277 G4RotationMatrix* pRotation = 0;
278 switch (axis) {
279 default:
280 case kXAxis:
281 translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0);
282 break;
283 case kYAxis:
284 translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0);
285 break;
286 case kZAxis:
287 translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width);
288 break;
289 case kRho:
290 if (pSol->GetEntityType() == "G4Tubs") {
291 ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset);
292 ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset);
293 } else {
294 if (fpMP->IsWarning())
295 G4cout <<
296 "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:"
297 "\n built-in replicated volumes replicated in radius for "
298 << pSol->GetEntityType() <<
299 "-type\n solids (your solid \""
300 << pSol->GetName() <<
301 "\") are not visualisable."
302 << G4endl;
303 visualisable = false;
304 }
305 break;
306 case kPhi:
307 rotation.rotateZ (-(offset+(n+0.5)*width));
308 // Minus Sign because for the physical volume we need the
309 // coordinate system rotation.
310 pRotation = &rotation;
311 break;
312 }
313 pVPV -> SetTranslation (translation);
314 pVPV -> SetRotation (pRotation);
315 pVPV -> SetCopyNo (n);
316 if (visualisable) {
317 DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
318 theAT, sceneHandler);
319 }
320 }
321 // Restore originals...
322 pVPV -> SetTranslation (originalTranslation);
323 pVPV -> SetRotation (pOriginalRotation);
324 if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
325 ((G4Tubs*)pSol)->SetInnerRadius(originalRMin);
326 ((G4Tubs*)pSol)->SetOuterRadius(originalRMax);
327 }
328 }
329 }
330
331 return;
332}
CLHEP::Hep3Vector G4ThreeVector
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
void DescribeAndDescend(G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
Definition: G4Tubs.hh:77
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
@ kRho
Definition: geomdefs.hh:54

Referenced by DescribeAndDescend(), and DescribeYourselfTo().

Member Data Documentation

◆ fBaseFullPVPath

std::vector<G4PhysicalVolumeNodeID> G4PhysicalVolumeModel::fBaseFullPVPath
protected

Definition at line 260 of file G4PhysicalVolumeModel.hh.

Referenced by DescribeYourselfTo(), and SetBaseFullPVPath().

◆ fClippingMode

ClippingMode G4PhysicalVolumeModel::fClippingMode
protected

Definition at line 265 of file G4PhysicalVolumeModel.hh.

Referenced by DescribeSolid(), and SetClippingMode().

◆ fCurrentDepth

G4int G4PhysicalVolumeModel::fCurrentDepth
protected

◆ fCurtailDescent

G4bool G4PhysicalVolumeModel::fCurtailDescent
protected

Definition at line 263 of file G4PhysicalVolumeModel.hh.

Referenced by CurtailDescent(), and DescribeAndDescend().

◆ fDrawnPVPath

std::vector<G4PhysicalVolumeNodeID> G4PhysicalVolumeModel::fDrawnPVPath
protected

◆ fFullPVPath

◆ fpClippingSolid

G4VSolid* G4PhysicalVolumeModel::fpClippingSolid
protected

◆ fpCurrentLV

G4LogicalVolume* G4PhysicalVolumeModel::fpCurrentLV
protected

◆ fpCurrentMaterial

G4Material* G4PhysicalVolumeModel::fpCurrentMaterial
protected

◆ fpCurrentPV

G4VPhysicalVolume* G4PhysicalVolumeModel::fpCurrentPV
protected

◆ fpCurrentTransform

G4Transform3D* G4PhysicalVolumeModel::fpCurrentTransform
protected

Definition at line 259 of file G4PhysicalVolumeModel.hh.

Referenced by CreateCurrentAttValues(), and DescribeAndDescend().

◆ fpTopPV

◆ fRequestedDepth

G4int G4PhysicalVolumeModel::fRequestedDepth
protected

◆ fTopPVCopyNo

G4int G4PhysicalVolumeModel::fTopPVCopyNo
protected

Definition at line 251 of file G4PhysicalVolumeModel.hh.

Referenced by G4PhysicalVolumeModel(), and Validate().

◆ fTopPVName

G4String G4PhysicalVolumeModel::fTopPVName
protected

Definition at line 250 of file G4PhysicalVolumeModel.hh.

Referenced by G4PhysicalVolumeModel(), and Validate().

◆ fUseFullExtent

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 254 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent().


The documentation for this class was generated from the following files: