Geant4  10.02.p03
G4PhysicalVolumeModel Class Reference

#include <G4PhysicalVolumeModel.hh>

Inheritance diagram for G4PhysicalVolumeModel:
Collaboration 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
 
G4Transform3DGetCurrentTransform () 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 Abort () const
 
void CurtailDescent () const
 
- Public Member Functions inherited from G4VModel
 G4VModel (const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
 
virtual ~G4VModel ()
 
const G4ModelingParametersGetModelingParameters () const
 
const G4StringGetType () 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 &)
 

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 fAbort
 
G4bool fCurtailDescent
 
G4VSolidfpClippingSolid
 
ClippingMode fClippingMode
 
- Protected Attributes inherited from G4VModel
G4String fType
 
G4String fGlobalTag
 
G4String fGlobalDescription
 
G4VisExtent fExtent
 
G4Transform3D fTransform
 
const G4ModelingParametersfpMP
 

Private Member Functions

 G4PhysicalVolumeModel (const G4PhysicalVolumeModel &)
 
G4PhysicalVolumeModeloperator= (const G4PhysicalVolumeModel &)
 

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() [1/2]

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

Definition at line 59 of file G4PhysicalVolumeModel.cc.

64 : G4VModel (modelTransformation, pMP)
65 , fpTopPV (pVPV)
66 , fTopPVCopyNo (0)
67 , fRequestedDepth (requestedDepth)
68 , fUseFullExtent (useFullExtent)
69 , fCurrentDepth (0)
70 , fpCurrentPV (0)
71 , fpCurrentLV (0)
74 , fAbort (false)
75 , fCurtailDescent (false)
76 , fpClippingSolid (0)
78 {
79  fType = "G4PhysicalVolumeModel";
80 
81  if (!fpTopPV) {
82 
83  // In some circumstances creating an "empty" G4PhysicalVolumeModel is
84  // allowed, so I have supressed the G4Exception below. If it proves to
85  // be a problem we might have to re-instate it, but it is unlikley to
86  // be used except by visualisation experts. See, for example, /vis/list,
87  // where it is used simply to get a list of G4AttDefs.
88  // G4Exception
89  // ("G4PhysicalVolumeModel::G4PhysicalVolumeModel",
90  // "modeling0010", FatalException, "Null G4PhysicalVolumeModel pointer.");
91 
92  fTopPVName = "NULL";
93  fGlobalTag = "Empty";
94  fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
95 
96  } else {
97 
98  fTopPVName = fpTopPV -> GetName ();
99  fTopPVCopyNo = fpTopPV -> GetCopyNo ();
100  std::ostringstream o;
101  o << fpTopPV -> GetCopyNo ();
102  fGlobalTag = fpTopPV -> GetName () + "." + o.str();
103  fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
104 
108  fpCurrentTransform = const_cast<G4Transform3D*>(&modelTransformation);
109 
110  CalculateExtent ();
111  }
112 }
G4String fType
Definition: G4VModel.hh:108
G4Material * GetMaterial() const
G4Transform3D * fpCurrentTransform
G4VPhysicalVolume * fpCurrentPV
G4String fGlobalTag
Definition: G4VModel.hh:109
G4VModel(const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0)
Definition: G4VModel.cc:38
G4String fGlobalDescription
Definition: G4VModel.hh:110
G4VPhysicalVolume * fpTopPV
G4LogicalVolume * GetLogicalVolume() const

◆ ~G4PhysicalVolumeModel()

G4PhysicalVolumeModel::~G4PhysicalVolumeModel ( )
virtual

Definition at line 114 of file G4PhysicalVolumeModel.cc.

115 {
116  delete fpClippingSolid;
117 }

◆ G4PhysicalVolumeModel() [2/2]

G4PhysicalVolumeModel::G4PhysicalVolumeModel ( const G4PhysicalVolumeModel )
private

Member Function Documentation

◆ Abort()

void G4PhysicalVolumeModel::Abort ( ) const
inline

Definition at line 225 of file G4PhysicalVolumeModel.hh.

225 {fAbort = true;}
Here is the caller graph for this function:

◆ CalculateExtent()

void G4PhysicalVolumeModel::CalculateExtent ( )
protected

Definition at line 119 of file G4PhysicalVolumeModel.cc.

120 {
121  if (fUseFullExtent) {
122  fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
123  }
124  else {
125  G4BoundingSphereScene bsScene(this);
126  const G4int tempRequestedDepth = fRequestedDepth;
127  fRequestedDepth = -1; // Always search to all depths to define extent.
128  const G4ModelingParameters* tempMP = fpMP;
129  G4ModelingParameters mParams
130  (0, // No default vis attributes needed.
131  G4ModelingParameters::wf, // wireframe (not relevant for this).
132  true, // Global culling.
133  true, // Cull invisible volumes.
134  false, // Density culling.
135  0., // Density (not relevant if density culling false).
136  true, // Cull daughters of opaque mothers.
137  24); // No of sides (not relevant for this operation).
138  fpMP = &mParams;
139  DescribeYourselfTo (bsScene);
140  G4double radius = bsScene.GetRadius();
141  if (radius < 0.) { // Nothing in the scene.
142  fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
143  } else {
144  // Transform back to coordinates relative to the top
145  // transformation, which is in G4VModel::fTransform. This makes
146  // it conform to all models, which are defined by a
147  // transformation and an extent relative to that
148  // transformation...
149  G4Point3D centre = bsScene.GetCentre();
150  centre.transform(fTransform.inverse());
151  fExtent = G4VisExtent(centre, radius);
152  }
153  fpMP = tempMP;
154  fRequestedDepth = tempRequestedDepth;
155  }
156 }
G4Transform3D fTransform
Definition: G4VModel.hh:112
const G4VisExtent & GetExtent() const
int G4int
Definition: G4Types.hh:78
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
G4VPhysicalVolume * fpTopPV
G4VisExtent fExtent
Definition: G4VModel.hh:111
double G4double
Definition: G4Types.hh:76
Transform3D inverse() const
Definition: Transform3D.cc:142
void DescribeYourselfTo(G4VGraphicsScene &)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateCurrentAttValues()

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

Definition at line 860 of file G4PhysicalVolumeModel.cc.

861 {
862  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
863  std::ostringstream oss;
864  for (size_t i = 0; i < fFullPVPath.size(); ++i) {
865  oss << fFullPVPath[i].GetPhysicalVolume()->GetName()
866  << ':' << fFullPVPath[i].GetCopyNo();
867  if (i != fFullPVPath.size() - 1) oss << '/';
868  }
869 
870  if (!fpCurrentLV) {
872  ("G4PhysicalVolumeModel::CreateCurrentAttValues",
873  "modeling0004",
874  JustWarning,
875  "Current logical volume not defined.");
876  return values;
877  }
878 
879  values->push_back(G4AttValue("PVPath", oss.str(),""));
880  values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
881  G4VSolid* pSol = fpCurrentLV->GetSolid();
882  values->push_back(G4AttValue("Solid", pSol->GetName(),""));
883  values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
884  oss.str(""); oss << '\n' << *pSol;
885  values->push_back(G4AttValue("DmpSol", oss.str(),""));
886  const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
887  const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
888  oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
889  values->push_back(G4AttValue("LocalTrans", oss.str(),""));
890  oss.str(""); oss << '\n' << *fpCurrentTransform;
891  values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
892  G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
893  values->push_back(G4AttValue("Material", matName,""));
895  values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
897  oss.str(""); oss << matState;
898  values->push_back(G4AttValue("State", oss.str(),""));
900  values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
901  G4Region* region = fpCurrentLV->GetRegion();
902  G4String regionName = region? region->GetName(): G4String("No region");
903  values->push_back(G4AttValue("Region", regionName,""));
904  oss.str(""); oss << fpCurrentLV->IsRootRegion();
905  values->push_back(G4AttValue("RootRegion", oss.str(),""));
906  return values;
907 }
G4State GetState() const
Definition: G4Material.hh:181
G4State
Definition: G4Material.hh:114
G4double GetDensity() const
Definition: G4Material.hh:180
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4Transform3D * fpCurrentTransform
virtual G4GeometryType GetEntityType() const =0
G4String GetName() const
G4VPhysicalVolume * fpCurrentPV
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
const G4String & GetName() const
G4bool IsRootRegion() const
HepGeom::Transform3D G4Transform3D
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
G4Region * GetRegion() const
G4RotationMatrix GetObjectRotationValue() const
double G4double
Definition: G4Types.hh:76
const G4ThreeVector & GetTranslation() const
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetRadlen() const
Definition: G4Material.hh:220
G4VSolid * GetSolid() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CurtailDescent()

void G4PhysicalVolumeModel::CurtailDescent ( ) const
inline

Definition at line 228 of file G4PhysicalVolumeModel.hh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ DescribeAndDescend()

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

Definition at line 358 of file G4PhysicalVolumeModel.cc.

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

◆ DescribeSolid()

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

Reimplemented in G4LogicalVolumeModel.

Definition at line 637 of file G4PhysicalVolumeModel.cc.

641 {
642  sceneHandler.PreAddSolid (theAT, *pVisAttribs);
643 
644  G4VSolid* pSectionSolid = fpMP->GetSectionSolid();
645  G4VSolid* pCutawaySolid = fpMP->GetCutawaySolid();
646 
647  if (!fpClippingSolid && !pSectionSolid && !pCutawaySolid) {
648 
649  pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
650 
651  } else {
652 
653  // Clipping, etc., performed by Boolean operations.
654 
655  // First, get polyhedron for current solid...
656  if (pVisAttribs->IsForceLineSegmentsPerCircle())
658  (pVisAttribs->GetForcedLineSegmentsPerCircle());
659  else
661  const G4Polyhedron* pOriginal = pSol->GetPolyhedron();
663 
664  if (!pOriginal) {
665 
666  if (fpMP->IsWarning())
667  G4cout <<
668  "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \""
669  << pSol->GetName() <<
670  "\" has no polyhedron. Cannot by clipped."
671  << G4endl;
672  pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment.
673 
674  } else {
675 
676  G4Polyhedron resultant(*pOriginal);
677  G4VisAttributes resultantVisAttribs(*pVisAttribs);
678  G4VSolid* resultantSolid = 0;
679 
680  if (fpClippingSolid) {
681  switch (fClippingMode) {
682  default:
683  case subtraction:
684  resultantSolid = new G4SubtractionSolid
685  ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
686  break;
687  case intersection:
688  resultantSolid = new G4IntersectionSolid
689  ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
690  break;
691  }
692  }
693 
694  if (pSectionSolid) {
695  resultantSolid = new G4IntersectionSolid
696  ("sectioned_solid", pSol, pSectionSolid, theAT.inverse());
697  }
698 
699  if (pCutawaySolid) {
700  resultantSolid = new G4SubtractionSolid
701  ("cutaway_solid", pSol, pCutawaySolid, theAT.inverse());
702  }
703 
704  G4Polyhedron* tmpResultant = resultantSolid->GetPolyhedron();
705  if (tmpResultant) resultant = *tmpResultant;
706  else {
707  if (fpMP->IsWarning())
708  G4cout <<
709  "WARNING: G4PhysicalVolumeModel::DescribeSolid: resultant polyhedron for"
710  "\n solid \"" << pSol->GetName() <<
711  "\" not defined due to error during Boolean processing."
712  "\n Original will be drawn in red."
713  << G4endl;
714  resultantVisAttribs.SetColour(G4Colour::Red());
715  }
716 
717  delete resultantSolid;
718 
719  // Finally, force polyhedron drawing...
720  resultant.SetVisAttributes(resultantVisAttribs);
721  sceneHandler.BeginPrimitives(theAT);
722  sceneHandler.AddPrimitive(resultant);
723  sceneHandler.EndPrimitives();
724  }
725  }
726  sceneHandler.PostAddSolid ();
727 }
virtual void PostAddSolid()=0
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
G4String GetName() const
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:644
G4int GetNoOfSides() const
virtual void AddPrimitive(const G4Polyline &)=0
G4VSolid * GetCutawaySolid() const
G4GLOB_DLL std::ostream G4cout
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &visAttribs)=0
G4bool IsForceLineSegmentsPerCircle() const
G4VSolid * GetSectionSolid() const
#define G4endl
Definition: G4ios.hh:61
static G4Colour Red()
Definition: G4Colour.hh:148
virtual void EndPrimitives()=0
static void SetNumberOfRotationSteps(G4int n)
G4bool IsWarning() const
Transform3D inverse() const
Definition: Transform3D.cc:142
G4int GetForcedLineSegmentsPerCircle() const
static void ResetNumberOfRotationSteps()
void DescribeYourselfTo(G4VGraphicsScene &)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DescribeYourselfTo()

void G4PhysicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene sceneHandler)
virtual

Implements G4VModel.

Definition at line 159 of file G4PhysicalVolumeModel.cc.

160 {
161  if (!fpTopPV) G4Exception
162  ("G4PhysicalVolumeModel::DescribeYourselfTo",
163  "modeling0012", FatalException, "No model.");
164 
165  if (!fpMP) G4Exception
166  ("G4PhysicalVolumeModel::DescribeYourselfTo",
167  "modeling0003", FatalException, "No modeling parameters.");
168 
169  // For safety...
170  fCurrentDepth = 0;
171 
172  G4Transform3D startingTransformation = fTransform;
173 
175  (fpTopPV,
177  startingTransformation,
178  sceneHandler);
179 
180  // Clear data...
181  fCurrentDepth = 0;
182  fpCurrentPV = 0;
183  fpCurrentLV = 0;
184  fpCurrentMaterial = 0;
185  if (fFullPVPath.size() != fBaseFullPVPath.size()) {
186  // They should be equal if pushing and popping is happening properly.
188  ("G4PhysicalVolumeModel::DescribeYourselfTo",
189  "modeling0013",
191  "Path at start of modeling not equal to base path. Something badly"
192  "\nwrong. Please contact visualisation coordinator.");
193  }
194  fDrawnPVPath.clear();
195  fAbort = false;
196  fCurtailDescent = false;
197 }
G4Transform3D fTransform
Definition: G4VModel.hh:112
G4VPhysicalVolume * fpCurrentPV
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4VPhysicalVolume * fpTopPV
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAttDefs()

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

Definition at line 783 of file G4PhysicalVolumeModel.cc.

784 {
785  G4bool isNew;
786  std::map<G4String,G4AttDef>* store
787  = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
788  if (isNew) {
789  (*store)["PVPath"] =
790  G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
791  (*store)["LVol"] =
792  G4AttDef("LVol","Logical Volume","Physics","","G4String");
793  (*store)["Solid"] =
794  G4AttDef("Solid","Solid Name","Physics","","G4String");
795  (*store)["EType"] =
796  G4AttDef("EType","Entity Type","Physics","","G4String");
797  (*store)["DmpSol"] =
798  G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
799  (*store)["LocalTrans"] =
800  G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
801  (*store)["GlobalTrans"] =
802  G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
803  (*store)["Material"] =
804  G4AttDef("Material","Material Name","Physics","","G4String");
805  (*store)["Density"] =
806  G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
807  (*store)["State"] =
808  G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
809  (*store)["Radlen"] =
810  G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
811  (*store)["Region"] =
812  G4AttDef("Region","Cuts Region","Physics","","G4String");
813  (*store)["RootRegion"] =
814  G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
815  }
816  return store;
817 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetClippingSolid()

const G4VSolid* G4PhysicalVolumeModel::GetClippingSolid ( ) const
inline

Definition at line 158 of file G4PhysicalVolumeModel.hh.

159  {return fpClippingSolid;}

◆ GetCurrentDepth()

G4int G4PhysicalVolumeModel::GetCurrentDepth ( ) const
inline

Definition at line 161 of file G4PhysicalVolumeModel.hh.

Here is the caller graph for this function:

◆ GetCurrentDescription()

G4String G4PhysicalVolumeModel::GetCurrentDescription ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 211 of file G4PhysicalVolumeModel.cc.

212 {
213  return "G4PhysicalVolumeModel " + GetCurrentTag ();
214 }
Here is the call graph for this function:

◆ GetCurrentLV()

G4LogicalVolume* G4PhysicalVolumeModel::GetCurrentLV ( ) const
inline

Definition at line 167 of file G4PhysicalVolumeModel.hh.

167 {return fpCurrentLV;}
Here is the caller graph for this function:

◆ GetCurrentMaterial()

G4Material* G4PhysicalVolumeModel::GetCurrentMaterial ( ) const
inline

Definition at line 170 of file G4PhysicalVolumeModel.hh.

170 {return fpCurrentMaterial;}
Here is the caller graph for this function:

◆ GetCurrentPV()

G4VPhysicalVolume* G4PhysicalVolumeModel::GetCurrentPV ( ) const
inline

Definition at line 164 of file G4PhysicalVolumeModel.hh.

164 {return fpCurrentPV;}
G4VPhysicalVolume * fpCurrentPV
Here is the caller graph for this function:

◆ GetCurrentTag()

G4String G4PhysicalVolumeModel::GetCurrentTag ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 199 of file G4PhysicalVolumeModel.cc.

200 {
201  if (fpCurrentPV) {
202  std::ostringstream o;
203  o << fpCurrentPV -> GetCopyNo ();
204  return fpCurrentPV -> GetName () + "." + o.str();
205  }
206  else {
207  return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
208  }
209 }
G4VPhysicalVolume * fpCurrentPV
G4String fGlobalTag
Definition: G4VModel.hh:109
Here is the caller graph for this function:

◆ GetCurrentTransform()

G4Transform3D* G4PhysicalVolumeModel::GetCurrentTransform ( ) const
inline

Definition at line 173 of file G4PhysicalVolumeModel.hh.

173 {return fpCurrentTransform;}
G4Transform3D * fpCurrentTransform
Here is the caller graph for this function:

◆ GetDrawnPVPath()

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

Definition at line 183 of file G4PhysicalVolumeModel.hh.

184  {return fDrawnPVPath;}
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetFullPVPath()

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

Definition at line 176 of file G4PhysicalVolumeModel.hh.

177  {return fFullPVPath;}
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
Here is the caller graph for this function:

◆ GetRequestedDepth()

G4int G4PhysicalVolumeModel::GetRequestedDepth ( ) const
inline

Definition at line 156 of file G4PhysicalVolumeModel.hh.

Here is the caller graph for this function:

◆ GetTopPhysicalVolume()

G4VPhysicalVolume* G4PhysicalVolumeModel::GetTopPhysicalVolume ( ) const
inline

Definition at line 154 of file G4PhysicalVolumeModel.hh.

154 {return fpTopPV;}
G4VPhysicalVolume * fpTopPV
Here is the caller graph for this function:

◆ operator=()

G4PhysicalVolumeModel& G4PhysicalVolumeModel::operator= ( const G4PhysicalVolumeModel )
private

◆ SetBaseFullPVPath()

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

Definition at line 204 of file G4PhysicalVolumeModel.hh.

205  {
206  fBaseFullPVPath = baseFullPVPath;
207  fFullPVPath = baseFullPVPath;
208  }
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
std::vector< G4PhysicalVolumeNodeID > fBaseFullPVPath
Here is the caller graph for this function:

◆ SetClippingMode()

void G4PhysicalVolumeModel::SetClippingMode ( ClippingMode  mode)
inline

Definition at line 218 of file G4PhysicalVolumeModel.hh.

218  {
219  fClippingMode = mode;
220  }
Here is the call graph for this function:

◆ SetClippingSolid()

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid pClippingSolid)
inline

Definition at line 214 of file G4PhysicalVolumeModel.hh.

214  {
215  fpClippingSolid = pClippingSolid;
216  }

◆ SetRequestedDepth()

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int  requestedDepth)
inline

Definition at line 210 of file G4PhysicalVolumeModel.hh.

210  {
211  fRequestedDepth = requestedDepth;
212  }

◆ Validate()

G4bool G4PhysicalVolumeModel::Validate ( G4bool  warn)
virtual

Reimplemented from G4VModel.

Definition at line 729 of file G4PhysicalVolumeModel.cc.

730 {
731  G4TransportationManager* transportationManager =
733 
734  size_t nWorlds = transportationManager->GetNoWorlds();
735 
736  G4bool found = false;
737 
738  std::vector<G4VPhysicalVolume*>::iterator iterWorld =
739  transportationManager->GetWorldsIterator();
740  for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
741  G4VPhysicalVolume* world = (*iterWorld);
742  // The idea now is to seek a PV with the same name and copy no
743  // in the hope it's the same one!!
744  G4PhysicalVolumeModel searchModel (world);
745  G4int verbosity = 0; // Suppress messages from G4PhysicalVolumeSearchScene.
746  G4PhysicalVolumeSearchScene searchScene
747  (&searchModel, fTopPVName, fTopPVCopyNo, verbosity);
748  G4ModelingParameters mp; // Default modeling parameters for this search.
750  searchModel.SetModelingParameters (&mp);
751  searchModel.DescribeYourselfTo (searchScene);
752  G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
753  if (foundVolume) {
754  if (foundVolume != fpTopPV && warn) {
755  G4cout <<
756  "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
757  "\n copy number (\""
758  << fTopPVName << "\", copy " << fTopPVCopyNo
759  << ") still exists and is being used."
760  "\n But it is not the same volume you originally specified"
761  "\n in /vis/scene/add/."
762  << G4endl;
763  }
764  fpTopPV = foundVolume;
765  CalculateExtent ();
766  found = true;
767  }
768  }
769  if (found) return true;
770  else {
771  if (warn) {
772  G4cout <<
773  "G4PhysicalVolumeModel::Validate(): No volume of name and"
774  "\n copy number (\""
775  << fTopPVName << "\", copy " << fTopPVCopyNo
776  << ") exists."
777  << G4endl;
778  }
779  return false;
780  }
781 }
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const
int G4int
Definition: G4Types.hh:78
const G4VisAttributes * GetDefaultVisAttributes() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
static G4TransportationManager * GetTransportationManager()
G4VPhysicalVolume * fpTopPV
#define G4endl
Definition: G4ios.hh:61
void SetDefaultVisAttributes(const G4VisAttributes *pDefaultVisAttributes)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ VisitGeometryAndGetVisReps()

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

Definition at line 217 of file G4PhysicalVolumeModel.cc.

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

Member Data Documentation

◆ fAbort

G4bool G4PhysicalVolumeModel::fAbort
mutableprotected

Definition at line 270 of file G4PhysicalVolumeModel.hh.

◆ fBaseFullPVPath

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

Definition at line 267 of file G4PhysicalVolumeModel.hh.

◆ fClippingMode

ClippingMode G4PhysicalVolumeModel::fClippingMode
protected

Definition at line 273 of file G4PhysicalVolumeModel.hh.

◆ fCurrentDepth

G4int G4PhysicalVolumeModel::fCurrentDepth
protected

Definition at line 262 of file G4PhysicalVolumeModel.hh.

◆ fCurtailDescent

G4bool G4PhysicalVolumeModel::fCurtailDescent
mutableprotected

Definition at line 271 of file G4PhysicalVolumeModel.hh.

◆ fDrawnPVPath

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

Definition at line 269 of file G4PhysicalVolumeModel.hh.

◆ fFullPVPath

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

Definition at line 268 of file G4PhysicalVolumeModel.hh.

◆ fpClippingSolid

G4VSolid* G4PhysicalVolumeModel::fpClippingSolid
protected

Definition at line 272 of file G4PhysicalVolumeModel.hh.

◆ fpCurrentLV

G4LogicalVolume* G4PhysicalVolumeModel::fpCurrentLV
protected

Definition at line 264 of file G4PhysicalVolumeModel.hh.

◆ fpCurrentMaterial

G4Material* G4PhysicalVolumeModel::fpCurrentMaterial
protected

Definition at line 265 of file G4PhysicalVolumeModel.hh.

◆ fpCurrentPV

G4VPhysicalVolume* G4PhysicalVolumeModel::fpCurrentPV
protected

Definition at line 263 of file G4PhysicalVolumeModel.hh.

◆ fpCurrentTransform

G4Transform3D* G4PhysicalVolumeModel::fpCurrentTransform
protected

Definition at line 266 of file G4PhysicalVolumeModel.hh.

◆ fpTopPV

G4VPhysicalVolume* G4PhysicalVolumeModel::fpTopPV
protected

Definition at line 256 of file G4PhysicalVolumeModel.hh.

◆ fRequestedDepth

G4int G4PhysicalVolumeModel::fRequestedDepth
protected

Definition at line 259 of file G4PhysicalVolumeModel.hh.

◆ fTopPVCopyNo

G4int G4PhysicalVolumeModel::fTopPVCopyNo
protected

Definition at line 258 of file G4PhysicalVolumeModel.hh.

◆ fTopPVName

G4String G4PhysicalVolumeModel::fTopPVName
protected

Definition at line 257 of file G4PhysicalVolumeModel.hh.

◆ fUseFullExtent

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 261 of file G4PhysicalVolumeModel.hh.


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