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 // G4PSPassageCellFlux 29 #include "G4PSPassageCellFlux.hh" 30 31 #include "G4SystemOfUnits.hh" 32 #include "G4StepStatus.hh" 33 #include "G4Track.hh" 34 #include "G4VSolid.hh" 35 #include "G4VPhysicalVolume.hh" 36 #include "G4VPVParameterisation.hh" 37 #include "G4UnitsTable.hh" 38 #include "G4VScoreHistFiller.hh" 39 40 //////////////////////////////////////////////////////////////////////////////// 41 // (Description) 42 // This is a primitive scorer class for scoring cell flux. 43 // The Cell Flux is defined by a track length divided by a geometry 44 // volume, where only tracks passing through the geometry are taken 45 // into account. e.g. the unit of Cell Flux is mm/mm3. 46 // 47 // If you want to score all tracks in the geometry volume, 48 // please use G4PSCellFlux. 49 // 50 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura. 51 // 2010-07-22 Introduce Unit specification. 52 // 2010-07-22 Add weighted option 53 // 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of kinetic energy (x) 54 // vs. cell flux * track weight (y) (Makoto Asai) 55 // 56 /////////////////////////////////////////////////////////////////////////////// 57 58 G4PSPassageCellFlux::G4PSPassageCellFlux(const G4String& name, G4int depth) 59 : G4PSPassageCellFlux(name, "percm2", depth) 60 {} 61 62 G4PSPassageCellFlux::G4PSPassageCellFlux(const G4String& name, const G4String& unit, 63 G4int depth) 64 : G4VPrimitivePlotter(name, depth) 65 , HCID(-1) 66 , fCurrentTrkID(-1) 67 , fCellFlux(0) 68 , EvtMap(nullptr) 69 , weighted(true) 70 { 71 DefineUnitAndCategory(); 72 SetUnit(unit); 73 } 74 75 G4bool G4PSPassageCellFlux::ProcessHits(G4Step* aStep, G4TouchableHistory*) 76 { 77 if(IsPassed(aStep)) 78 { 79 G4int idx = 80 ((G4TouchableHistory*) (aStep->GetPreStepPoint()->GetTouchable())) 81 ->GetReplicaNumber(indexDepth); 82 G4double cubicVolume = ComputeVolume(aStep, idx); 83 84 fCellFlux /= cubicVolume; 85 G4int index = GetIndex(aStep); 86 EvtMap->add(index, fCellFlux); 87 88 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend()) 89 { 90 auto filler = G4VScoreHistFiller::Instance(); 91 if(filler == nullptr) 92 { 93 G4Exception( 94 "G4PSPassageCellFlux::ProcessHits", "SCORER0123", JustWarning, 95 "G4TScoreHistFiller is not instantiated!! Histogram is not filled."); 96 } 97 else 98 { 99 filler->FillH1(hitIDMap[index], 100 aStep->GetPreStepPoint()->GetKineticEnergy(), fCellFlux); 101 } 102 } 103 } 104 105 return true; 106 } 107 108 G4bool G4PSPassageCellFlux::IsPassed(G4Step* aStep) 109 { 110 G4bool Passed = false; 111 112 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary; 113 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary; 114 115 G4int trkid = aStep->GetTrack()->GetTrackID(); 116 G4double trklength = aStep->GetStepLength(); 117 if(weighted) 118 trklength *= aStep->GetPreStepPoint()->GetWeight(); 119 120 if(IsEnter && IsExit) 121 { // Passed at one step 122 fCellFlux = trklength; // Track length is absolutely given. 123 Passed = true; 124 } 125 else if(IsEnter) 126 { // Enter a new geometry 127 fCurrentTrkID = trkid; // Resetting the current track. 128 fCellFlux = trklength; 129 } 130 else if(IsExit) 131 { // Exit a current geometry 132 if(fCurrentTrkID == trkid) 133 { // if the track is same as entered, 134 fCellFlux += trklength; // add the track length to current one. 135 Passed = true; 136 } 137 } 138 else 139 { // Inside geometry 140 if(fCurrentTrkID == trkid) 141 { // if the track is same as entered, 142 fCellFlux += trklength; // adding the track length to current one. 143 } 144 } 145 146 return Passed; 147 } 148 149 void G4PSPassageCellFlux::Initialize(G4HCofThisEvent* HCE) 150 { 151 fCurrentTrkID = -1; 152 153 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName()); 154 if(HCID < 0) 155 HCID = GetCollectionID(0); 156 HCE->AddHitsCollection(HCID, EvtMap); 157 } 158 159 void G4PSPassageCellFlux::clear() { EvtMap->clear(); } 160 161 void G4PSPassageCellFlux::PrintAll() 162 { 163 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl; 164 G4cout << " PrimitiveScorer " << GetName() << G4endl; 165 G4cout << " Number of entries " << EvtMap->entries() << G4endl; 166 for(const auto& [copy, flux] : *(EvtMap->GetMap())) 167 { 168 G4cout << " copy no.: " << copy 169 << " cell flux : " << *(flux) / GetUnitValue() << " [" 170 << GetUnit() << G4endl; 171 } 172 } 173 174 void G4PSPassageCellFlux::SetUnit(const G4String& unit) 175 { 176 CheckAndSetUnit(unit, "Per Unit Surface"); 177 } 178 179 void G4PSPassageCellFlux::DefineUnitAndCategory() 180 { 181 // Per Unit Surface 182 new G4UnitDefinition("percentimeter2", "percm2", "Per Unit Surface", 183 (1. / cm2)); 184 new G4UnitDefinition("permillimeter2", "permm2", "Per Unit Surface", 185 (1. / mm2)); 186 new G4UnitDefinition("permeter2", "perm2", "Per Unit Surface", (1. / m2)); 187 } 188 189 G4double G4PSPassageCellFlux::ComputeVolume(G4Step* aStep, G4int idx) 190 { 191 return ComputeSolid(aStep, idx)->GetCubicVolume(); 192 } 193