Geant4 Cross Reference |
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 // 28 // 29 // John Allison 10th August 1998. 30 // An artificial scene to find physical volumes. 31 32 #include "G4PhysicalVolumeMassScene.hh" 33 34 #include "G4VSolid.hh" 35 #include "G4Vector3D.hh" 36 #include "G4PhysicalVolumeModel.hh" 37 #include "G4LogicalVolume.hh" 38 #include "G4Polyhedron.hh" 39 #include "G4Material.hh" 40 #include "G4VPVParameterisation.hh" 41 #include "G4UnitsTable.hh" 42 43 #define G4warn G4cout 44 45 G4PhysicalVolumeMassScene::G4PhysicalVolumeMassScene 46 (G4PhysicalVolumeModel* pPVModel): 47 fpPVModel (pPVModel), 48 fVolume (0.), 49 fMass (0.), 50 fpLastPV (0), 51 fPVPCount (0), 52 fLastDepth (0), 53 fLastDensity (0.) 54 {} 55 56 G4PhysicalVolumeMassScene::~G4PhysicalVolumeMassScene () {} 57 58 void G4PhysicalVolumeMassScene::Reset () 59 { 60 fVolume = 0.; 61 fMass = 0.; 62 fpLastPV = 0; 63 fPVPCount = 0; 64 fLastDepth = 0; 65 fLastDensity = 0.; 66 fDensityStack.clear(); 67 } 68 69 void G4PhysicalVolumeMassScene::ProcessVolume (const G4VSolid& solid) 70 { 71 G4int currentDepth = fpPVModel->GetCurrentDepth(); 72 G4VPhysicalVolume* pCurrentPV = fpPVModel->GetCurrentPV(); 73 //G4LogicalVolume* pCurrentLV = fpPVModel->GetCurrentLV(); 74 G4Material* pCurrentMaterial = fpPVModel->GetCurrentMaterial(); 75 76 if (pCurrentPV != fpLastPV) { 77 fpLastPV = pCurrentPV; 78 fPVPCount = 0; 79 } 80 81 G4double currentVolume = ((G4VSolid&)solid).GetCubicVolume(); 82 G4double currentDensity = pCurrentMaterial? pCurrentMaterial->GetDensity() : 0.; 83 /* Using G4Polyhedron... (gives slightly different answers on Tubs, e.g.). 84 G4Polyhedron* pPolyhedron = solid.GetPolyhedron(); 85 if (!pPolyhedron) { 86 G4cout << 87 "G4PhysicalVolumeMassScene::AccrueMass: WARNING:" 88 "\n No G4Polyhedron for" << solid.GetEntityType() << 89 ". \"" << solid.GetName() << "\" will not be accounted." 90 "\n It will be as though not there, i.e., the density as its mother." 91 "\n Its daughters will still be found and accounted." 92 << G4endl; 93 currentVolume = 0.; 94 currentDensity = 0.; 95 } 96 */ 97 98 if (currentDepth == 0) fVolume = currentVolume; 99 100 if (currentDepth > fLastDepth) { 101 fDensityStack.push_back (fLastDensity); 102 } else if (currentDepth < fLastDepth) { 103 fDensityStack.pop_back(); 104 } 105 fLastDepth = currentDepth; 106 fLastDensity = currentDensity; 107 G4double motherDensity = 0.; 108 if (currentDepth > 0) motherDensity = fDensityStack.back(); 109 110 G4double subtractedMass = currentVolume * motherDensity; 111 G4double addedMass = currentVolume * currentDensity; 112 fMass -= subtractedMass; 113 fMass += addedMass; 114 /* Debug 115 G4cout << "current vol = " 116 << G4BestUnit (currentVolume,"Volume") 117 << ", current density = " 118 << G4BestUnit (currentDensity, "Volumic Mass") 119 << ", mother density = " 120 << G4BestUnit (motherDensity, "Volumic Mass") 121 << G4endl; 122 G4cout << "Subtracted mass = " << G4BestUnit (subtractedMass, "Mass") 123 << ", added mass = " << G4BestUnit (addedMass, "Mass") 124 << ", new mass = " << G4BestUnit (fMass, "Mass") 125 << G4endl; 126 */ 127 if (fMass < 0.) { 128 G4warn << 129 "G4PhysicalVolumeMassScene::AccrueMass: WARNING:" 130 "\n Mass going negative for \"" 131 << pCurrentPV->GetName() << 132 "\", copy " 133 << pCurrentPV->GetCopyNo() << 134 ". Larger than mother?" 135 << G4endl; 136 } 137 } 138