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 // G4PSPassageTrackLength 29 #include "G4PSPassageTrackLength.hh" 30 #include "G4StepStatus.hh" 31 #include "G4Track.hh" 32 #include "G4UnitsTable.hh" 33 #include "G4VScoreHistFiller.hh" 34 35 //////////////////////////////////////////////////////////////////////////////// 36 // (Description) 37 // This is a primitive scorer class for scoring only track length. 38 // The tracks which passed a geometry is taken into account. 39 // 40 // 41 // Created: 2005-11-14 Tsukasa ASO, Akinori Kimura. 42 // 2010-07-22 Introduce Unit specification. 43 // 2020-10-06 Use G4VPrimitivePlotter and fill 1-D histo of track length 44 // vs. number of tracks (not weighted) (Makoto Asai) 45 // 46 /////////////////////////////////////////////////////////////////////////////// 47 48 G4PSPassageTrackLength::G4PSPassageTrackLength(const G4String& name, G4int depth) 49 : G4PSPassageTrackLength(name, "mm", depth) 50 {} 51 52 G4PSPassageTrackLength::G4PSPassageTrackLength(const G4String& name, 53 const G4String& unit, 54 G4int depth) 55 : G4VPrimitivePlotter(name, depth) 56 , HCID(-1) 57 , fCurrentTrkID(-1) 58 , fTrackLength(0.) 59 , EvtMap(nullptr) 60 , weighted(false) 61 { 62 SetUnit(unit); 63 } 64 65 G4bool G4PSPassageTrackLength::ProcessHits(G4Step* aStep, G4TouchableHistory*) 66 { 67 if(IsPassed(aStep)) 68 { 69 G4int index = GetIndex(aStep); 70 EvtMap->add(index, fTrackLength); 71 72 if(!hitIDMap.empty() && hitIDMap.find(index) != hitIDMap.cend()) 73 { 74 auto filler = G4VScoreHistFiller::Instance(); 75 if(filler == nullptr) 76 { 77 G4Exception( 78 "G4PSPassageTrackLength::ProcessHits", "SCORER0123", JustWarning, 79 "G4TScoreHistFiller is not instantiated!! Histogram is not filled."); 80 } 81 else 82 { 83 filler->FillH1(hitIDMap[index], fTrackLength, 1.); 84 } 85 } 86 } 87 88 return true; 89 } 90 91 G4bool G4PSPassageTrackLength::IsPassed(G4Step* aStep) 92 { 93 G4bool Passed = false; 94 95 G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary; 96 G4bool IsExit = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary; 97 98 G4int trkid = aStep->GetTrack()->GetTrackID(); 99 G4double trklength = aStep->GetStepLength(); 100 if(weighted) 101 trklength *= aStep->GetPreStepPoint()->GetWeight(); 102 103 if(IsEnter && IsExit) 104 { // Passed at one step 105 fTrackLength = trklength; // Track length is absolutely given. 106 Passed = true; 107 } 108 else if(IsEnter) 109 { // Enter a new geometry 110 fCurrentTrkID = trkid; // Resetting the current track. 111 fTrackLength = trklength; 112 } 113 else if(IsExit) 114 { // Exit a current geometry 115 if(fCurrentTrkID == trkid) 116 { 117 fTrackLength += trklength; // Adding the track length to current one, 118 Passed = true; // if the track is same as entered. 119 } 120 } 121 else 122 { // Inside geometry 123 if(fCurrentTrkID == trkid) 124 { // Adding the track length to current one , 125 fTrackLength += trklength; // if the track is same as entered. 126 } 127 } 128 129 return Passed; 130 } 131 132 void G4PSPassageTrackLength::Initialize(G4HCofThisEvent* HCE) 133 { 134 fCurrentTrkID = -1; 135 136 EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName()); 137 if(HCID < 0) 138 HCID = GetCollectionID(0); 139 HCE->AddHitsCollection(HCID, EvtMap); 140 } 141 142 void G4PSPassageTrackLength::clear() { EvtMap->clear(); } 143 144 void G4PSPassageTrackLength::PrintAll() 145 { 146 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl; 147 G4cout << " PrimitiveSenstivity " << GetName() << G4endl; 148 G4cout << " Number of entries " << EvtMap->entries() << G4endl; 149 for(const auto& [copy, length] : *(EvtMap->GetMap())) 150 { 151 G4cout << " copy no.: " << copy 152 << " track length : " << *(length) / GetUnitValue() << " [" 153 << GetUnit() << "]" << G4endl; 154 } 155 } 156 157 void G4PSPassageTrackLength::SetUnit(const G4String& unit) 158 { 159 CheckAndSetUnit(unit, "Length"); 160 } 161