G4PhysicalVolumeModel.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // 
00030 // John Allison  31st December 1997.
00031 // Model for physical volumes.
00032 
00033 #include "G4PhysicalVolumeModel.hh"
00034 
00035 #include "G4ModelingParameters.hh"
00036 #include "G4VGraphicsScene.hh"
00037 #include "G4VPhysicalVolume.hh"
00038 #include "G4VPVParameterisation.hh"
00039 #include "G4LogicalVolume.hh"
00040 #include "G4VSolid.hh"
00041 #include "G4SubtractionSolid.hh"
00042 #include "G4IntersectionSolid.hh"
00043 #include "G4Material.hh"
00044 #include "G4VisAttributes.hh"
00045 #include "G4BoundingSphereScene.hh"
00046 #include "G4PhysicalVolumeSearchScene.hh"
00047 #include "G4TransportationManager.hh"
00048 #include "G4Polyhedron.hh"
00049 #include "HepPolyhedronProcessor.h"
00050 #include "G4AttDefStore.hh"
00051 #include "G4AttDef.hh"
00052 #include "G4AttValue.hh"
00053 #include "G4UnitsTable.hh"
00054 #include "G4Vector3D.hh"
00055 
00056 #include <sstream>
00057 
00058 G4PhysicalVolumeModel::G4PhysicalVolumeModel
00059 (G4VPhysicalVolume*          pVPV,
00060  G4int                       requestedDepth,
00061  const G4Transform3D& modelTransformation,
00062  const G4ModelingParameters* pMP,
00063  G4bool useFullExtent):
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 }
00090 
00091 G4PhysicalVolumeModel::~G4PhysicalVolumeModel ()
00092 {
00093   delete fpClippingSolid;
00094 }
00095 
00096 void G4PhysicalVolumeModel::CalculateExtent ()
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 }
00134 
00135 void G4PhysicalVolumeModel::DescribeYourselfTo
00136 (G4VGraphicsScene& sceneHandler)
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 }
00173 
00174 G4String G4PhysicalVolumeModel::GetCurrentTag () const
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 }
00185 
00186 G4String G4PhysicalVolumeModel::GetCurrentDescription () const
00187 {
00188   return "G4PhysicalVolumeModel " + GetCurrentTag ();
00189 }
00190 
00191 void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps
00192 (G4VPhysicalVolume* pVPV,
00193  G4int requestedDepth,
00194  const G4Transform3D& theAT,
00195  G4VGraphicsScene& sceneHandler)
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 }
00333 
00334 void G4PhysicalVolumeModel::DescribeAndDescend
00335 (G4VPhysicalVolume* pVPV,
00336  G4int requestedDepth,
00337  G4LogicalVolume* pLV,
00338  G4VSolid* pSol,
00339  G4Material* pMaterial,
00340  const G4Transform3D& theAT,
00341  G4VGraphicsScene& sceneHandler)
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 }
00607 
00608 void G4PhysicalVolumeModel::DescribeSolid
00609 (const G4Transform3D& theAT,
00610  G4VSolid* pSol,
00611  const G4VisAttributes* pVisAttribs,
00612  G4VGraphicsScene& sceneHandler)
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 }
00700 
00701 G4bool G4PhysicalVolumeModel::Validate (G4bool warn)
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 }
00753 
00754 const std::map<G4String,G4AttDef>* G4PhysicalVolumeModel::GetAttDefs() const
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 }
00787 
00788 #include <iomanip>
00789 
00790 static std::ostream& operator<< (std::ostream& o, const G4Transform3D t)
00791 {
00792   using namespace std;
00793 
00794   G4Scale3D sc;
00795   G4Rotate3D r;
00796   G4Translate3D tl;
00797   t.getDecomposition(sc, r, tl);
00798 
00799   const int w = 10;
00800 
00801   // Transformation itself
00802   o << setw(w) << t.xx() << setw(w) << t.xy() << setw(w) << t.xz() << setw(w) << t.dx() << endl;
00803   o << setw(w) << t.yx() << setw(w) << t.yy() << setw(w) << t.yz() << setw(w) << t.dy() << endl;
00804   o << setw(w) << t.zx() << setw(w) << t.zy() << setw(w) << t.zz() << setw(w) << t.dz() << endl;
00805 
00806   // Translation
00807   o << "= translation:" << endl;
00808   o << setw(w) << tl.dx() << setw(w) << tl.dy() << setw(w) << tl.dz() << endl;
00809 
00810   // Rotation
00811   o << "* rotation:" << endl;
00812   o << setw(w) << r.xx() << setw(w) << r.xy() << setw(w) << r.xz() << endl;
00813   o << setw(w) << r.yx() << setw(w) << r.yy() << setw(w) << r.yz() << endl;
00814   o << setw(w) << r.zx() << setw(w) << r.zy() << setw(w) << r.zz() << endl;
00815 
00816   // Scale
00817   o << "* scale:" << endl;
00818   o << setw(w) << sc.xx() << setw(w) << sc.yy() << setw(w) << sc.zz() << endl;
00819 
00820   // Transformed axes
00821   o << "Transformed axes:" << endl;
00822   o << "x': " << r * G4Vector3D(1., 0., 0.) << endl;
00823   o << "y': " << r * G4Vector3D(0., 1., 0.) << endl;
00824   o << "z': " << r * G4Vector3D(0., 0., 1.) << endl;
00825 
00826   return o;
00827 }
00828 
00829 std::vector<G4AttValue>* G4PhysicalVolumeModel::CreateCurrentAttValues() const
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 }
00873 
00874 G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator<
00875   (const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID& right) const
00876 {
00877   if (fpPV < right.fpPV) return true;
00878   if (fpPV == right.fpPV) {
00879     if (fCopyNo < right.fCopyNo) return true;
00880     if (fCopyNo == right.fCopyNo)
00881       return fNonCulledDepth < right.fNonCulledDepth;
00882   }
00883   return false;
00884 }
00885 
00886 std::ostream& operator<<
00887   (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID node)
00888 {
00889   G4VPhysicalVolume* pPV = node.GetPhysicalVolume();
00890   if (pPV) {
00891     os << pPV->GetName()
00892        << ':' << node.GetCopyNo()
00893        << '[' << node.GetNonCulledDepth() << ']'
00894        << ':' << node.GetTransform();
00895     if (!node.GetDrawn()) os << "  Not "; os << "drawn";
00896   } else {
00897     os << "Null node";
00898   }
00899   return os;
00900 }
00901 
00902 G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::G4PhysicalVolumeModelTouchable
00903 (const std::vector<G4PhysicalVolumeNodeID>& fullPVPath):
00904   fFullPVPath(fullPVPath) {}
00905 
00906 const G4ThreeVector& G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetTranslation(G4int depth) const
00907 {
00908   size_t i = fFullPVPath.size() - depth - 1;
00909   if (i >= fFullPVPath.size()) {
00910     G4Exception("G4PhysicalVolumeModelTouchable::GetTranslation",
00911                 "modeling0005",
00912                 FatalErrorInArgument,
00913                 "Index out of range. Asking for non-existent depth");
00914   }
00915   static G4ThreeVector tempTranslation;
00916   tempTranslation = fFullPVPath[i].GetTransform().getTranslation();
00917   return tempTranslation;
00918 }
00919 
00920 const G4RotationMatrix* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetRotation(G4int depth) const
00921 {
00922   size_t i = fFullPVPath.size() - depth - 1;
00923   if (i >= fFullPVPath.size()) {
00924     G4Exception("G4PhysicalVolumeModelTouchable::GetRotation",
00925                 "modeling0006",
00926                 FatalErrorInArgument,
00927                 "Index out of range. Asking for non-existent depth");
00928   }
00929   static G4RotationMatrix tempRotation;
00930   tempRotation = fFullPVPath[i].GetTransform().getRotation();
00931   return &tempRotation;
00932 }
00933 
00934 G4VPhysicalVolume* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetVolume(G4int depth) const
00935 {
00936   size_t i = fFullPVPath.size() - depth - 1;
00937   if (i >= fFullPVPath.size()) {
00938     G4Exception("G4PhysicalVolumeModelTouchable::GetVolume",
00939                 "modeling0007",
00940                 FatalErrorInArgument,
00941                 "Index out of range. Asking for non-existent depth");
00942   }
00943   return fFullPVPath[i].GetPhysicalVolume();
00944 }
00945 
00946 G4VSolid* G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetSolid(G4int depth) const
00947 {
00948   size_t i = fFullPVPath.size() - depth - 1;
00949   if (i >= fFullPVPath.size()) {
00950     G4Exception("G4PhysicalVolumeModelTouchable::GetSolid",
00951                 "modeling0008",
00952                 FatalErrorInArgument,
00953                 "Index out of range. Asking for non-existent depth");
00954   }
00955   return fFullPVPath[i].GetPhysicalVolume()->GetLogicalVolume()->GetSolid();
00956 }
00957 
00958 G4int G4PhysicalVolumeModel::G4PhysicalVolumeModelTouchable::GetReplicaNumber(G4int depth) const
00959 {
00960   size_t i = fFullPVPath.size() - depth - 1;
00961   if (i >= fFullPVPath.size()) {
00962     G4Exception("G4PhysicalVolumeModelTouchable::GetReplicaNumber",
00963                 "modeling0009",
00964                 FatalErrorInArgument,
00965                 "Index out of range. Asking for non-existent depth");
00966   }
00967   return fFullPVPath[i].GetCopyNo();
00968 }

Generated on Mon May 27 17:49:19 2013 for Geant4 by  doxygen 1.4.7