G4PhysicalVolumeModel Class Reference

#include <G4PhysicalVolumeModel.hh>

Inheritance diagram for G4PhysicalVolumeModel:

G4VModel G4FlavoredParallelWorldModel G4LogicalVolumeModel

Public Types

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

Public Member Functions

 G4PhysicalVolumeModel (G4VPhysicalVolume *=0, G4int requestedDepth=UNLIMITED, const G4Transform3D &modelTransformation=G4Transform3D(), const G4ModelingParameters *=0, G4bool useFullExtent=false)
virtual ~G4PhysicalVolumeModel ()
void DescribeYourselfTo (G4VGraphicsScene &)
G4String GetCurrentDescription () const
G4String GetCurrentTag () const
G4VPhysicalVolumeGetTopPhysicalVolume () const
G4int GetRequestedDepth () const
const G4VSolidGetClippingSolid () const
G4int GetCurrentDepth () const
G4VPhysicalVolumeGetCurrentPV () const
G4LogicalVolumeGetCurrentLV () const
G4MaterialGetCurrentMaterial () const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath () const
const std::vector< G4PhysicalVolumeNodeID > & GetDrawnPVPath () const
const std::map< G4String,
G4AttDef > * 
GetAttDefs () const
std::vector< G4AttValue > * CreateCurrentAttValues () const
void SetBaseFullPVPath (const std::vector< G4PhysicalVolumeNodeID > &baseFullPVPath)
void SetRequestedDepth (G4int requestedDepth)
void SetClippingSolid (G4VSolid *pClippingSolid)
void SetClippingMode (ClippingMode mode)
G4bool Validate (G4bool warn)
void CurtailDescent ()

Protected Member Functions

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

Protected Attributes

G4VPhysicalVolumefpTopPV
G4String fTopPVName
G4int fTopPVCopyNo
G4int fRequestedDepth
G4bool fUseFullExtent
G4int fCurrentDepth
G4VPhysicalVolumefpCurrentPV
G4LogicalVolumefpCurrentLV
G4MaterialfpCurrentMaterial
G4Transform3DfpCurrentTransform
std::vector< G4PhysicalVolumeNodeIDfBaseFullPVPath
std::vector< G4PhysicalVolumeNodeIDfFullPVPath
std::vector< G4PhysicalVolumeNodeIDfDrawnPVPath
G4bool fCurtailDescent
G4VSolidfpClippingSolid
ClippingMode fClippingMode

Data Structures

class  G4PhysicalVolumeModelTouchable
class  G4PhysicalVolumeNodeID

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.

00086 {UNLIMITED = -1};

enum G4PhysicalVolumeModel::ClippingMode

Enumerator:
subtraction 
intersection 

Definition at line 88 of file G4PhysicalVolumeModel.hh.


Constructor & Destructor Documentation

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

Definition at line 59 of file G4PhysicalVolumeModel.cc.

00063                       :
00064   G4VModel        (modelTransformation, pMP),
00065   fpTopPV         (pVPV),
00066   fRequestedDepth (requestedDepth),
00067   fUseFullExtent  (useFullExtent),
00068   fCurrentDepth   (0),
00069   fpCurrentPV     (0),
00070   fpCurrentLV     (0),
00071   fpCurrentMaterial (0),
00072   fpCurrentTransform (0),
00073   fCurtailDescent (false),
00074   fpClippingSolid (0),
00075   fClippingMode   (subtraction)
00076 {
00077   if (fpTopPV) {
00078     
00079     fType = "G4PhysicalVolumeModel";
00080     fTopPVName = fpTopPV -> GetName ();
00081     fTopPVCopyNo = fpTopPV -> GetCopyNo ();
00082     std::ostringstream o;
00083     o << fpTopPV -> GetCopyNo ();
00084     fGlobalTag = fpTopPV -> GetName () + "." + o.str();
00085     fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag;
00086 
00087     CalculateExtent ();
00088   }
00089 }

G4PhysicalVolumeModel::~G4PhysicalVolumeModel (  )  [virtual]

Definition at line 91 of file G4PhysicalVolumeModel.cc.

References fpClippingSolid.

00092 {
00093   delete fpClippingSolid;
00094 }


Member Function Documentation

void G4PhysicalVolumeModel::CalculateExtent (  )  [protected]

Definition at line 96 of file G4PhysicalVolumeModel.cc.

References DescribeYourselfTo(), G4VModel::fExtent, G4VModel::fpMP, fpTopPV, fRequestedDepth, G4VModel::fTransform, fUseFullExtent, G4BoundingSphereScene::GetCentre(), G4VModel::GetExtent(), G4BoundingSphereScene::GetRadius(), and G4ModelingParameters::wf.

Referenced by Validate().

00097 {
00098   if (fUseFullExtent) {
00099     fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
00100   }
00101   else {
00102     G4BoundingSphereScene bsScene(this);
00103     const G4int tempRequestedDepth = fRequestedDepth;
00104     fRequestedDepth = -1;  // Always search to all depths to define extent.
00105     const G4ModelingParameters* tempMP = fpMP;
00106     G4ModelingParameters mParams
00107       (0,      // No default vis attributes needed.
00108        G4ModelingParameters::wf,  // wireframe (not relevant for this).
00109        true,   // Global culling.
00110        true,   // Cull invisible volumes.
00111        false,  // Density culling.
00112        0.,     // Density (not relevant if density culling false).
00113        true,   // Cull daughters of opaque mothers.
00114        24);    // No of sides (not relevant for this operation).
00115     fpMP = &mParams;
00116     DescribeYourselfTo (bsScene);
00117     G4double radius = bsScene.GetRadius();
00118     if (radius < 0.) {  // Nothing in the scene.
00119       fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent ();
00120     } else {
00121       // Transform back to coordinates relative to the top
00122       // transformation, which is in G4VModel::fTransform.  This makes
00123       // it conform to all models, which are defined by a
00124       // transformation and an extent relative to that
00125       // transformation...
00126       G4Point3D centre = bsScene.GetCentre();
00127       centre.transform(fTransform.inverse());
00128       fExtent = G4VisExtent(centre, radius);
00129     }
00130     fpMP = tempMP;
00131     fRequestedDepth = tempRequestedDepth;
00132   }
00133 }

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

Definition at line 829 of file G4PhysicalVolumeModel.cc.

References fFullPVPath, fpCurrentLV, fpCurrentMaterial, fpCurrentTransform, G4BestUnit, G4Exception(), G4Material::GetDensity(), G4VSolid::GetEntityType(), G4Region::GetName(), G4Material::GetName(), G4VSolid::GetName(), G4LogicalVolume::GetName(), G4Material::GetRadlen(), G4LogicalVolume::GetRegion(), G4LogicalVolume::GetSolid(), G4Material::GetState(), G4LogicalVolume::IsRootRegion(), JustWarning, and kStateUndefined.

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

00830 {
00831   std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
00832   std::ostringstream oss;
00833   for (size_t i = 0; i < fFullPVPath.size(); ++i) {
00834     oss << fFullPVPath[i].GetPhysicalVolume()->GetName()
00835         << ':' << fFullPVPath[i].GetCopyNo();
00836     if (i != fFullPVPath.size() - 1) oss << '/';
00837   }
00838 
00839   if (!fpCurrentLV) {
00840      G4Exception
00841         ("G4PhysicalVolumeModel::CreateCurrentAttValues",
00842          "modeling0004",
00843          JustWarning,
00844          "Current logical volume not defined.");
00845      return values;
00846   }
00847 
00848   values->push_back(G4AttValue("PVPath", oss.str(),""));
00849   values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),""));
00850   G4VSolid* pSol = fpCurrentLV->GetSolid();
00851   values->push_back(G4AttValue("Solid", pSol->GetName(),""));
00852   values->push_back(G4AttValue("EType", pSol->GetEntityType(),""));
00853   oss.str(""); oss << '\n' << *pSol;
00854   values->push_back(G4AttValue("DmpSol", oss.str(),""));
00855   oss.str(""); oss << '\n' << *fpCurrentTransform;
00856   values->push_back(G4AttValue("Trans", oss.str(),""));
00857   G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material");
00858   values->push_back(G4AttValue("Material", matName,""));
00859   G4double matDensity = fpCurrentMaterial? fpCurrentMaterial->GetDensity(): 0.;
00860   values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),""));
00861   G4State matState = fpCurrentMaterial? fpCurrentMaterial->GetState(): kStateUndefined;
00862   oss.str(""); oss << matState;
00863   values->push_back(G4AttValue("State", oss.str(),""));
00864   G4double matRadlen = fpCurrentMaterial? fpCurrentMaterial->GetRadlen(): 0.;
00865   values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),""));
00866   G4Region* region = fpCurrentLV->GetRegion();
00867   G4String regionName = region? region->GetName(): G4String("No region");
00868   values->push_back(G4AttValue("Region", regionName,""));
00869   oss.str(""); oss << fpCurrentLV->IsRootRegion();
00870   values->push_back(G4AttValue("RootRegion", oss.str(),""));
00871   return values;
00872 }

void G4PhysicalVolumeModel::CurtailDescent (  )  [inline]

Definition at line 222 of file G4PhysicalVolumeModel.hh.

References fCurtailDescent.

Referenced by G4ASCIITreeSceneHandler::RequestPrimitives().

00222 {fCurtailDescent = true;}

void G4PhysicalVolumeModel::DescribeAndDescend ( G4VPhysicalVolume ,
G4int  requestedDepth,
G4LogicalVolume ,
G4VSolid ,
G4Material ,
const G4Transform3D ,
G4VGraphicsScene  
) [protected]

Definition at line 335 of file G4PhysicalVolumeModel.cc.

References G4Colour::GetAlpha(), G4VisAttributes::GetColour(), G4Material::GetDensity(), G4VisAttributes::GetForcedDrawingStyle(), G4VisAttributes::GetForcedLineSegmentsPerCircle(), G4VisAttributes::GetLineStyle(), G4VisAttributes::GetLineWidth(), G4LogicalVolume::GetNoDaughters(), G4LogicalVolume::GetVisAttributes(), G4ModelingParameters::hlhsr, G4ModelingParameters::hsr, G4VisAttributes::IsDaughtersInvisible(), G4VisAttributes::IsForceAuxEdgeVisible(), G4VisAttributes::IsForceDrawingStyle(), G4VisAttributes::IsVisible(), G4VisAttributes::SetColour(), G4VisAttributes::SetDaughtersInvisible(), G4VisAttributes::SetForceAuxEdgeVisible(), G4VisAttributes::SetForceLineSegmentsPerCircle(), G4VisAttributes::SetForceSolid(), G4VisAttributes::SetForceWireframe(), G4VisAttributes::SetLineStyle(), G4VisAttributes::SetLineWidth(), G4VisAttributes::SetVisibility(), G4VisAttributes::solid, G4ModelingParameters::VASColour, G4ModelingParameters::VASDaughtersInvisible, G4ModelingParameters::VASForceAuxEdgeVisible, G4ModelingParameters::VASForceLineSegmentsPerCircle, G4ModelingParameters::VASForceSolid, G4ModelingParameters::VASForceWireframe, G4ModelingParameters::VASLineStyle, G4ModelingParameters::VASLineWidth, G4ModelingParameters::VASVisibility, and G4VisAttributes::wireframe.

00342 {
00343   // Maintain useful data members...
00344   fpCurrentPV = pVPV;
00345   fpCurrentLV = pLV;
00346   fpCurrentMaterial = pMaterial;
00347 
00348   const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue ();
00349   const G4ThreeVector&  translation     = pVPV -> GetTranslation ();
00350   G4Transform3D theLT (G4Transform3D (objectRotation, translation));
00351 
00352   // Compute the accumulated transformation...
00353   // Note that top volume's transformation relative to the world
00354   // coordinate system is specified in theAT == startingTransformation
00355   // = fTransform (see DescribeYourselfTo), so first time through the
00356   // volume's own transformation, which is only relative to its
00357   // mother, i.e., not relative to the world coordinate system, should
00358   // not be accumulated.
00359   G4Transform3D theNewAT (theAT);
00360   if (fCurrentDepth != 0) theNewAT = theAT * theLT;
00361   fpCurrentTransform = &theNewAT;
00362 
00363   const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes();
00364   if (!pVisAttribs) pVisAttribs = fpMP->GetDefaultVisAttributes();
00365   // Beware - pVisAttribs might still be zero - create a temporary default one...
00366   G4bool visAttsCreated = false;
00367   if (!pVisAttribs) {
00368     pVisAttribs = new G4VisAttributes;
00369     visAttsCreated = true;
00370   }
00371   // From here, can assume pVisAttribs is a valid pointer.
00372   
00373   // Make decision to draw...
00374   G4bool thisToBeDrawn = true;
00375   
00376   // Update full path of physical volumes...
00377   G4int copyNo = fpCurrentPV->GetCopyNo();
00378   fFullPVPath.push_back
00379   (G4PhysicalVolumeNodeID
00380    (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform));
00381   
00382   // In case we need to copy the vis atts for modification...
00383   G4bool copyForVAM = false;
00384   const G4VisAttributes* pUnmodifiedVisAtts = pVisAttribs;
00385   G4VisAttributes* pModifiedVisAtts = 0;
00386   
00387   // Check if vis attributes are to be modified by a /vis/touchable/set/ command.
00388   const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
00389     fpMP->GetVisAttributesModifiers();
00390   std::vector<G4ModelingParameters::VisAttributesModifier>::const_iterator
00391     iModifier;
00392   for (iModifier = vams.begin();
00393        iModifier != vams.end();
00394        ++iModifier) {
00395     const G4ModelingParameters::PVNameCopyNoPath& vamPath =
00396       iModifier->GetPVNameCopyNoPath();
00397     if (vamPath.size() == fFullPVPath.size()) {
00398       // OK - there's a size match.  Check it out.
00399 //      G4cout << "Size match" << G4endl;
00400       G4ModelingParameters::PVNameCopyNoPathConstIterator iVAMNameCopyNo;
00401       std::vector<G4PhysicalVolumeNodeID>::const_iterator iPVNodeId;
00402       for (iVAMNameCopyNo = vamPath.begin(), iPVNodeId = fFullPVPath.begin();
00403            iVAMNameCopyNo != vamPath.end();
00404            ++iVAMNameCopyNo, ++iPVNodeId) {
00405 //        G4cout
00406 //        << iVAMNameCopyNo->fName
00407 //        << ',' << iVAMNameCopyNo->fCopyNo
00408 //        << "; " << iPVNodeId->GetPhysicalVolume()->GetName()
00409 //        << ','  << iPVNodeId->GetPhysicalVolume()->GetCopyNo()
00410 //        << G4endl;
00411         if (!(
00412               iVAMNameCopyNo->GetName() ==
00413               iPVNodeId->GetPhysicalVolume()->GetName() &&
00414               iVAMNameCopyNo->GetCopyNo() ==
00415               iPVNodeId->GetPhysicalVolume()->GetCopyNo()
00416               )) {
00417           break;
00418         }
00419       }
00420       if (iVAMNameCopyNo == vamPath.end()) {
00421 //        G4cout << "Match found" << G4endl;
00422         if (!copyForVAM) {
00423           pModifiedVisAtts = new G4VisAttributes(*pUnmodifiedVisAtts);
00424           pVisAttribs = pModifiedVisAtts;
00425           copyForVAM = true;
00426         }
00427         const G4VisAttributes& transVisAtts = iModifier->GetVisAttributes();
00428         switch (iModifier->GetVisAttributesSignifier()) {
00429           case G4ModelingParameters::VASVisibility:
00430             pModifiedVisAtts->SetVisibility(transVisAtts.IsVisible());
00431             break;
00432           case G4ModelingParameters::VASDaughtersInvisible:
00433             pModifiedVisAtts->SetDaughtersInvisible
00434             (transVisAtts.IsDaughtersInvisible());
00435             break;
00436           case G4ModelingParameters::VASColour:
00437             pModifiedVisAtts->SetColour(transVisAtts.GetColour());
00438             break;
00439           case G4ModelingParameters::VASLineStyle:
00440             pModifiedVisAtts->SetLineStyle(transVisAtts.GetLineStyle());
00441             break;
00442           case G4ModelingParameters::VASLineWidth:
00443             pModifiedVisAtts->SetLineWidth(transVisAtts.GetLineWidth());
00444             break;
00445           case G4ModelingParameters::VASForceWireframe:
00446             if (transVisAtts.GetForcedDrawingStyle() ==
00447                 G4VisAttributes::wireframe) {
00448               pModifiedVisAtts->SetForceWireframe
00449               (transVisAtts.IsForceDrawingStyle());
00450             }
00451             break;
00452           case G4ModelingParameters::VASForceSolid:
00453             if (transVisAtts.GetForcedDrawingStyle() ==
00454                 G4VisAttributes::solid) {
00455               pModifiedVisAtts->SetForceSolid
00456               (transVisAtts.IsForceDrawingStyle());
00457             }
00458             break;
00459           case G4ModelingParameters::VASForceAuxEdgeVisible:
00460             pModifiedVisAtts->SetForceAuxEdgeVisible
00461             (transVisAtts.IsForceAuxEdgeVisible());
00462             break;
00463           case G4ModelingParameters::VASForceLineSegmentsPerCircle:
00464             pModifiedVisAtts->SetForceLineSegmentsPerCircle
00465             (transVisAtts.GetForcedLineSegmentsPerCircle());
00466             break;
00467         }
00468       }
00469     }
00470   }
00471 
00472   // There are various reasons why this volume
00473   // might not be drawn...
00474   G4bool culling = fpMP->IsCulling();
00475   G4bool cullingInvisible = fpMP->IsCullingInvisible();
00476   G4bool markedVisible = pVisAttribs->IsVisible();
00477   G4bool cullingLowDensity = fpMP->IsDensityCulling();
00478   G4double density = pMaterial? pMaterial->GetDensity(): 0;
00479   G4double densityCut = fpMP -> GetVisibleDensity ();
00480 
00481   // 1) Global culling is on....
00482   if (culling) {
00483     // 2) Culling of invisible volumes is on...
00484     if (cullingInvisible) {
00485       // 3) ...and the volume is marked not visible...
00486       if (!markedVisible) thisToBeDrawn = false;
00487     }
00488     // 4) Or culling of low density volumes is on...
00489     if (cullingLowDensity) {
00490       // 5) ...and density is less than cut value...
00491       if (density < densityCut) thisToBeDrawn = false;
00492     }
00493   }
00494   
00495   // Record thisToBeDrawn in path...
00496   fFullPVPath.back().SetDrawn(thisToBeDrawn);
00497 
00498   if (thisToBeDrawn) {
00499 
00500     // Update path of drawn physical volumes...
00501     fDrawnPVPath.push_back
00502       (G4PhysicalVolumeNodeID
00503        (fpCurrentPV,copyNo,fCurrentDepth,*fpCurrentTransform,thisToBeDrawn));
00504 
00505     if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) {
00506       // For top-level drawn volumes, explode along radius...
00507       G4Transform3D centering = G4Translate3D(fpMP->GetExplodeCentre());
00508       G4Transform3D centred = centering.inverse() * theNewAT;
00509       G4Scale3D oldScale;
00510       G4Rotate3D oldRotation;
00511       G4Translate3D oldTranslation;
00512       centred.getDecomposition(oldScale, oldRotation, oldTranslation);
00513       G4double explodeFactor = fpMP->GetExplodeFactor();
00514       G4Translate3D newTranslation =
00515         G4Translate3D(explodeFactor * oldTranslation.dx(),
00516                       explodeFactor * oldTranslation.dy(),
00517                       explodeFactor * oldTranslation.dz());
00518       theNewAT = centering * newTranslation * oldRotation * oldScale;
00519     }
00520 
00521     DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler);
00522 
00523   }
00524 
00525   // Make decision to draw daughters, if any.  There are various
00526   // reasons why daughters might not be drawn...
00527 
00528   // First, reasons that do not depend on culling policy...
00529   G4int nDaughters = pLV->GetNoDaughters();
00530   G4bool daughtersToBeDrawn = true;
00531   // 1) There are no daughters...
00532   if (!nDaughters) daughtersToBeDrawn = false;
00533   // 2) We are at the limit if requested depth...
00534   else if (requestedDepth == 0) daughtersToBeDrawn = false;
00535   // 3) The user has asked that the descent be curtailed...
00536   else if (fCurtailDescent) daughtersToBeDrawn = false;
00537 
00538   // Now, reasons that depend on culling policy...
00539   else {
00540     G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible();
00541     // Culling of covered daughters request.  This is computed in
00542     // G4VSceneHandler::CreateModelingParameters() depending on view
00543     // parameters...
00544     G4bool cullingCovered = fpMP->IsCullingCovered();
00545     G4bool surfaceDrawing =
00546       fpMP->GetDrawingStyle() == G4ModelingParameters::hsr ||
00547       fpMP->GetDrawingStyle() == G4ModelingParameters::hlhsr;    
00548     if (pVisAttribs->IsForceDrawingStyle()) {
00549       switch (pVisAttribs->GetForcedDrawingStyle()) {
00550       default:
00551       case G4VisAttributes::wireframe: surfaceDrawing = false; break;
00552       case G4VisAttributes::solid: surfaceDrawing = true; break;
00553       }
00554     }
00555     G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.;
00556     // 4) Global culling is on....
00557     if (culling) {
00558       // 5) ..and culling of invisible volumes is on...
00559       if (cullingInvisible) {
00560         // 6) ...and the mother requests daughters invisible
00561         if (daughtersInvisible) daughtersToBeDrawn = false;
00562       }
00563       // 7) Or culling of covered daughters is requested...
00564       if (cullingCovered) {
00565         // 8) ...and surface drawing is operating...
00566         if (surfaceDrawing) {
00567           // 9) ...but only if mother is visible...
00568           if (thisToBeDrawn) {
00569             // 10) ...and opaque...
00570               if (opaque) daughtersToBeDrawn = false;
00571           }
00572         }
00573       }
00574     }
00575   }
00576 
00577   // Delete modified vis atts if created...
00578   if (copyForVAM) {
00579     delete pModifiedVisAtts;
00580     pVisAttribs = pUnmodifiedVisAtts;
00581   }
00582   
00583   // Vis atts for this volume no longer needed if created...
00584   if (visAttsCreated) delete pVisAttribs;
00585 
00586   if (daughtersToBeDrawn) {
00587     for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
00588       // Store daughter pVPV in local variable ready for recursion...
00589       G4VPhysicalVolume* pDaughterVPV = pLV -> GetDaughter (iDaughter);
00590       // Descend the geometry structure recursively...
00591       fCurrentDepth++;
00592       VisitGeometryAndGetVisReps
00593         (pDaughterVPV, requestedDepth - 1, theNewAT, sceneHandler);
00594       fCurrentDepth--;
00595     }
00596   }
00597 
00598   // Reset for normal descending of next volume at this level...
00599   fCurtailDescent = false;
00600 
00601   // Pop item from paths physical volumes...
00602   fFullPVPath.pop_back();
00603   if (thisToBeDrawn) {
00604     fDrawnPVPath.pop_back();
00605   }
00606 }

void G4PhysicalVolumeModel::DescribeSolid ( const G4Transform3D theAT,
G4VSolid pSol,
const G4VisAttributes pVisAttribs,
G4VGraphicsScene sceneHandler 
) [protected, virtual]

Reimplemented in G4LogicalVolumeModel.

Definition at line 609 of file G4PhysicalVolumeModel.cc.

References G4VGraphicsScene::AddPrimitive(), G4VGraphicsScene::BeginPrimitives(), G4VGraphicsScene::EndPrimitives(), G4cout, G4endl, G4VisAttributes::GetForcedLineSegmentsPerCircle(), G4VSolid::GetName(), G4VSolid::GetPolyhedron(), G4VisAttributes::IsForceLineSegmentsPerCircle(), G4VGraphicsScene::PostAddSolid(), G4VGraphicsScene::PreAddSolid(), G4Colour::Red(), G4VisAttributes::SetColour(), and G4Visible::SetVisAttributes().

00613 {
00614   sceneHandler.PreAddSolid (theAT, *pVisAttribs);
00615 
00616   G4VSolid* pSectionSolid = fpMP->GetSectionSolid();
00617   G4VSolid* pCutawaySolid = fpMP->GetCutawaySolid();
00618 
00619   if (!fpClippingSolid && !pSectionSolid && !pCutawaySolid) {
00620 
00621     pSol -> DescribeYourselfTo (sceneHandler);  // Standard treatment.
00622 
00623   } else {
00624 
00625     // Clipping, etc., performed by Boolean operations.
00626 
00627     // First, get polyhedron for current solid...
00628     if (pVisAttribs->IsForceLineSegmentsPerCircle())
00629       G4Polyhedron::SetNumberOfRotationSteps
00630         (pVisAttribs->GetForcedLineSegmentsPerCircle());
00631     else
00632       G4Polyhedron::SetNumberOfRotationSteps(fpMP->GetNoOfSides());
00633     const G4Polyhedron* pOriginal = pSol->GetPolyhedron();
00634     G4Polyhedron::ResetNumberOfRotationSteps();
00635 
00636     if (!pOriginal) {
00637 
00638       if (fpMP->IsWarning())
00639         G4cout <<
00640           "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n  \""
00641                << pSol->GetName() <<
00642           "\" has no polyhedron.  Cannot by clipped."
00643                << G4endl;
00644       pSol -> DescribeYourselfTo (sceneHandler);  // Standard treatment.
00645 
00646     } else {
00647 
00648       G4Polyhedron resultant(*pOriginal);
00649       G4VisAttributes resultantVisAttribs(*pVisAttribs);
00650       G4VSolid* resultantSolid = 0;
00651 
00652       if (fpClippingSolid) {
00653         switch (fClippingMode) {
00654         default:
00655         case subtraction:
00656           resultantSolid = new G4SubtractionSolid
00657             ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
00658           break;
00659         case intersection:
00660           resultantSolid = new G4IntersectionSolid
00661             ("resultant_solid", pSol, fpClippingSolid, theAT.inverse());
00662           break;
00663         }
00664       }
00665 
00666       if (pSectionSolid) {
00667         resultantSolid = new G4IntersectionSolid
00668           ("sectioned_solid", pSol, pSectionSolid, theAT.inverse());
00669       }
00670 
00671       if (pCutawaySolid) {
00672         resultantSolid = new G4SubtractionSolid
00673           ("cutaway_solid", pSol, pCutawaySolid, theAT.inverse());
00674       }
00675 
00676       G4Polyhedron* tmpResultant = resultantSolid->GetPolyhedron();
00677       if (tmpResultant) resultant = *tmpResultant;
00678       else {
00679         if (fpMP->IsWarning())
00680           G4cout <<
00681             "WARNING: G4PhysicalVolumeModel::DescribeSolid: resultant polyhedron for"
00682             "\n  solid \"" << pSol->GetName() <<
00683             "\" not defined due to error during Boolean processing."
00684             "\n  Original will be drawn in red."
00685                  << G4endl;
00686         resultantVisAttribs.SetColour(G4Colour::Red());
00687       }
00688 
00689       delete resultantSolid;
00690 
00691       // Finally, force polyhedron drawing...
00692       resultant.SetVisAttributes(resultantVisAttribs);
00693       sceneHandler.BeginPrimitives(theAT);
00694       sceneHandler.AddPrimitive(resultant);
00695       sceneHandler.EndPrimitives();
00696     }
00697   }
00698   sceneHandler.PostAddSolid ();
00699 }

void G4PhysicalVolumeModel::DescribeYourselfTo ( G4VGraphicsScene  )  [virtual]

Implements G4VModel.

Reimplemented in G4LogicalVolumeModel.

Definition at line 136 of file G4PhysicalVolumeModel.cc.

References FatalException, and G4Exception().

Referenced by CalculateExtent(), G4LogicalVolumeModel::DescribeYourselfTo(), and Validate().

00137 {
00138   if (!fpTopPV) G4Exception
00139     ("G4PhysicalVolumeModel::DescribeYourselfTo",
00140      "modeling0012", FatalException, "No model.");
00141 
00142   if (!fpMP) G4Exception
00143     ("G4PhysicalVolumeModel::DescribeYourselfTo",
00144      "modeling0003", FatalException, "No modeling parameters.");
00145 
00146   // For safety...
00147   fCurrentDepth = 0;
00148 
00149   G4Transform3D startingTransformation = fTransform;
00150 
00151   VisitGeometryAndGetVisReps
00152     (fpTopPV,
00153      fRequestedDepth,
00154      startingTransformation,
00155      sceneHandler);
00156 
00157   // Clear data...
00158   fCurrentDepth     = 0;
00159   fpCurrentPV       = 0;
00160   fpCurrentLV       = 0;
00161   fpCurrentMaterial = 0;
00162   if (fFullPVPath.size() != fBaseFullPVPath.size()) {
00163     // They should be equal if pushing and popping is happening properly.
00164     G4Exception
00165     ("G4PhysicalVolumeModel::DescribeYourselfTo",
00166      "modeling0013",
00167      FatalException,
00168      "Path at start of modeling not equal to base path.  Something badly"
00169      "\nwrong.  Please contact visualisation coordinator.");
00170   }
00171   fDrawnPVPath.clear();
00172 }

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

Definition at line 754 of file G4PhysicalVolumeModel.cc.

References G4AttDefStore::GetInstance().

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

00755 {
00756     G4bool isNew;
00757     std::map<G4String,G4AttDef>* store
00758       = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew);
00759     if (isNew) {
00760       (*store)["PVPath"] =
00761         G4AttDef("PVPath","Physical Volume Path","Physics","","G4String");
00762       (*store)["LVol"] =
00763         G4AttDef("LVol","Logical Volume","Physics","","G4String");
00764       (*store)["Solid"] =
00765         G4AttDef("Solid","Solid Name","Physics","","G4String");
00766       (*store)["EType"] =
00767         G4AttDef("EType","Entity Type","Physics","","G4String");
00768       (*store)["DmpSol"] =
00769         G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String");
00770       (*store)["Trans"] =
00771         G4AttDef("Trans","Transformation of volume","Physics","","G4String");
00772       (*store)["Material"] =
00773         G4AttDef("Material","Material Name","Physics","","G4String");
00774       (*store)["Density"] =
00775         G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double");
00776       (*store)["State"] =
00777         G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String");
00778       (*store)["Radlen"] =
00779         G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double");
00780       (*store)["Region"] =
00781         G4AttDef("Region","Cuts Region","Physics","","G4String");
00782       (*store)["RootRegion"] =
00783         G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool");
00784     }
00785   return store;
00786 }

const G4VSolid* G4PhysicalVolumeModel::GetClippingSolid (  )  const [inline]

Definition at line 158 of file G4PhysicalVolumeModel.hh.

References fpClippingSolid.

00159   {return fpClippingSolid;}

G4int G4PhysicalVolumeModel::GetCurrentDepth (  )  const [inline]

Definition at line 161 of file G4PhysicalVolumeModel.hh.

References fCurrentDepth.

Referenced by G4HepRepSceneHandler::AddPrimitive(), G4HepRepSceneHandler::AddSolid(), and G4GMocrenFileSceneHandler::AddSolid().

00161 {return fCurrentDepth;}

G4String G4PhysicalVolumeModel::GetCurrentDescription (  )  const [virtual]

Reimplemented from G4VModel.

Definition at line 186 of file G4PhysicalVolumeModel.cc.

References GetCurrentTag().

00187 {
00188   return "G4PhysicalVolumeModel " + GetCurrentTag ();
00189 }

G4LogicalVolume* G4PhysicalVolumeModel::GetCurrentLV (  )  const [inline]

Definition at line 167 of file G4PhysicalVolumeModel.hh.

References fpCurrentLV.

Referenced by G4GMocrenFileSceneHandler::AddPrimitive(), G4HepRepSceneHandler::AddSolid(), and G4ASCIITreeSceneHandler::RequestPrimitives().

00167 {return fpCurrentLV;}

G4Material* G4PhysicalVolumeModel::GetCurrentMaterial (  )  const [inline]

Definition at line 170 of file G4PhysicalVolumeModel.hh.

References fpCurrentMaterial.

Referenced by G4HepRepSceneHandler::AddSolid(), G4GMocrenFileSceneHandler::AddSolid(), and G4ASCIITreeSceneHandler::RequestPrimitives().

00170 {return fpCurrentMaterial;}

G4VPhysicalVolume* G4PhysicalVolumeModel::GetCurrentPV (  )  const [inline]

Definition at line 164 of file G4PhysicalVolumeModel.hh.

References fpCurrentPV.

Referenced by G4GMocrenFileSceneHandler::AddPrimitive(), G4GMocrenFileSceneHandler::AddSolid(), and G4ASCIITreeSceneHandler::RequestPrimitives().

00164 {return fpCurrentPV;}

G4String G4PhysicalVolumeModel::GetCurrentTag (  )  const [virtual]

Reimplemented from G4VModel.

Definition at line 174 of file G4PhysicalVolumeModel.cc.

References G4VModel::fGlobalTag, and fpCurrentPV.

Referenced by GetCurrentDescription().

00175 {
00176   if (fpCurrentPV) {
00177     std::ostringstream o;
00178     o << fpCurrentPV -> GetCopyNo ();
00179     return fpCurrentPV -> GetName () + "." + o.str();
00180   }
00181   else {
00182     return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag;
00183   }
00184 }

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

Definition at line 180 of file G4PhysicalVolumeModel.hh.

References fDrawnPVPath.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4XXXSGSceneHandler::CreateCurrentItem(), G4VTreeSceneHandler::PreAddSolid(), and G4ASCIITreeSceneHandler::RequestPrimitives().

00181   {return fDrawnPVPath;}

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

Definition at line 173 of file G4PhysicalVolumeModel.hh.

References fFullPVPath.

00174   {return fFullPVPath;}

G4int G4PhysicalVolumeModel::GetRequestedDepth (  )  const [inline]

Definition at line 156 of file G4PhysicalVolumeModel.hh.

References fRequestedDepth.

00156 {return fRequestedDepth;}

G4VPhysicalVolume* G4PhysicalVolumeModel::GetTopPhysicalVolume (  )  const [inline]

Definition at line 154 of file G4PhysicalVolumeModel.hh.

References fpTopPV.

Referenced by G4GMocrenFileSceneHandler::AddSolid().

00154 {return fpTopPV;}

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

Definition at line 201 of file G4PhysicalVolumeModel.hh.

References fBaseFullPVPath, and fFullPVPath.

00202                      {
00203       fBaseFullPVPath = baseFullPVPath;
00204       fFullPVPath = baseFullPVPath;
00205   }

void G4PhysicalVolumeModel::SetClippingMode ( ClippingMode  mode  )  [inline]

Definition at line 215 of file G4PhysicalVolumeModel.hh.

References fClippingMode.

00215                                            {
00216     fClippingMode = mode;
00217   }

void G4PhysicalVolumeModel::SetClippingSolid ( G4VSolid pClippingSolid  )  [inline]

Definition at line 211 of file G4PhysicalVolumeModel.hh.

References fpClippingSolid.

00211                                                    {
00212     fpClippingSolid = pClippingSolid;
00213   }

void G4PhysicalVolumeModel::SetRequestedDepth ( G4int  requestedDepth  )  [inline]

Definition at line 207 of file G4PhysicalVolumeModel.hh.

References fRequestedDepth.

00207                                                 {
00208     fRequestedDepth = requestedDepth;
00209   }

G4bool G4PhysicalVolumeModel::Validate ( G4bool  warn  )  [virtual]

Reimplemented from G4VModel.

Reimplemented in G4LogicalVolumeModel.

Definition at line 701 of file G4PhysicalVolumeModel.cc.

References CalculateExtent(), DescribeYourselfTo(), G4VModel::fpMP, fpTopPV, fTopPVCopyNo, fTopPVName, G4cout, G4endl, G4ModelingParameters::GetDefaultVisAttributes(), G4PhysicalVolumeSearchScene::GetFoundVolume(), G4TransportationManager::GetNoWorlds(), G4TransportationManager::GetTransportationManager(), G4TransportationManager::GetWorldsIterator(), G4ModelingParameters::SetDefaultVisAttributes(), and G4VModel::SetModelingParameters().

00702 {
00703   G4TransportationManager* transportationManager =
00704     G4TransportationManager::GetTransportationManager ();
00705 
00706   size_t nWorlds = transportationManager->GetNoWorlds();
00707 
00708   G4bool found = false;
00709 
00710   std::vector<G4VPhysicalVolume*>::iterator iterWorld =
00711     transportationManager->GetWorldsIterator();
00712   for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
00713     G4VPhysicalVolume* world = (*iterWorld);
00714     // The idea now is to seek a PV with the same name and copy no
00715     // in the hope it's the same one!!
00716     G4PhysicalVolumeModel searchModel (world);
00717     G4PhysicalVolumeSearchScene searchScene
00718       (&searchModel, fTopPVName, fTopPVCopyNo);
00719     G4ModelingParameters mp;  // Default modeling parameters for this search.
00720     mp.SetDefaultVisAttributes(fpMP? fpMP->GetDefaultVisAttributes(): 0);
00721     searchModel.SetModelingParameters (&mp);
00722     searchModel.DescribeYourselfTo (searchScene);
00723     G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume ();
00724     if (foundVolume) {
00725       if (foundVolume != fpTopPV && warn) {
00726         G4cout <<
00727           "G4PhysicalVolumeModel::Validate(): A volume of the same name and"
00728           "\n  copy number (\""
00729                << fTopPVName << "\", copy " << fTopPVCopyNo
00730                << ") still exists and is being used."
00731           "\n  But it is not the same volume you originally specified"
00732           "\n  in /vis/scene/add/."
00733                << G4endl;
00734       }
00735       fpTopPV = foundVolume;
00736       CalculateExtent ();
00737       found = true;
00738     }
00739   }
00740   if (found) return true;
00741   else {
00742     if (warn) {
00743       G4cout <<
00744         "G4PhysicalVolumeModel::Validate(): No volume of name and"
00745         "\n  copy number (\""
00746              << fTopPVName << "\", copy " << fTopPVCopyNo
00747              << ") exists."
00748              << G4endl;
00749     }
00750     return false;
00751   }
00752 }

void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps ( G4VPhysicalVolume ,
G4int  requestedDepth,
const G4Transform3D ,
G4VGraphicsScene  
) [protected]

Definition at line 192 of file G4PhysicalVolumeModel.cc.

References G4cout, G4endl, G4VSolid::GetEntityType(), G4VSolid::GetName(), kPhi, kRho, kXAxis, kYAxis, kZAxis, and CLHEP::detail::n.

00196 {
00197   // Visits geometry structure to a given depth (requestedDepth), starting
00198   //   at given physical volume with given starting transformation and
00199   //   describes volumes to the scene handler.
00200   // requestedDepth < 0 (default) implies full visit.
00201   // theAT is the Accumulated Transformation.
00202 
00203   // Find corresponding logical volume and (later) solid, storing in
00204   // local variables to preserve re-entrancy.
00205   G4LogicalVolume* pLV  = pVPV -> GetLogicalVolume ();
00206 
00207   G4VSolid* pSol;
00208   G4Material* pMaterial;
00209 
00210   if (!(pVPV -> IsReplicated ())) {
00211     // Non-replicated physical volume.
00212     pSol = pLV -> GetSolid ();
00213     pMaterial = pLV -> GetMaterial ();
00214     DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
00215                         theAT, sceneHandler);
00216   }
00217   else {
00218     // Replicated or parametrised physical volume.
00219     EAxis axis;
00220     G4int nReplicas;
00221     G4double width;
00222     G4double offset;
00223     G4bool consuming;
00224     pVPV -> GetReplicationData (axis, nReplicas, width,  offset, consuming);
00225     if (fCurrentDepth == 0) nReplicas = 1;  // Just draw first
00226     G4VPVParameterisation* pP = pVPV -> GetParameterisation ();
00227     if (pP) {  // Parametrised volume.
00228       for (int n = 0; n < nReplicas; n++) {
00229         pSol = pP -> ComputeSolid (n, pVPV);
00230         pP -> ComputeTransformation (n, pVPV);
00231         pSol -> ComputeDimensions (pP, n, pVPV);
00232         pVPV -> SetCopyNo (n);
00233         // Create a touchable of current parent for ComputeMaterial.
00234         // fFullPVPath has not been updated yet so at this point it
00235         // corresponds to the parent.
00236         G4PhysicalVolumeModelTouchable parentTouchable(fFullPVPath);
00237         pMaterial = pP -> ComputeMaterial (n, pVPV, &parentTouchable);
00238         DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
00239                             theAT, sceneHandler);
00240       }
00241     }
00242     else {  // Plain replicated volume.  From geometry_guide.txt...
00243       // The replica's positions are claculated by means of a linear formula.
00244       // Replication may occur along:
00245       // 
00246       // o Cartesian axes (kXAxis,kYAxis,kZAxis)
00247       // 
00248       //   The replications, of specified width have coordinates of
00249       //   form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1
00250       //   for the case of kXAxis, and are unrotated.
00251       // 
00252       // o Radial axis (cylindrical polar) (kRho)
00253       // 
00254       //   The replications are cons/tubs sections, centred on the origin
00255       //   and are unrotated.
00256       //   They have radii of width*n+offset to width*(n+1)+offset
00257       //                      where n=0..nReplicas-1
00258       // 
00259       // o Phi axis (cylindrical polar) (kPhi)
00260       //   The replications are `phi sections' or wedges, and of cons/tubs form
00261       //   They have phi of offset+n*width to offset+(n+1)*width where
00262       //   n=0..nReplicas-1
00263       // 
00264       pSol = pLV -> GetSolid ();
00265       pMaterial = pLV -> GetMaterial ();
00266       G4ThreeVector originalTranslation = pVPV -> GetTranslation ();
00267       G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation ();
00268       G4double originalRMin = 0., originalRMax = 0.;
00269       if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
00270         originalRMin = ((G4Tubs*)pSol)->GetInnerRadius();
00271         originalRMax = ((G4Tubs*)pSol)->GetOuterRadius();
00272       }
00273       G4bool visualisable = true;
00274       for (int n = 0; n < nReplicas; n++) {
00275         G4ThreeVector translation;  // Null.
00276         G4RotationMatrix rotation;  // Null - life long enough for visualizing.
00277         G4RotationMatrix* pRotation = 0;
00278         switch (axis) {
00279         default:
00280         case kXAxis:
00281           translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0);
00282           break;
00283         case kYAxis:
00284           translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0);
00285           break;
00286         case kZAxis:
00287           translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width);
00288           break;
00289         case kRho:
00290           if (pSol->GetEntityType() == "G4Tubs") {
00291             ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset);
00292             ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset);
00293           } else {
00294             if (fpMP->IsWarning())
00295               G4cout <<
00296                 "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:"
00297                 "\n  built-in replicated volumes replicated in radius for "
00298                      << pSol->GetEntityType() <<
00299                 "-type\n  solids (your solid \""
00300                      << pSol->GetName() <<
00301                 "\") are not visualisable."
00302                      << G4endl;
00303             visualisable = false;
00304           }
00305           break;
00306         case kPhi:
00307           rotation.rotateZ (-(offset+(n+0.5)*width));
00308           // Minus Sign because for the physical volume we need the
00309           // coordinate system rotation.
00310           pRotation = &rotation;
00311           break;
00312         } 
00313         pVPV -> SetTranslation (translation);
00314         pVPV -> SetRotation    (pRotation);
00315         pVPV -> SetCopyNo (n);
00316         if (visualisable) {
00317           DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial,
00318                             theAT, sceneHandler);
00319         }
00320       }
00321       // Restore originals...
00322       pVPV -> SetTranslation (originalTranslation);
00323       pVPV -> SetRotation    (pOriginalRotation);
00324       if (axis == kRho && pSol->GetEntityType() == "G4Tubs") {
00325         ((G4Tubs*)pSol)->SetInnerRadius(originalRMin);
00326         ((G4Tubs*)pSol)->SetOuterRadius(originalRMax);
00327       }
00328     }
00329   }
00330 
00331   return;
00332 }


Field Documentation

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

Definition at line 260 of file G4PhysicalVolumeModel.hh.

Referenced by SetBaseFullPVPath().

ClippingMode G4PhysicalVolumeModel::fClippingMode [protected]

Definition at line 265 of file G4PhysicalVolumeModel.hh.

Referenced by SetClippingMode().

G4int G4PhysicalVolumeModel::fCurrentDepth [protected]

Definition at line 255 of file G4PhysicalVolumeModel.hh.

Referenced by GetCurrentDepth().

G4bool G4PhysicalVolumeModel::fCurtailDescent [protected]

Definition at line 263 of file G4PhysicalVolumeModel.hh.

Referenced by CurtailDescent().

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

Definition at line 262 of file G4PhysicalVolumeModel.hh.

Referenced by GetDrawnPVPath().

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

Definition at line 261 of file G4PhysicalVolumeModel.hh.

Referenced by CreateCurrentAttValues(), GetFullPVPath(), and SetBaseFullPVPath().

G4VSolid* G4PhysicalVolumeModel::fpClippingSolid [protected]

Definition at line 264 of file G4PhysicalVolumeModel.hh.

Referenced by GetClippingSolid(), SetClippingSolid(), and ~G4PhysicalVolumeModel().

G4LogicalVolume* G4PhysicalVolumeModel::fpCurrentLV [protected]

Definition at line 257 of file G4PhysicalVolumeModel.hh.

Referenced by CreateCurrentAttValues(), and GetCurrentLV().

G4Material* G4PhysicalVolumeModel::fpCurrentMaterial [protected]

Definition at line 258 of file G4PhysicalVolumeModel.hh.

Referenced by CreateCurrentAttValues(), and GetCurrentMaterial().

G4VPhysicalVolume* G4PhysicalVolumeModel::fpCurrentPV [protected]

Definition at line 256 of file G4PhysicalVolumeModel.hh.

Referenced by GetCurrentPV(), and GetCurrentTag().

G4Transform3D* G4PhysicalVolumeModel::fpCurrentTransform [protected]

Definition at line 259 of file G4PhysicalVolumeModel.hh.

Referenced by CreateCurrentAttValues().

G4VPhysicalVolume* G4PhysicalVolumeModel::fpTopPV [protected]

Definition at line 249 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent(), GetTopPhysicalVolume(), and Validate().

G4int G4PhysicalVolumeModel::fRequestedDepth [protected]

Definition at line 252 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent(), GetRequestedDepth(), and SetRequestedDepth().

G4int G4PhysicalVolumeModel::fTopPVCopyNo [protected]

Definition at line 251 of file G4PhysicalVolumeModel.hh.

Referenced by Validate().

G4String G4PhysicalVolumeModel::fTopPVName [protected]

Definition at line 250 of file G4PhysicalVolumeModel.hh.

Referenced by Validate().

G4bool G4PhysicalVolumeModel::fUseFullExtent [protected]

Definition at line 254 of file G4PhysicalVolumeModel.hh.

Referenced by CalculateExtent().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:53 2013 for Geant4 by  doxygen 1.4.7