Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
< G4PhysicalVolumeNodeID
fBaseFullPVPath
 
std::vector
< G4PhysicalVolumeNodeID
fFullPVPath
 
std::vector
< G4PhysicalVolumeNodeID
fDrawnPVPath
 
G4bool fAbort
 
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
Enumerator
UNLIMITED 

Definition at line 86 of file G4PhysicalVolumeModel.hh.

Constructor & Destructor Documentation

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 }
G4Material * GetMaterial() const
G4String fType
Definition: G4VModel.hh:108
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 ( )
virtual

Definition at line 114 of file G4PhysicalVolumeModel.cc.

115 {
116  delete fpClippingSolid;
117 }

Member Function Documentation

void G4PhysicalVolumeModel::Abort ( ) const
inline

Definition at line 225 of file G4PhysicalVolumeModel.hh.

225 {fAbort = true;}
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  72); // 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
Transform3D inverse() const
Definition: Transform3D.cc:142
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
void DescribeYourselfTo(G4VGraphicsScene &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 841 of file G4PhysicalVolumeModel.cc.

842 {
843  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
844  std::ostringstream oss;
845  for (size_t i = 0; i < fFullPVPath.size(); ++i) {
846  oss << fFullPVPath[i].GetPhysicalVolume()->GetName()
847  << ':' << fFullPVPath[i].GetCopyNo();
848  if (i != fFullPVPath.size() - 1) oss << '/';
849  }
850 
851  if (!fpCurrentLV) {
853  ("G4PhysicalVolumeModel::CreateCurrentAttValues",
854  "modeling0004",
855  JustWarning,
856  "Current logical volume not defined.");
857  return values;
858  }
859 
860  values->push_back(G4AttValue("PVPath", oss.str(),""));
861  values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
862  G4VSolid* pSol = fpCurrentLV->GetSolid();
863  values->push_back(G4AttValue("Solid", pSol->GetName(),""));
864  values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
865  oss.str(""); oss << '\n' << *pSol;
866  values->push_back(G4AttValue("DmpSol", oss.str(),""));
867  const G4RotationMatrix localRotation = fpCurrentPV->GetObjectRotationValue();
868  const G4ThreeVector& localTranslation = fpCurrentPV->GetTranslation();
869  oss.str(""); oss << '\n' << G4Transform3D(localRotation,localTranslation);
870  values->push_back(G4AttValue("LocalTrans", oss.str(),""));
871  oss.str(""); oss << '\n' << *fpCurrentTransform;
872  values->push_back(G4AttValue("GlobalTrans", oss.str(),""));
873  G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
874  values->push_back(G4AttValue("Material", matName,""));
876  values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
878  oss.str(""); oss << matState;
879  values->push_back(G4AttValue("State", oss.str(),""));
881  values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
882  G4Region* region = fpCurrentLV->GetRegion();
883  G4String regionName = region? region->GetName(): G4String("No region");
884  values->push_back(G4AttValue("Region", regionName,""));
885  oss.str(""); oss << fpCurrentLV->IsRootRegion();
886  values->push_back(G4AttValue("RootRegion", oss.str(),""));
887  return values;
888 }
G4String GetName() const
const G4String & GetName() const
G4State
Definition: G4Material.hh:114
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetDensity() const
Definition: G4Material.hh:180
G4VSolid * GetSolid() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4Region * GetRegion() const
G4Transform3D * fpCurrentTransform
virtual G4GeometryType GetEntityType() const =0
G4VPhysicalVolume * fpCurrentPV
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
G4bool IsRootRegion() const
G4RotationMatrix GetObjectRotationValue() const
HepGeom::Transform3D G4Transform3D
G4double GetRadlen() const
Definition: G4Material.hh:220
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4ThreeVector & GetTranslation() const
G4State GetState() const
Definition: G4Material.hh:181
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4PhysicalVolumeModel::CurtailDescent ( ) const
inline

Definition at line 228 of file G4PhysicalVolumeModel.hh.

Here is the caller graph for this function:

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 - probably will, since that's
389  // the default for G4ModelingParameters. So create one if necessary...
390  if (!pVisAttribs) {
391  static G4VisAttributes defaultVisAttribs;
392  pVisAttribs = &defaultVisAttribs;
393  }
394  // From here, can assume pVisAttribs is a valid pointer. This is necessary
395  // because PreAddSolid needs a vis attributes object.
396 
397  // Make decision to draw...
398  G4bool thisToBeDrawn = true;
399 
400  // Update full path of physical volumes...
401  G4int copyNo = fpCurrentPV->GetCopyNo();
402  fFullPVPath.push_back
403  (G4PhysicalVolumeNodeID
405 
406  // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
407  const auto& vams = fpMP->GetVisAttributesModifiers();
408  if (vams.size()) {
409  // OK, we have some VAMs (Vis Attributes Modifiers).
410  for (const auto& vam: vams) {
411  const auto& vamPath = vam.GetPVNameCopyNoPath();
412  if (vamPath.size() == fFullPVPath.size()) {
413  // OK, we have a size match.
414  // Check the volume name/copy number path.
415  auto iVAMNameCopyNo = vamPath.begin();
416  auto iPVNodeId = fFullPVPath.begin();
417  for (; iVAMNameCopyNo != vamPath.end(); ++iVAMNameCopyNo, ++iPVNodeId) {
418  if (!(
419  iVAMNameCopyNo->GetName() ==
420  iPVNodeId->GetPhysicalVolume()->GetName() &&
421  iVAMNameCopyNo->GetCopyNo() ==
422  iPVNodeId->GetPhysicalVolume()->GetCopyNo()
423  )) {
424  // This path element does NOT match.
425  break;
426  }
427  }
428  if (iVAMNameCopyNo == vamPath.end()) {
429  // OK, the paths match (the above loop terminated normally).
430  // Create a vis atts object for the modified vis atts.
431  // It is static so that we may return a reliable pointer to it.
432  static G4VisAttributes modifiedVisAtts;
433  // Initialise it with the current vis atts and reset the pointer.
434  modifiedVisAtts = *pVisAttribs;
435  pVisAttribs = &modifiedVisAtts;
436  const G4VisAttributes& transVisAtts = vam.GetVisAttributes();
437  switch (vam.GetVisAttributesSignifier()) {
439  modifiedVisAtts.SetVisibility(transVisAtts.IsVisible());
440  break;
442  modifiedVisAtts.SetDaughtersInvisible
443  (transVisAtts.IsDaughtersInvisible());
444  break;
446  modifiedVisAtts.SetColour(transVisAtts.GetColour());
447  break;
449  modifiedVisAtts.SetLineStyle(transVisAtts.GetLineStyle());
450  break;
452  modifiedVisAtts.SetLineWidth(transVisAtts.GetLineWidth());
453  break;
455  if (transVisAtts.IsForceDrawingStyle()) {
456  if (transVisAtts.GetForcedDrawingStyle() ==
458  modifiedVisAtts.SetForceWireframe(true);
459  }
460  }
461  break;
463  if (transVisAtts.IsForceDrawingStyle()) {
464  if (transVisAtts.GetForcedDrawingStyle() ==
466  modifiedVisAtts.SetForceSolid(true);
467  }
468  }
469  break;
471  if (transVisAtts.IsForceAuxEdgeVisible()) {
472  modifiedVisAtts.SetForceAuxEdgeVisible
473  (transVisAtts.IsForcedAuxEdgeVisible());
474  }
475  break;
477  modifiedVisAtts.SetForceLineSegmentsPerCircle
478  (transVisAtts.GetForcedLineSegmentsPerCircle());
479  break;
480  }
481  }
482  }
483  }
484  }
485 
486  // There are various reasons why this volume
487  // might not be drawn...
488  G4bool culling = fpMP->IsCulling();
489  G4bool cullingInvisible = fpMP->IsCullingInvisible();
490  G4bool markedVisible = pVisAttribs->IsVisible();
491  G4bool cullingLowDensity = fpMP->IsDensityCulling();
492  G4double density = pMaterial? pMaterial->GetDensity(): 0;
493  G4double densityCut = fpMP -> GetVisibleDensity ();
494 
495  // 1) Global culling is on....
496  if (culling) {
497  // 2) Culling of invisible volumes is on...
498  if (cullingInvisible) {
499  // 3) ...and the volume is marked not visible...
500  if (!markedVisible) thisToBeDrawn = false;
501  }
502  // 4) Or culling of low density volumes is on...
503  if (cullingLowDensity) {
504  // 5) ...and density is less than cut value...
505  if (density < densityCut) thisToBeDrawn = false;
506  }
507  }
508  // 6) The user has asked for all further traversing to be aborted...
509  if (fAbort) thisToBeDrawn = false;
510 
511  // Record thisToBeDrawn in path...
512  fFullPVPath.back().SetDrawn(thisToBeDrawn);
513 
514  if (thisToBeDrawn) {
515 
516  // Update path of drawn physical volumes...
517  fDrawnPVPath.push_back
518  (G4PhysicalVolumeNodeID
519  (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform,thisToBeDrawn));
520 
521  if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
522  // For top-level drawn volumes, explode along radius...
524  G4Transform3D centred = centering.inverse() * theNewAT;
525  G4Scale3D oldScale;
526  G4Rotate3D oldRotation;
527  G4Translate3D oldTranslation;
528  centred.getDecomposition(oldScale, oldRotation, oldTranslation);
529  G4double explodeFactor = fpMP->GetExplodeFactor();
530  G4Translate3D newTranslation =
531  G4Translate3D(explodeFactor * oldTranslation.dx(),
532  explodeFactor * oldTranslation.dy(),
533  explodeFactor * oldTranslation.dz());
534  theNewAT = centering * newTranslation * oldRotation * oldScale;
535  }
536 
537  DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
538 
539  }
540 
541  // Make decision to draw daughters, if any. There are various
542  // reasons why daughters might not be drawn...
543 
544  // First, reasons that do not depend on culling policy...
545  G4int nDaughters = pLV->GetNoDaughters();
546  G4bool daughtersToBeDrawn = true;
547  // 1) There are no daughters...
548  if (!nDaughters) daughtersToBeDrawn = false;
549  // 2) We are at the limit if requested depth...
550  else if (requestedDepth == 0) daughtersToBeDrawn = false;
551  // 3) The user has asked for all further traversing to be aborted...
552  else if (fAbort) daughtersToBeDrawn = false;
553  // 4) The user has asked that the descent be curtailed...
554  else if (fCurtailDescent) daughtersToBeDrawn = false;
555 
556  // Now, reasons that depend on culling policy...
557  else {
558  G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
559  // Culling of covered daughters request. This is computed in
560  // G4VSceneHandler::CreateModelingParameters() depending on view
561  // parameters...
562  G4bool cullingCovered = fpMP->IsCullingCovered();
563  G4bool surfaceDrawing =
566  if (pVisAttribs->IsForceDrawingStyle()) {
567  switch (pVisAttribs->GetForcedDrawingStyle()) {
568  default:
569  case G4VisAttributes::wireframe: surfaceDrawing = false; break;
570  case G4VisAttributes::solid: surfaceDrawing = true; break;
571  }
572  }
573  G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
574  // 5) Global culling is on....
575  if (culling) {
576  // 6) ..and culling of invisible volumes is on...
577  if (cullingInvisible) {
578  // 7) ...and the mother requests daughters invisible
579  if (daughtersInvisible) daughtersToBeDrawn = false;
580  }
581  // 8) Or culling of covered daughters is requested...
582  if (cullingCovered) {
583  // 9) ...and surface drawing is operating...
584  if (surfaceDrawing) {
585  // 10) ...but only if mother is visible...
586  if (thisToBeDrawn) {
587  // 11) ...and opaque...
588  if (opaque) daughtersToBeDrawn = false;
589  }
590  }
591  }
592  }
593  }
594 
595  if (daughtersToBeDrawn) {
596  for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
597  // Store daughter pVPV in local variable ready for recursion...
598  G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
599  // Descend the geometry structure recursively...
600  fCurrentDepth++;
602  (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
603  fCurrentDepth--;
604  }
605  }
606 
607  // Reset for normal descending of next volume at this level...
608  fCurtailDescent = false;
609 
610  // Pop item from paths physical volumes...
611  fFullPVPath.pop_back();
612  if (thisToBeDrawn) {
613  fDrawnPVPath.pop_back();
614  }
615 }
G4bool IsForceAuxEdgeVisible() const
void SetColour(const G4Colour &)
G4double GetAlpha() const
Definition: G4Colour.hh:142
const G4Point3D & GetExplodeCentre() const
const G4VisAttributes * GetDefaultVisAttributes() const
G4double GetLineWidth() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:174
void SetLineStyle(LineStyle)
void SetLineWidth(G4double)
G4double GetDensity() const
Definition: G4Material.hh:180
G4bool IsVisible() const
const G4Colour & GetColour() const
Transform3D inverse() const
Definition: Transform3D.cc:142
double dx() const
Definition: Transform3D.h:279
G4bool IsExplode() const
G4bool IsDensityCulling() const
G4Transform3D * fpCurrentTransform
int G4int
Definition: G4Types.hh:78
G4VPhysicalVolume * fpCurrentPV
void SetForceSolid(G4bool=true)
LineStyle GetLineStyle() const
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
G4bool IsDaughtersInvisible() const
void SetVisibility(G4bool=true)
bool G4bool
Definition: G4Types.hh:79
G4bool IsCullingCovered() const
const G4VisAttributes * GetVisAttributes() const
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
G4bool IsForcedAuxEdgeVisible() const
virtual void DescribeSolid(const G4Transform3D &theAT, G4VSolid *pSol, const G4VisAttributes *pVisAttribs, G4VGraphicsScene &sceneHandler)
G4bool IsCullingInvisible() const
double dy() const
Definition: Transform3D.h:282
G4int GetNoDaughters() const
void SetDaughtersInvisible(G4bool=true)
G4int GetForcedLineSegmentsPerCircle() const
double dz() const
Definition: Transform3D.h:285
G4double GetExplodeFactor() const
G4bool IsCulling() const
G4bool IsForceDrawingStyle() const
virtual G4int GetCopyNo() const =0
HepGeom::Translate3D G4Translate3D
void SetForceAuxEdgeVisible(G4bool=true)
std::vector< G4PhysicalVolumeNodeID > fDrawnPVPath
double G4double
Definition: G4Types.hh:76
void VisitGeometryAndGetVisReps(G4VPhysicalVolume *, G4int requestedDepth, const G4Transform3D &, G4VGraphicsScene &)
void SetForceWireframe(G4bool=true)
ForcedDrawingStyle GetForcedDrawingStyle() const
DrawingStyle GetDrawingStyle() const
void SetForceLineSegmentsPerCircle(G4int nSegments)
const std::vector< VisAttributesModifier > & GetVisAttributesModifiers() const

Here is the call graph for this function:

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

Reimplemented in G4LogicalVolumeModel.

Definition at line 618 of file G4PhysicalVolumeModel.cc.

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

Here is the call graph for this function:

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:

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

Definition at line 764 of file G4PhysicalVolumeModel.cc.

765 {
766  G4bool isNew;
767  std::map<G4String,G4AttDef>* store
768  = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
769  if (isNew) {
770  (*store)["PVPath"] =
771  G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
772  (*store)["LVol"] =
773  G4AttDef("LVol","Logical Volume","Physics","","G4String");
774  (*store)["Solid"] =
775  G4AttDef("Solid","Solid Name","Physics","","G4String");
776  (*store)["EType"] =
777  G4AttDef("EType","Entity Type","Physics","","G4String");
778  (*store)["DmpSol"] =
779  G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
780  (*store)["LocalTrans"] =
781  G4AttDef("LocalTrans","Local transformation of volume","Physics","","G4String");
782  (*store)["GlobalTrans"] =
783  G4AttDef("GlobalTrans","Global transformation of volume","Physics","","G4String");
784  (*store)["Material"] =
785  G4AttDef("Material","Material Name","Physics","","G4String");
786  (*store)["Density"] =
787  G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
788  (*store)["State"] =
789  G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
790  (*store)["Radlen"] =
791  G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
792  (*store)["Region"] =
793  G4AttDef("Region","Cuts Region","Physics","","G4String");
794  (*store)["RootRegion"] =
795  G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
796  }
797  return store;
798 }
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:

const G4VSolid* G4PhysicalVolumeModel::GetClippingSolid ( ) const
inline

Definition at line 158 of file G4PhysicalVolumeModel.hh.

159  {return fpClippingSolid;}
G4int G4PhysicalVolumeModel::GetCurrentDepth ( ) const
inline

Definition at line 161 of file G4PhysicalVolumeModel.hh.

Here is the caller graph for this function:

G4String G4PhysicalVolumeModel::GetCurrentDescription ( ) const
virtual

Reimplemented from G4VModel.

Definition at line 211 of file G4PhysicalVolumeModel.cc.

212 {
213  return "G4PhysicalVolumeModel " + GetCurrentTag ();
214 }
G4String GetCurrentTag() const

Here is the call graph for this function:

G4LogicalVolume* G4PhysicalVolumeModel::GetCurrentLV ( ) const
inline

Definition at line 167 of file G4PhysicalVolumeModel.hh.

167 {return fpCurrentLV;}

Here is the caller graph for this function:

G4Material* G4PhysicalVolumeModel::GetCurrentMaterial ( ) const
inline

Definition at line 170 of file G4PhysicalVolumeModel.hh.

170 {return fpCurrentMaterial;}

Here is the caller graph for this function:

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:

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:

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:

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 caller graph for this function:

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:

G4int G4PhysicalVolumeModel::GetRequestedDepth ( ) const
inline

Definition at line 156 of file G4PhysicalVolumeModel.hh.

Here is the caller graph for this function:

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:

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:

void G4PhysicalVolumeModel::SetClippingMode ( ClippingMode  mode)
inline

Definition at line 218 of file G4PhysicalVolumeModel.hh.

218  {
219  fClippingMode = mode;
220  }
void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid pClippingSolid)
inline

Definition at line 214 of file G4PhysicalVolumeModel.hh.

214  {
215  fpClippingSolid = pClippingSolid;
216  }
void G4PhysicalVolumeModel::SetRequestedDepth ( G4int  requestedDepth)
inline

Definition at line 210 of file G4PhysicalVolumeModel.hh.

210  {
211  fRequestedDepth = requestedDepth;
212  }
G4bool G4PhysicalVolumeModel::Validate ( G4bool  warn)
virtual

Reimplemented from G4VModel.

Definition at line 710 of file G4PhysicalVolumeModel.cc.

711 {
712  G4TransportationManager* transportationManager =
714 
715  size_t nWorlds = transportationManager->GetNoWorlds();
716 
717  G4bool found = false;
718 
719  std::vector<G4VPhysicalVolume*>::iterator iterWorld =
720  transportationManager->GetWorldsIterator();
721  for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
722  G4VPhysicalVolume* world = (*iterWorld);
723  // The idea now is to seek a PV with the same name and copy no
724  // in the hope it's the same one!!
725  G4PhysicalVolumeModel searchModel (world);
726  G4int verbosity = 0; // Suppress messages from G4PhysicalVolumeSearchScene.
727  G4PhysicalVolumeSearchScene searchScene
728  (&searchModel, fTopPVName, fTopPVCopyNo, verbosity);
729  G4ModelingParameters mp; // Default modeling parameters for this search.
731  searchModel.SetModelingParameters (&mp);
732  searchModel.DescribeYourselfTo (searchScene);
733  G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
734  if (foundVolume) {
735  if (foundVolume != fpTopPV && warn) {
736  G4cout <<
737  "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
738  "\n copy number (\""
739  << fTopPVName << "\", copy " << fTopPVCopyNo
740  << ") still exists and is being used."
741  "\n But it is not the same volume you originally specified"
742  "\n in /vis/scene/add/."
743  << G4endl;
744  }
745  fpTopPV = foundVolume;
746  CalculateExtent ();
747  found = true;
748  }
749  }
750  if (found) return true;
751  else {
752  if (warn) {
753  G4cout <<
754  "G4PhysicalVolumeModel::Validate(): No volume of name and"
755  "\n copy number (\""
756  << fTopPVName << "\", copy " << fTopPVCopyNo
757  << ") exists."
758  << G4endl;
759  }
760  return false;
761  }
762 }
const G4VisAttributes * GetDefaultVisAttributes() const
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4ModelingParameters * fpMP
Definition: G4VModel.hh:113
static G4TransportationManager * GetTransportationManager()
G4VPhysicalVolume * fpTopPV
size_t GetNoWorlds() const
#define G4endl
Definition: G4ios.hh:61
void SetDefaultVisAttributes(const G4VisAttributes *pDefaultVisAttributes)

Here is the call graph for this function:

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
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
Definition: G4Tubs.hh:85
virtual G4GeometryType GetEntityType() const =0
int G4int
Definition: G4Types.hh:78
std::vector< G4PhysicalVolumeNodeID > fFullPVPath
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
void DescribeAndDescend(G4VPhysicalVolume *, G4int requestedDepth, G4LogicalVolume *, G4VSolid *, G4Material *, const G4Transform3D &, G4VGraphicsScene &)
G4bool IsWarning() const

Here is the call graph for this function:

Member Data Documentation

G4bool G4PhysicalVolumeModel::fAbort
mutableprotected

Definition at line 270 of file G4PhysicalVolumeModel.hh.

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

Definition at line 267 of file G4PhysicalVolumeModel.hh.

ClippingMode G4PhysicalVolumeModel::fClippingMode
protected

Definition at line 273 of file G4PhysicalVolumeModel.hh.

G4int G4PhysicalVolumeModel::fCurrentDepth
protected

Definition at line 262 of file G4PhysicalVolumeModel.hh.

G4bool G4PhysicalVolumeModel::fCurtailDescent
mutableprotected

Definition at line 271 of file G4PhysicalVolumeModel.hh.

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

Definition at line 269 of file G4PhysicalVolumeModel.hh.

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

Definition at line 268 of file G4PhysicalVolumeModel.hh.

G4VSolid* G4PhysicalVolumeModel::fpClippingSolid
protected

Definition at line 272 of file G4PhysicalVolumeModel.hh.

G4LogicalVolume* G4PhysicalVolumeModel::fpCurrentLV
protected

Definition at line 264 of file G4PhysicalVolumeModel.hh.

G4Material* G4PhysicalVolumeModel::fpCurrentMaterial
protected

Definition at line 265 of file G4PhysicalVolumeModel.hh.

G4VPhysicalVolume* G4PhysicalVolumeModel::fpCurrentPV
protected

Definition at line 263 of file G4PhysicalVolumeModel.hh.

G4Transform3D* G4PhysicalVolumeModel::fpCurrentTransform
protected

Definition at line 266 of file G4PhysicalVolumeModel.hh.

G4VPhysicalVolume* G4PhysicalVolumeModel::fpTopPV
protected

Definition at line 256 of file G4PhysicalVolumeModel.hh.

G4int G4PhysicalVolumeModel::fRequestedDepth
protected

Definition at line 259 of file G4PhysicalVolumeModel.hh.

G4int G4PhysicalVolumeModel::fTopPVCopyNo
protected

Definition at line 258 of file G4PhysicalVolumeModel.hh.

G4String G4PhysicalVolumeModel::fTopPVName
protected

Definition at line 257 of file G4PhysicalVolumeModel.hh.

G4bool G4PhysicalVolumeModel::fUseFullExtent
protected

Definition at line 261 of file G4PhysicalVolumeModel.hh.


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