[834] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // |
---|
| 27 | // $Id: G4PhysicalVolumeModel.cc,v 1.63 2007/11/10 14:56:36 allison Exp $ |
---|
[850] | 28 | // GEANT4 tag $Name: HEAD $ |
---|
[834] | 29 | // |
---|
| 30 | // |
---|
| 31 | // John Allison 31st December 1997. |
---|
| 32 | // Model for physical volumes. |
---|
| 33 | |
---|
| 34 | #include "G4PhysicalVolumeModel.hh" |
---|
| 35 | |
---|
| 36 | #include "G4ModelingParameters.hh" |
---|
| 37 | #include "G4VGraphicsScene.hh" |
---|
| 38 | #include "G4VPhysicalVolume.hh" |
---|
| 39 | #include "G4VPVParameterisation.hh" |
---|
| 40 | #include "G4LogicalVolume.hh" |
---|
| 41 | #include "G4VSolid.hh" |
---|
| 42 | #include "G4Material.hh" |
---|
| 43 | #include "G4VisAttributes.hh" |
---|
| 44 | #include "G4BoundingSphereScene.hh" |
---|
| 45 | #include "G4PhysicalVolumeSearchScene.hh" |
---|
| 46 | #include "G4TransportationManager.hh" |
---|
| 47 | #include "G4Polyhedron.hh" |
---|
| 48 | #include "G4AttDefStore.hh" |
---|
| 49 | #include "G4AttDef.hh" |
---|
| 50 | #include "G4AttValue.hh" |
---|
| 51 | #include "G4UnitsTable.hh" |
---|
| 52 | #include "G4Vector3D.hh" |
---|
| 53 | |
---|
| 54 | #include <sstream> |
---|
| 55 | |
---|
| 56 | G4bool G4PhysicalVolumeModel::G4PhysicalVolumeNodeID::operator< |
---|
| 57 | (const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID& right) const |
---|
| 58 | { |
---|
| 59 | if (fpPV < right.fpPV) return true; |
---|
| 60 | if (fpPV == right.fpPV) { |
---|
| 61 | if (fCopyNo < right.fCopyNo) return true; |
---|
| 62 | if (fCopyNo == right.fCopyNo) |
---|
| 63 | return fNonCulledDepth < right.fNonCulledDepth; |
---|
| 64 | } |
---|
| 65 | return false; |
---|
| 66 | } |
---|
| 67 | |
---|
| 68 | std::ostream& operator<< |
---|
| 69 | (std::ostream& os, const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID node) |
---|
| 70 | { |
---|
| 71 | G4VPhysicalVolume* pPV = node.GetPhysicalVolume(); |
---|
| 72 | if (pPV) { |
---|
| 73 | os << pPV->GetName() |
---|
| 74 | << ':' << node.GetCopyNo() |
---|
| 75 | << '[' << node.GetNonCulledDepth() << ']'; |
---|
| 76 | } else { |
---|
| 77 | os << "Null node"; |
---|
| 78 | } |
---|
| 79 | return os; |
---|
| 80 | } |
---|
| 81 | |
---|
| 82 | G4PhysicalVolumeModel::G4PhysicalVolumeModel |
---|
| 83 | (G4VPhysicalVolume* pVPV, |
---|
| 84 | G4int requestedDepth, |
---|
| 85 | const G4Transform3D& modelTransformation, |
---|
| 86 | const G4ModelingParameters* pMP, |
---|
| 87 | G4bool useFullExtent): |
---|
| 88 | G4VModel (modelTransformation, pMP), |
---|
| 89 | fpTopPV (pVPV), |
---|
| 90 | fTopPVName (pVPV -> GetName ()), |
---|
| 91 | fTopPVCopyNo (pVPV -> GetCopyNo ()), |
---|
| 92 | fRequestedDepth (requestedDepth), |
---|
| 93 | fUseFullExtent (useFullExtent), |
---|
| 94 | fCurrentDepth (0), |
---|
| 95 | fpCurrentPV (0), |
---|
| 96 | fpCurrentLV (0), |
---|
| 97 | fpCurrentMaterial (0), |
---|
| 98 | fpCurrentTransform (0), |
---|
| 99 | fCurtailDescent (false), |
---|
| 100 | fpClippingPolyhedron (0), |
---|
| 101 | fClippingMode (subtraction) |
---|
| 102 | { |
---|
| 103 | std::ostringstream o; |
---|
| 104 | o << fpTopPV -> GetCopyNo (); |
---|
| 105 | fGlobalTag = fpTopPV -> GetName () + "." + o.str(); |
---|
| 106 | fGlobalDescription = "G4PhysicalVolumeModel " + fGlobalTag; |
---|
| 107 | |
---|
| 108 | CalculateExtent (); |
---|
| 109 | } |
---|
| 110 | |
---|
| 111 | G4PhysicalVolumeModel::~G4PhysicalVolumeModel () |
---|
| 112 | { |
---|
| 113 | delete fpClippingPolyhedron; |
---|
| 114 | } |
---|
| 115 | |
---|
| 116 | void G4PhysicalVolumeModel::CalculateExtent () |
---|
| 117 | { |
---|
| 118 | if (fUseFullExtent) { |
---|
| 119 | fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent (); |
---|
| 120 | } |
---|
| 121 | else { |
---|
| 122 | G4BoundingSphereScene bsScene(this); |
---|
| 123 | const G4int tempRequestedDepth = fRequestedDepth; |
---|
| 124 | fRequestedDepth = -1; // Always search to all depths to define extent. |
---|
| 125 | const G4ModelingParameters* tempMP = fpMP; |
---|
| 126 | G4ModelingParameters mParams |
---|
| 127 | (0, // No default vis attributes needed. |
---|
| 128 | G4ModelingParameters::wf, // wireframe (not relevant for this). |
---|
| 129 | true, // Global culling. |
---|
| 130 | true, // Cull invisible volumes. |
---|
| 131 | false, // Density culling. |
---|
| 132 | 0., // Density (not relevant if density culling false). |
---|
| 133 | true, // Cull daughters of opaque mothers. |
---|
| 134 | 24); // No of sides (not relevant for this operation). |
---|
| 135 | fpMP = &mParams; |
---|
| 136 | DescribeYourselfTo (bsScene); |
---|
| 137 | G4double radius = bsScene.GetRadius(); |
---|
| 138 | if (radius < 0.) { // Nothing in the scene. |
---|
| 139 | fExtent = fpTopPV -> GetLogicalVolume () -> GetSolid () -> GetExtent (); |
---|
| 140 | } else { |
---|
| 141 | // Transform back to coordinates relative to the top |
---|
| 142 | // transformation, which is in G4VModel::fTransform. This makes |
---|
| 143 | // it conform to all models, which are defined by a |
---|
| 144 | // transformation and an extent relative to that |
---|
| 145 | // transformation... |
---|
| 146 | G4Point3D centre = bsScene.GetCentre(); |
---|
| 147 | centre.transform(fTransform.inverse()); |
---|
| 148 | fExtent = G4VisExtent(centre, radius); |
---|
| 149 | } |
---|
| 150 | fpMP = tempMP; |
---|
| 151 | fRequestedDepth = tempRequestedDepth; |
---|
| 152 | } |
---|
| 153 | } |
---|
| 154 | |
---|
| 155 | void G4PhysicalVolumeModel::DescribeYourselfTo |
---|
| 156 | (G4VGraphicsScene& sceneHandler) |
---|
| 157 | { |
---|
| 158 | if (!fpMP) G4Exception |
---|
| 159 | ("G4PhysicalVolumeModel::DescribeYourselfTo: No modeling parameters."); |
---|
| 160 | |
---|
| 161 | // For safety... |
---|
| 162 | fCurrentDepth = 0; |
---|
| 163 | |
---|
| 164 | G4Transform3D startingTransformation = fTransform; |
---|
| 165 | |
---|
| 166 | VisitGeometryAndGetVisReps |
---|
| 167 | (fpTopPV, |
---|
| 168 | fRequestedDepth, |
---|
| 169 | startingTransformation, |
---|
| 170 | sceneHandler); |
---|
| 171 | |
---|
| 172 | // Clear data... |
---|
| 173 | fCurrentDepth = 0; |
---|
| 174 | fpCurrentPV = 0; |
---|
| 175 | fpCurrentLV = 0; |
---|
| 176 | fpCurrentMaterial = 0; |
---|
| 177 | fFullPVPath.clear(); |
---|
| 178 | fDrawnPVPath.clear(); |
---|
| 179 | } |
---|
| 180 | |
---|
| 181 | G4String G4PhysicalVolumeModel::GetCurrentTag () const |
---|
| 182 | { |
---|
| 183 | if (fpCurrentPV) { |
---|
| 184 | std::ostringstream o; |
---|
| 185 | o << fpCurrentPV -> GetCopyNo (); |
---|
| 186 | return fpCurrentPV -> GetName () + "." + o.str(); |
---|
| 187 | } |
---|
| 188 | else { |
---|
| 189 | return "WARNING: NO CURRENT VOLUME - global tag is " + fGlobalTag; |
---|
| 190 | } |
---|
| 191 | } |
---|
| 192 | |
---|
| 193 | G4String G4PhysicalVolumeModel::GetCurrentDescription () const |
---|
| 194 | { |
---|
| 195 | return "G4PhysicalVolumeModel " + GetCurrentTag (); |
---|
| 196 | } |
---|
| 197 | |
---|
| 198 | void G4PhysicalVolumeModel::VisitGeometryAndGetVisReps |
---|
| 199 | (G4VPhysicalVolume* pVPV, |
---|
| 200 | G4int requestedDepth, |
---|
| 201 | const G4Transform3D& theAT, |
---|
| 202 | G4VGraphicsScene& sceneHandler) |
---|
| 203 | { |
---|
| 204 | // Visits geometry structure to a given depth (requestedDepth), starting |
---|
| 205 | // at given physical volume with given starting transformation and |
---|
| 206 | // describes volumes to the scene handler. |
---|
| 207 | // requestedDepth < 0 (default) implies full visit. |
---|
| 208 | // theAT is the Accumulated Transformation. |
---|
| 209 | |
---|
| 210 | // Find corresponding logical volume and (later) solid, storing in |
---|
| 211 | // local variables to preserve re-entrancy. |
---|
| 212 | G4LogicalVolume* pLV = pVPV -> GetLogicalVolume (); |
---|
| 213 | |
---|
| 214 | G4VSolid* pSol; |
---|
| 215 | G4Material* pMaterial; |
---|
| 216 | |
---|
| 217 | if (!(pVPV -> IsReplicated ())) { |
---|
| 218 | // Non-replicated physical volume. |
---|
| 219 | pSol = pLV -> GetSolid (); |
---|
| 220 | pMaterial = pLV -> GetMaterial (); |
---|
| 221 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
| 222 | theAT, sceneHandler); |
---|
| 223 | } |
---|
| 224 | else { |
---|
| 225 | // Replicated or parametrised physical volume. |
---|
| 226 | EAxis axis; |
---|
| 227 | G4int nReplicas; |
---|
| 228 | G4double width; |
---|
| 229 | G4double offset; |
---|
| 230 | G4bool consuming; |
---|
| 231 | pVPV -> GetReplicationData (axis, nReplicas, width, offset, consuming); |
---|
| 232 | G4VPVParameterisation* pP = pVPV -> GetParameterisation (); |
---|
| 233 | if (pP) { // Parametrised volume. |
---|
| 234 | for (int n = 0; n < nReplicas; n++) { |
---|
| 235 | pSol = pP -> ComputeSolid (n, pVPV); |
---|
| 236 | pMaterial = pP -> ComputeMaterial (n, pVPV); |
---|
| 237 | pP -> ComputeTransformation (n, pVPV); |
---|
| 238 | pSol -> ComputeDimensions (pP, n, pVPV); |
---|
| 239 | pVPV -> SetCopyNo (n); |
---|
| 240 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
| 241 | theAT, sceneHandler); |
---|
| 242 | } |
---|
| 243 | } |
---|
| 244 | else { // Plain replicated volume. From geometry_guide.txt... |
---|
| 245 | // The replica's positions are claculated by means of a linear formula. |
---|
| 246 | // Replication may occur along: |
---|
| 247 | // |
---|
| 248 | // o Cartesian axes (kXAxis,kYAxis,kZAxis) |
---|
| 249 | // |
---|
| 250 | // The replications, of specified width have coordinates of |
---|
| 251 | // form (-width*(nReplicas-1)*0.5+n*width,0,0) where n=0.. nReplicas-1 |
---|
| 252 | // for the case of kXAxis, and are unrotated. |
---|
| 253 | // |
---|
| 254 | // o Radial axis (cylindrical polar) (kRho) |
---|
| 255 | // |
---|
| 256 | // The replications are cons/tubs sections, centred on the origin |
---|
| 257 | // and are unrotated. |
---|
| 258 | // They have radii of width*n+offset to width*(n+1)+offset |
---|
| 259 | // where n=0..nReplicas-1 |
---|
| 260 | // |
---|
| 261 | // o Phi axis (cylindrical polar) (kPhi) |
---|
| 262 | // The replications are `phi sections' or wedges, and of cons/tubs form |
---|
| 263 | // They have phi of offset+n*width to offset+(n+1)*width where |
---|
| 264 | // n=0..nReplicas-1 |
---|
| 265 | // |
---|
| 266 | pSol = pLV -> GetSolid (); |
---|
| 267 | pMaterial = pLV -> GetMaterial (); |
---|
| 268 | G4ThreeVector originalTranslation = pVPV -> GetTranslation (); |
---|
| 269 | G4RotationMatrix* pOriginalRotation = pVPV -> GetRotation (); |
---|
| 270 | G4double originalRMin = 0., originalRMax = 0.; |
---|
| 271 | if (axis == kRho && pSol->GetEntityType() == "G4Tubs") { |
---|
| 272 | originalRMin = ((G4Tubs*)pSol)->GetInnerRadius(); |
---|
| 273 | originalRMax = ((G4Tubs*)pSol)->GetOuterRadius(); |
---|
| 274 | } |
---|
| 275 | G4bool visualisable = true; |
---|
| 276 | for (int n = 0; n < nReplicas; n++) { |
---|
| 277 | G4ThreeVector translation; // Null. |
---|
| 278 | G4RotationMatrix rotation; // Null - life long enough for visualizing. |
---|
| 279 | G4RotationMatrix* pRotation = 0; |
---|
| 280 | switch (axis) { |
---|
| 281 | default: |
---|
| 282 | case kXAxis: |
---|
| 283 | translation = G4ThreeVector (-width*(nReplicas-1)*0.5+n*width,0,0); |
---|
| 284 | break; |
---|
| 285 | case kYAxis: |
---|
| 286 | translation = G4ThreeVector (0,-width*(nReplicas-1)*0.5+n*width,0); |
---|
| 287 | break; |
---|
| 288 | case kZAxis: |
---|
| 289 | translation = G4ThreeVector (0,0,-width*(nReplicas-1)*0.5+n*width); |
---|
| 290 | break; |
---|
| 291 | case kRho: |
---|
| 292 | if (pSol->GetEntityType() == "G4Tubs") { |
---|
| 293 | ((G4Tubs*)pSol)->SetInnerRadius(width*n+offset); |
---|
| 294 | ((G4Tubs*)pSol)->SetOuterRadius(width*(n+1)+offset); |
---|
| 295 | } else { |
---|
| 296 | if (fpMP->IsWarning()) |
---|
| 297 | G4cout << |
---|
| 298 | "G4PhysicalVolumeModel::VisitGeometryAndGetVisReps: WARNING:" |
---|
| 299 | "\n built-in replicated volumes replicated in radius for " |
---|
| 300 | << pSol->GetEntityType() << |
---|
| 301 | "-type\n solids (your solid \"" |
---|
| 302 | << pSol->GetName() << |
---|
| 303 | "\") are not visualisable." |
---|
| 304 | << G4endl; |
---|
| 305 | visualisable = false; |
---|
| 306 | } |
---|
| 307 | break; |
---|
| 308 | case kPhi: |
---|
| 309 | rotation.rotateZ (-(offset+(n+0.5)*width)); |
---|
| 310 | // Minus Sign because for the physical volume we need the |
---|
| 311 | // coordinate system rotation. |
---|
| 312 | pRotation = &rotation; |
---|
| 313 | break; |
---|
| 314 | } |
---|
| 315 | pVPV -> SetTranslation (translation); |
---|
| 316 | pVPV -> SetRotation (pRotation); |
---|
| 317 | pVPV -> SetCopyNo (n); |
---|
| 318 | if (visualisable) { |
---|
| 319 | DescribeAndDescend (pVPV, requestedDepth, pLV, pSol, pMaterial, |
---|
| 320 | theAT, sceneHandler); |
---|
| 321 | } |
---|
| 322 | } |
---|
| 323 | // Restore originals... |
---|
| 324 | pVPV -> SetTranslation (originalTranslation); |
---|
| 325 | pVPV -> SetRotation (pOriginalRotation); |
---|
| 326 | if (axis == kRho && pSol->GetEntityType() == "G4Tubs") { |
---|
| 327 | ((G4Tubs*)pSol)->SetInnerRadius(originalRMin); |
---|
| 328 | ((G4Tubs*)pSol)->SetOuterRadius(originalRMax); |
---|
| 329 | } |
---|
| 330 | } |
---|
| 331 | } |
---|
| 332 | |
---|
| 333 | return; |
---|
| 334 | } |
---|
| 335 | |
---|
| 336 | void G4PhysicalVolumeModel::DescribeAndDescend |
---|
| 337 | (G4VPhysicalVolume* pVPV, |
---|
| 338 | G4int requestedDepth, |
---|
| 339 | G4LogicalVolume* pLV, |
---|
| 340 | G4VSolid* pSol, |
---|
| 341 | G4Material* pMaterial, |
---|
| 342 | const G4Transform3D& theAT, |
---|
| 343 | G4VGraphicsScene& sceneHandler) |
---|
| 344 | { |
---|
| 345 | // Maintain useful data members... |
---|
| 346 | fpCurrentPV = pVPV; |
---|
| 347 | fpCurrentLV = pLV; |
---|
| 348 | fpCurrentMaterial = pMaterial; |
---|
| 349 | |
---|
| 350 | const G4RotationMatrix objectRotation = pVPV -> GetObjectRotationValue (); |
---|
| 351 | const G4ThreeVector& translation = pVPV -> GetTranslation (); |
---|
| 352 | G4Transform3D theLT (G4Transform3D (objectRotation, translation)); |
---|
| 353 | |
---|
| 354 | // Compute the accumulated transformation... |
---|
| 355 | // Note that top volume's transformation relative to the world |
---|
| 356 | // coordinate system is specified in theAT == startingTransformation |
---|
| 357 | // = fTransform (see DescribeYourselfTo), so first time through the |
---|
| 358 | // volume's own transformation, which is only relative to its |
---|
| 359 | // mother, i.e., not relative to the world coordinate system, should |
---|
| 360 | // not be accumulated. |
---|
| 361 | G4Transform3D theNewAT (theAT); |
---|
| 362 | if (fCurrentDepth != 0) theNewAT = theAT * theLT; |
---|
| 363 | fpCurrentTransform = &theNewAT; |
---|
| 364 | |
---|
| 365 | /******************************************************** |
---|
| 366 | G4cout << "G4PhysicalVolumeModel::DescribeAndDescend: " |
---|
| 367 | << pVPV -> GetName () << "." << pVPV -> GetCopyNo (); |
---|
| 368 | G4cout << "\n theAT: "; |
---|
| 369 | G4cout << "\n Rotation: "; |
---|
| 370 | G4RotationMatrix rotation = theAT.getRotation (); |
---|
| 371 | G4cout << rotation.thetaX() << ", " |
---|
| 372 | << rotation.phiX() << ", " |
---|
| 373 | << rotation.thetaY() << ", " |
---|
| 374 | << rotation.phiY() << ", " |
---|
| 375 | << rotation.thetaZ() << ", " |
---|
| 376 | << rotation.phiZ(); |
---|
| 377 | G4cout << "\n Translation: " << theAT.getTranslation(); |
---|
| 378 | G4cout << "\n theNewAT: "; |
---|
| 379 | G4cout << "\n Rotation: "; |
---|
| 380 | rotation = theNewAT.getRotation (); |
---|
| 381 | G4cout << rotation.thetaX() << ", " |
---|
| 382 | << rotation.phiX() << ", " |
---|
| 383 | << rotation.thetaY() << ", " |
---|
| 384 | << rotation.phiY() << ", " |
---|
| 385 | << rotation.thetaZ() << ", " |
---|
| 386 | << rotation.phiZ(); |
---|
| 387 | G4cout << "\n Translation: " << theNewAT.getTranslation(); |
---|
| 388 | G4cout << G4endl; |
---|
| 389 | **********************************************************/ |
---|
| 390 | |
---|
| 391 | // Make decision to draw... |
---|
| 392 | const G4VisAttributes* pVisAttribs = pLV->GetVisAttributes(); |
---|
| 393 | if (!pVisAttribs) pVisAttribs = fpMP->GetDefaultVisAttributes(); |
---|
| 394 | // Beware - pVisAttribs might still be zero - create a temporary default one... |
---|
| 395 | G4bool visAttsCreated = false; |
---|
| 396 | if (!pVisAttribs) { |
---|
| 397 | pVisAttribs = new G4VisAttributes; |
---|
| 398 | visAttsCreated = true; |
---|
| 399 | } |
---|
| 400 | |
---|
| 401 | // From here, can assume pVisAttribs is a valid pointer. |
---|
| 402 | |
---|
| 403 | G4bool thisToBeDrawn = true; |
---|
| 404 | |
---|
| 405 | // There are various reasons why this volume |
---|
| 406 | // might not be drawn... |
---|
| 407 | G4bool culling = fpMP->IsCulling(); |
---|
| 408 | G4bool cullingInvisible = fpMP->IsCullingInvisible(); |
---|
| 409 | G4bool markedVisible = pVisAttribs->IsVisible(); |
---|
| 410 | G4bool cullingLowDensity = fpMP->IsDensityCulling(); |
---|
| 411 | G4double density = pMaterial? pMaterial->GetDensity(): 0; |
---|
| 412 | G4double densityCut = fpMP -> GetVisibleDensity (); |
---|
| 413 | |
---|
| 414 | // 1) Global culling is on.... |
---|
| 415 | if (culling) { |
---|
| 416 | // 2) Culling of invisible volumes is on... |
---|
| 417 | if (cullingInvisible) { |
---|
| 418 | // 3) ...and the volume is marked not visible... |
---|
| 419 | if (!markedVisible) thisToBeDrawn = false; |
---|
| 420 | } |
---|
| 421 | // 4) Or culling of low density volumes is on... |
---|
| 422 | if (cullingLowDensity) { |
---|
| 423 | // 5) ...and density is less than cut value... |
---|
| 424 | if (density < densityCut) thisToBeDrawn = false; |
---|
| 425 | } |
---|
| 426 | } |
---|
| 427 | |
---|
| 428 | // Update full path of physical volumes... |
---|
| 429 | G4int copyNo = fpCurrentPV->GetCopyNo(); |
---|
| 430 | fFullPVPath.push_back |
---|
| 431 | (G4PhysicalVolumeNodeID(fpCurrentPV,copyNo,fCurrentDepth)); |
---|
| 432 | |
---|
| 433 | if (thisToBeDrawn) { |
---|
| 434 | |
---|
| 435 | // Update path of drawn physical volumes... |
---|
| 436 | G4int copyNo = fpCurrentPV->GetCopyNo(); |
---|
| 437 | fDrawnPVPath.push_back |
---|
| 438 | (G4PhysicalVolumeNodeID(fpCurrentPV,copyNo,fCurrentDepth)); |
---|
| 439 | |
---|
| 440 | if (fpMP->IsExplode() && fDrawnPVPath.size() == 1) { |
---|
| 441 | // For top-level drawn volumes, explode along radius... |
---|
| 442 | G4Transform3D centering = G4Translate3D(fpMP->GetExplodeCentre()); |
---|
| 443 | G4Transform3D centred = centering.inverse() * theNewAT; |
---|
| 444 | G4Scale3D scale; |
---|
| 445 | G4Rotate3D rotation; |
---|
| 446 | G4Translate3D translation; |
---|
| 447 | centred.getDecomposition(scale, rotation, translation); |
---|
| 448 | G4double explodeFactor = fpMP->GetExplodeFactor(); |
---|
| 449 | G4Translate3D newTranslation = |
---|
| 450 | G4Translate3D(explodeFactor * translation.dx(), |
---|
| 451 | explodeFactor * translation.dy(), |
---|
| 452 | explodeFactor * translation.dz()); |
---|
| 453 | theNewAT = centering * newTranslation * rotation * scale; |
---|
| 454 | } |
---|
| 455 | |
---|
| 456 | DescribeSolid (theNewAT, pSol, pVisAttribs, sceneHandler); |
---|
| 457 | |
---|
| 458 | } |
---|
| 459 | |
---|
| 460 | // Make decision to draw daughters, if any. There are various |
---|
| 461 | // reasons why daughters might not be drawn... |
---|
| 462 | |
---|
| 463 | // First, reasons that do not depend on culling policy... |
---|
| 464 | G4int nDaughters = pLV->GetNoDaughters(); |
---|
| 465 | G4bool daughtersToBeDrawn = true; |
---|
| 466 | // 1) There are no daughters... |
---|
| 467 | if (!nDaughters) daughtersToBeDrawn = false; |
---|
| 468 | // 2) We are at the limit if requested depth... |
---|
| 469 | else if (requestedDepth == 0) daughtersToBeDrawn = false; |
---|
| 470 | // 3) The user has asked that the descent be curtailed... |
---|
| 471 | else if (fCurtailDescent) daughtersToBeDrawn = false; |
---|
| 472 | |
---|
| 473 | // Now, reasons that depend on culling policy... |
---|
| 474 | else { |
---|
| 475 | G4bool culling = fpMP->IsCulling(); |
---|
| 476 | G4bool cullingInvisible = fpMP->IsCullingInvisible(); |
---|
| 477 | G4bool daughtersInvisible = pVisAttribs->IsDaughtersInvisible(); |
---|
| 478 | // Culling of covered daughters request. This is computed in |
---|
| 479 | // G4VSceneHandler::CreateModelingParameters() depending on view |
---|
| 480 | // parameters... |
---|
| 481 | G4bool cullingCovered = fpMP->IsCullingCovered(); |
---|
| 482 | G4bool surfaceDrawing = |
---|
| 483 | fpMP->GetDrawingStyle() == G4ModelingParameters::hsr || |
---|
| 484 | fpMP->GetDrawingStyle() == G4ModelingParameters::hlhsr; |
---|
| 485 | if (pVisAttribs->IsForceDrawingStyle()) { |
---|
| 486 | switch (pVisAttribs->GetForcedDrawingStyle()) { |
---|
| 487 | default: |
---|
| 488 | case G4VisAttributes::wireframe: surfaceDrawing = false; break; |
---|
| 489 | case G4VisAttributes::solid: surfaceDrawing = true; break; |
---|
| 490 | } |
---|
| 491 | } |
---|
| 492 | G4bool opaque = pVisAttribs->GetColour().GetAlpha() >= 1.; |
---|
| 493 | // 4) Global culling is on.... |
---|
| 494 | if (culling) { |
---|
| 495 | // 5) ..and culling of invisible volumes is on... |
---|
| 496 | if (cullingInvisible) { |
---|
| 497 | // 6) ...and the mother requests daughters invisible |
---|
| 498 | if (daughtersInvisible) daughtersToBeDrawn = false; |
---|
| 499 | } |
---|
| 500 | // 7) Or culling of covered daughters is requested... |
---|
| 501 | if (cullingCovered) { |
---|
| 502 | // 8) ...and surface drawing is operating... |
---|
| 503 | if (surfaceDrawing) { |
---|
| 504 | // 9) ...but only if mother is visible... |
---|
| 505 | if (thisToBeDrawn) { |
---|
| 506 | // 10) ...and opaque... |
---|
| 507 | if (opaque) daughtersToBeDrawn = false; |
---|
| 508 | } |
---|
| 509 | } |
---|
| 510 | } |
---|
| 511 | } |
---|
| 512 | } |
---|
| 513 | |
---|
| 514 | // Vis atts for this volume no longer needed if created... |
---|
| 515 | if (visAttsCreated) delete pVisAttribs; |
---|
| 516 | |
---|
| 517 | if (daughtersToBeDrawn) { |
---|
| 518 | for (G4int iDaughter = 0; iDaughter < nDaughters; iDaughter++) { |
---|
| 519 | G4VPhysicalVolume* pVPV = pLV -> GetDaughter (iDaughter); |
---|
| 520 | // Descend the geometry structure recursively... |
---|
| 521 | fCurrentDepth++; |
---|
| 522 | VisitGeometryAndGetVisReps |
---|
| 523 | (pVPV, requestedDepth - 1, theNewAT, sceneHandler); |
---|
| 524 | fCurrentDepth--; |
---|
| 525 | } |
---|
| 526 | } |
---|
| 527 | |
---|
| 528 | // Reset for normal descending of next volume at this level... |
---|
| 529 | fCurtailDescent = false; |
---|
| 530 | |
---|
| 531 | // Pop item from paths physical volumes... |
---|
| 532 | fFullPVPath.pop_back(); |
---|
| 533 | if (thisToBeDrawn) { |
---|
| 534 | fDrawnPVPath.pop_back(); |
---|
| 535 | } |
---|
| 536 | } |
---|
| 537 | |
---|
| 538 | void G4PhysicalVolumeModel::DescribeSolid |
---|
| 539 | (const G4Transform3D& theAT, |
---|
| 540 | G4VSolid* pSol, |
---|
| 541 | const G4VisAttributes* pVisAttribs, |
---|
| 542 | G4VGraphicsScene& sceneHandler) |
---|
| 543 | { |
---|
| 544 | sceneHandler.PreAddSolid (theAT, *pVisAttribs); |
---|
| 545 | |
---|
| 546 | const G4Polyhedron* pSectionPolyhedron = fpMP->GetSectionPolyhedron(); |
---|
| 547 | const G4Polyhedron* pCutawayPolyhedron = fpMP->GetCutawayPolyhedron(); |
---|
| 548 | |
---|
| 549 | if (!fpClippingPolyhedron && !pSectionPolyhedron && !pCutawayPolyhedron) { |
---|
| 550 | |
---|
| 551 | pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment. |
---|
| 552 | |
---|
| 553 | } else { |
---|
| 554 | |
---|
| 555 | // Clipping, etc., performed by Boolean operations on polyhedron objects. |
---|
| 556 | |
---|
| 557 | // First, get polyhedron for current solid... |
---|
| 558 | if (pVisAttribs->IsForceLineSegmentsPerCircle()) |
---|
| 559 | G4Polyhedron::SetNumberOfRotationSteps |
---|
| 560 | (pVisAttribs->GetForcedLineSegmentsPerCircle()); |
---|
| 561 | else |
---|
| 562 | G4Polyhedron::SetNumberOfRotationSteps(fpMP->GetNoOfSides()); |
---|
| 563 | G4Polyhedron* pOriginal = pSol->GetPolyhedron(); |
---|
| 564 | G4Polyhedron::ResetNumberOfRotationSteps(); |
---|
| 565 | if (!pOriginal) { |
---|
| 566 | if (fpMP->IsWarning()) |
---|
| 567 | G4cout << |
---|
| 568 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: solid\n \"" |
---|
| 569 | << pSol->GetName() << |
---|
| 570 | "\" has no polyhedron. Cannot by clipped." |
---|
| 571 | << G4endl; |
---|
| 572 | pSol -> DescribeYourselfTo (sceneHandler); // Standard treatment. |
---|
| 573 | } else { |
---|
| 574 | |
---|
| 575 | G4Polyhedron resultant = *pOriginal; |
---|
| 576 | |
---|
| 577 | if (fpClippingPolyhedron) { |
---|
| 578 | G4Polyhedron clipper = *fpClippingPolyhedron; // Local copy. |
---|
| 579 | clipper.Transform(theAT.inverse()); |
---|
| 580 | switch (fClippingMode) { |
---|
| 581 | default: |
---|
| 582 | case subtraction: resultant = resultant.subtract(clipper); break; |
---|
| 583 | case intersection: resultant = resultant.intersect(clipper); break; |
---|
| 584 | } |
---|
| 585 | if(resultant.IsErrorBooleanProcess()) { |
---|
| 586 | if (fpMP->IsWarning()) |
---|
| 587 | G4cout << |
---|
| 588 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: clipped polyhedron for" |
---|
| 589 | "\n solid \"" << pSol->GetName() << |
---|
| 590 | "\" not defined due to error during Boolean processing." |
---|
| 591 | << G4endl; |
---|
| 592 | // Nevertheless, keep resultant. |
---|
| 593 | } |
---|
| 594 | } |
---|
| 595 | |
---|
| 596 | if (pSectionPolyhedron) { |
---|
| 597 | G4Polyhedron sectioner = *pSectionPolyhedron; // Local copy. |
---|
| 598 | sectioner.Transform(theAT.inverse()); |
---|
| 599 | resultant = resultant.intersect(sectioner); |
---|
| 600 | if(resultant.IsErrorBooleanProcess()) { |
---|
| 601 | if (fpMP->IsWarning()) |
---|
| 602 | G4cout << |
---|
| 603 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: sectioned polyhedron for" |
---|
| 604 | "\n solid \"" << pSol->GetName() << |
---|
| 605 | "\" not defined due to error during Boolean processing." |
---|
| 606 | << G4endl; |
---|
| 607 | // Nevertheless, keep resultant. |
---|
| 608 | } |
---|
| 609 | } |
---|
| 610 | |
---|
| 611 | if (pCutawayPolyhedron) { |
---|
| 612 | G4Polyhedron cutter = *pCutawayPolyhedron; // Local copy. |
---|
| 613 | cutter.Transform(theAT.inverse()); |
---|
| 614 | resultant = resultant.subtract(cutter); |
---|
| 615 | if(resultant.IsErrorBooleanProcess()) { |
---|
| 616 | if (fpMP->IsWarning()) |
---|
| 617 | G4cout << |
---|
| 618 | "WARNING: G4PhysicalVolumeModel::DescribeSolid: cutaway polyhedron for" |
---|
| 619 | "\n solid \"" << pSol->GetName() << |
---|
| 620 | "\" not defined due to error during Boolean processing." |
---|
| 621 | << G4endl; |
---|
| 622 | // Nevertheless, keep resultant. |
---|
| 623 | } |
---|
| 624 | } |
---|
| 625 | |
---|
| 626 | // Finally, force polyhedron drawing... |
---|
| 627 | resultant.SetVisAttributes(pVisAttribs); |
---|
| 628 | sceneHandler.BeginPrimitives(theAT); |
---|
| 629 | sceneHandler.AddPrimitive(resultant); |
---|
| 630 | sceneHandler.EndPrimitives(); |
---|
| 631 | } |
---|
| 632 | } |
---|
| 633 | sceneHandler.PostAddSolid (); |
---|
| 634 | } |
---|
| 635 | |
---|
| 636 | G4bool G4PhysicalVolumeModel::Validate (G4bool warn) |
---|
| 637 | { |
---|
| 638 | G4VPhysicalVolume* world = |
---|
| 639 | G4TransportationManager::GetTransportationManager () |
---|
| 640 | -> GetNavigatorForTracking () -> GetWorldVolume (); |
---|
| 641 | // The idea now is to seek a PV with the same name and copy no |
---|
| 642 | // in the hope it's the same one!! |
---|
| 643 | if (warn) { |
---|
| 644 | G4cout << "G4PhysicalVolumeModel::Validate() called." << G4endl; |
---|
| 645 | } |
---|
| 646 | G4PhysicalVolumeModel searchModel (world); |
---|
| 647 | G4PhysicalVolumeSearchScene searchScene |
---|
| 648 | (&searchModel, fTopPVName, fTopPVCopyNo); |
---|
| 649 | G4ModelingParameters mp; // Default modeling parameters for this search. |
---|
| 650 | mp.SetDefaultVisAttributes(fpMP? fpMP->GetDefaultVisAttributes(): 0); |
---|
| 651 | searchModel.SetModelingParameters (&mp); |
---|
| 652 | searchModel.DescribeYourselfTo (searchScene); |
---|
| 653 | G4VPhysicalVolume* foundVolume = searchScene.GetFoundVolume (); |
---|
| 654 | if (foundVolume) { |
---|
| 655 | if (warn) { |
---|
| 656 | G4cout << " Volume of the same name and copy number (\"" |
---|
| 657 | << fTopPVName << "\", copy " << fTopPVCopyNo |
---|
| 658 | << ") still exists and is being used." |
---|
| 659 | "\n WARNING: This does not necessarily guarantee it's the same" |
---|
| 660 | "\n volume you originally specified in /vis/scene/add/." |
---|
| 661 | << G4endl; |
---|
| 662 | } |
---|
| 663 | fpTopPV = foundVolume; |
---|
| 664 | CalculateExtent (); |
---|
| 665 | return true; |
---|
| 666 | } |
---|
| 667 | else { |
---|
| 668 | if (warn) { |
---|
| 669 | G4cout << " A volume of the same name and copy number (\"" |
---|
| 670 | << fTopPVName << "\", copy " << fTopPVCopyNo |
---|
| 671 | << ") no longer exists." |
---|
| 672 | << G4endl; |
---|
| 673 | } |
---|
| 674 | return false; |
---|
| 675 | } |
---|
| 676 | } |
---|
| 677 | |
---|
| 678 | const std::map<G4String,G4AttDef>* G4PhysicalVolumeModel::GetAttDefs() const |
---|
| 679 | { |
---|
| 680 | G4bool isNew; |
---|
| 681 | std::map<G4String,G4AttDef>* store |
---|
| 682 | = G4AttDefStore::GetInstance("G4PhysicalVolumeModel", isNew); |
---|
| 683 | if (isNew) { |
---|
| 684 | (*store)["PVPath"] = |
---|
| 685 | G4AttDef("PVPath","Physical Volume Path","Physics","","G4String"); |
---|
| 686 | (*store)["LVol"] = |
---|
| 687 | G4AttDef("LVol","Logical Volume","Physics","","G4String"); |
---|
| 688 | (*store)["Solid"] = |
---|
| 689 | G4AttDef("Solid","Solid Name","Physics","","G4String"); |
---|
| 690 | (*store)["EType"] = |
---|
| 691 | G4AttDef("EType","Entity Type","Physics","","G4String"); |
---|
| 692 | (*store)["DmpSol"] = |
---|
| 693 | G4AttDef("DmpSol","Dump of Solid properties","Physics","","G4String"); |
---|
| 694 | (*store)["Trans"] = |
---|
| 695 | G4AttDef("Trans","Transformation of volume","Physics","","G4String"); |
---|
| 696 | (*store)["Material"] = |
---|
| 697 | G4AttDef("Material","Material Name","Physics","","G4String"); |
---|
| 698 | (*store)["Density"] = |
---|
| 699 | G4AttDef("Density","Material Density","Physics","G4BestUnit","G4double"); |
---|
| 700 | (*store)["State"] = |
---|
| 701 | G4AttDef("State","Material State (enum undefined,solid,liquid,gas)","Physics","","G4String"); |
---|
| 702 | (*store)["Radlen"] = |
---|
| 703 | G4AttDef("Radlen","Material Radiation Length","Physics","G4BestUnit","G4double"); |
---|
| 704 | } |
---|
| 705 | (*store)["Region"] = |
---|
| 706 | G4AttDef("Region","Cuts Region","Physics","","G4String"); |
---|
| 707 | (*store)["RootRegion"] = |
---|
| 708 | G4AttDef("RootRegion","Root Region (0/1 = false/true)","Physics","","G4bool"); |
---|
| 709 | return store; |
---|
| 710 | } |
---|
| 711 | |
---|
| 712 | #include <iomanip> |
---|
| 713 | |
---|
| 714 | static std::ostream& operator<< (std::ostream& o, const G4Transform3D t) |
---|
| 715 | { |
---|
| 716 | using namespace std; |
---|
| 717 | |
---|
| 718 | G4Scale3D s; |
---|
| 719 | G4Rotate3D r; |
---|
| 720 | G4Translate3D tl; |
---|
| 721 | t.getDecomposition(s, r, tl); |
---|
| 722 | |
---|
| 723 | const int w = 10; |
---|
| 724 | |
---|
| 725 | // Transformation itself |
---|
| 726 | o << setw(w) << t.xx() << setw(w) << t.xy() << setw(w) << t.xz() << setw(w) << t.dx() << endl; |
---|
| 727 | o << setw(w) << t.yx() << setw(w) << t.yy() << setw(w) << t.yz() << setw(w) << t.dy() << endl; |
---|
| 728 | o << setw(w) << t.zx() << setw(w) << t.zy() << setw(w) << t.zz() << setw(w) << t.dz() << endl; |
---|
| 729 | |
---|
| 730 | // Translation |
---|
| 731 | o << "= translation:" << endl; |
---|
| 732 | o << setw(w) << tl.dx() << setw(w) << tl.dy() << setw(w) << tl.dz() << endl; |
---|
| 733 | |
---|
| 734 | // Rotation |
---|
| 735 | o << "* rotation:" << endl; |
---|
| 736 | o << setw(w) << r.xx() << setw(w) << r.xy() << setw(w) << r.xz() << endl; |
---|
| 737 | o << setw(w) << r.yx() << setw(w) << r.yy() << setw(w) << r.yz() << endl; |
---|
| 738 | o << setw(w) << r.zx() << setw(w) << r.zy() << setw(w) << r.zz() << endl; |
---|
| 739 | |
---|
| 740 | // Scale |
---|
| 741 | o << "* scale:" << endl; |
---|
| 742 | o << setw(w) << s.xx() << setw(w) << s.yy() << setw(w) << s.zz() << endl; |
---|
| 743 | |
---|
| 744 | // Transformed axes |
---|
| 745 | o << "Transformed axes:" << endl; |
---|
| 746 | o << "x': " << r * G4Vector3D(1., 0., 0.) << endl; |
---|
| 747 | o << "y': " << r * G4Vector3D(0., 1., 0.) << endl; |
---|
| 748 | o << "z': " << r * G4Vector3D(0., 0., 1.) << endl; |
---|
| 749 | |
---|
| 750 | return o; |
---|
| 751 | } |
---|
| 752 | |
---|
| 753 | std::vector<G4AttValue>* G4PhysicalVolumeModel::CreateCurrentAttValues() const |
---|
| 754 | { |
---|
| 755 | std::vector<G4AttValue>* values = new std::vector<G4AttValue>; |
---|
| 756 | std::ostringstream oss; |
---|
| 757 | for (size_t i = 0; i < fFullPVPath.size(); ++i) { |
---|
| 758 | oss << fFullPVPath[i].GetPhysicalVolume()->GetName() |
---|
| 759 | << ':' << fFullPVPath[i].GetCopyNo(); |
---|
| 760 | if (i != fFullPVPath.size() - 1) oss << '/'; |
---|
| 761 | } |
---|
| 762 | values->push_back(G4AttValue("PVPath", oss.str(),"")); |
---|
| 763 | values->push_back(G4AttValue("LVol", fpCurrentLV->GetName(),"")); |
---|
| 764 | G4VSolid* pSol = fpCurrentLV->GetSolid(); |
---|
| 765 | values->push_back(G4AttValue("Solid", pSol->GetName(),"")); |
---|
| 766 | values->push_back(G4AttValue("EType", pSol->GetEntityType(),"")); |
---|
| 767 | oss.str(""); oss << '\n' << *pSol; |
---|
| 768 | values->push_back(G4AttValue("DmpSol", oss.str(),"")); |
---|
| 769 | oss.str(""); oss << '\n' << *fpCurrentTransform; |
---|
| 770 | values->push_back(G4AttValue("Trans", oss.str(),"")); |
---|
| 771 | G4String matName = fpCurrentMaterial? fpCurrentMaterial->GetName(): G4String("No material"); |
---|
| 772 | values->push_back(G4AttValue("Material", matName,"")); |
---|
| 773 | G4double matDensity = fpCurrentMaterial? fpCurrentMaterial->GetDensity(): 0.; |
---|
| 774 | values->push_back(G4AttValue("Density", G4BestUnit(matDensity,"Volumic Mass"),"")); |
---|
| 775 | G4State matState = fpCurrentMaterial? fpCurrentMaterial->GetState(): kStateUndefined; |
---|
| 776 | oss.str(""); oss << matState; |
---|
| 777 | values->push_back(G4AttValue("State", oss.str(),"")); |
---|
| 778 | G4double matRadlen = fpCurrentMaterial? fpCurrentMaterial->GetRadlen(): 0.; |
---|
| 779 | values->push_back(G4AttValue("Radlen", G4BestUnit(matRadlen,"Length"),"")); |
---|
| 780 | G4Region* region = fpCurrentLV->GetRegion(); |
---|
| 781 | G4String regionName = region? region->GetName(): G4String("No region"); |
---|
| 782 | values->push_back(G4AttValue("Region", regionName,"")); |
---|
| 783 | oss.str(""); oss << fpCurrentLV->IsRootRegion(); |
---|
| 784 | values->push_back(G4AttValue("RootRegion", oss.str(),"")); |
---|
| 785 | return values; |
---|
| 786 | } |
---|