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 // This example is provided by the Geant4-DNA collaboration 27 // Any report or published results obtained using the Geant4-DNA software 28 // shall cite the following Geant4-DNA collaboration publications: 29 // Med. Phys. 45 (2018) e722-e739 30 // Phys. Med. 31 (2015) 861-874 31 // Med. Phys. 37 (2010) 4692-4708 32 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 33 // 34 // The Geant4-DNA web site is available at http://geant4-dna.org 35 // 36 /// \file TrackerSD.cc 37 /// \brief Implementation of the TrackerSD class 38 39 #include "TrackerSD.hh" 40 41 #include "G4AnalysisManager.hh" 42 #include "G4SDManager.hh" 43 #include "G4SystemOfUnits.hh" 44 #include "Randomize.hh" 45 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 47 48 TrackerSD::TrackerSD(const G4String& name, const G4String& hitsCollectionName) 49 : G4VSensitiveDetector(name), fHitsCollection(NULL) 50 { 51 collectionName.insert(hitsCollectionName); 52 } 53 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 55 56 TrackerSD::~TrackerSD() {} 57 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 59 60 void TrackerSD::Initialize(G4HCofThisEvent* hce) 61 { 62 // Create hits collection 63 fHitsCollection = new TrackerHitsCollection(SensitiveDetectorName, collectionName[0]); 64 65 // Add this collection in hce 66 G4int hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); 67 68 hce->AddHitsCollection(hcID, fHitsCollection); 69 } 70 71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 72 73 G4bool TrackerSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) 74 { 75 // Energy deposit 76 G4double edep = aStep->GetTotalEnergyDeposit(); 77 78 if (edep == 0.) return false; 79 80 TrackerHit* newHit = new TrackerHit(); 81 82 newHit->SetTrackID(aStep->GetTrack()->GetTrackID()); 83 newHit->SetEdep(edep); 84 newHit->SetPos(aStep->GetPostStepPoint()->GetPosition()); 85 86 if (aStep->GetTrack()->GetTrackID() == 1 && aStep->GetTrack()->GetParentID() == 0) 87 newHit->SetIncidentEnergy(aStep->GetTrack()->GetVertexKineticEnergy()); 88 89 fHitsCollection->insert(newHit); 90 91 // Print 92 // newHit->Print(); 93 94 return true; 95 } 96 97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 98 void TrackerSD::EndOfEvent(G4HCofThisEvent*) 99 { 100 G4int nofHits = fHitsCollection->entries(); 101 102 G4double Einc = 0; 103 104 /* 105 G4cout << G4endl 106 << "-------->Hits Collection: in this event they are " 107 << nofHits 108 << " hits in the target volume " << G4endl; 109 */ 110 111 // PROCESSING OF PROXIMITY FUNCTION t(x) 112 113 // ************************************* 114 // Please select herebelow : 115 // - the minimum value of x radius 116 // - the maximum value of x radius 117 // - the number of steps in radius 118 119 G4double minRadius = 0.1 * nm; 120 121 G4double maxRadius = 10000 * nm; 122 G4int nRadiusSteps = 101; 123 124 // 125 126 auto analysisManager = G4AnalysisManager::Instance(); 127 128 G4double radius(minRadius); 129 G4double stpRadius(std::pow(maxRadius / radius, 1. / static_cast<G4double>(nRadiusSteps - 1))); 130 G4int step(nRadiusSteps); 131 G4int noRadius(0); 132 133 // 1) loop on radius 134 135 while (step > 0) { 136 step--; 137 noRadius = nRadiusSteps - step; 138 139 // G4cout << "---radius/nm=" << radius/nm << G4endl; 140 141 // Computation of t(x) 142 143 G4double tNum = 0.; 144 G4double tDenom = 0.; 145 G4int nbEdep = 0; 146 147 // 2) loop on hits 148 149 for (G4int k = 0; k < nofHits; k++) { 150 G4ThreeVector hitPos = (*fHitsCollection)[k]->GetPos(); 151 G4double hitNrj = (*fHitsCollection)[k]->GetEdep(); 152 153 // G4cout << "======> hit position x/nm =" << hitPos.x()/nm << G4endl; 154 // G4cout << "======> hit position y/nm =" << hitPos.y()/nm << G4endl; 155 // G4cout << "======> hit position z/nm =" << hitPos.z()/nm << G4endl; 156 157 // 3) loop on all other hits located within shell 158 159 G4double localSum = 0.; 160 161 for (G4int i = 0; i < nofHits; i++) { 162 if ((*fHitsCollection)[i]->GetIncidentEnergy() > 0) 163 Einc = (*fHitsCollection)[i]->GetIncidentEnergy(); 164 165 G4ThreeVector localPosi = (*fHitsCollection)[i]->GetPos(); 166 167 if (((localPosi.x() - hitPos.x()) * (localPosi.x() - hitPos.x()) 168 + (localPosi.y() - hitPos.y()) * (localPosi.y() - hitPos.y()) 169 + (localPosi.z() - hitPos.z()) * (localPosi.z() - hitPos.z()) 170 < radius * stpRadius * radius * stpRadius) 171 && ((localPosi.x() - hitPos.x()) * (localPosi.x() - hitPos.x()) 172 + (localPosi.y() - hitPos.y()) * (localPosi.y() - hitPos.y()) 173 + (localPosi.z() - hitPos.z()) * (localPosi.z() - hitPos.z()) 174 >= radius * radius)) 175 176 { 177 localSum = localSum + (*fHitsCollection)[i]->GetEdep(); 178 nbEdep = nbEdep + 1; 179 } 180 } 181 182 tNum = tNum + localSum * hitNrj; 183 tDenom = tDenom + hitNrj; 184 185 } // loop on hits 186 187 // fill ntuple including weighting 188 // does not work with ntuple merging... 189 190 analysisManager->FillNtupleDColumn(0, 0, radius / nm); 191 analysisManager->FillNtupleIColumn(0, 1, noRadius); 192 analysisManager->FillNtupleDColumn(0, 2, nofHits); 193 analysisManager->FillNtupleDColumn(0, 3, nbEdep); 194 analysisManager->FillNtupleDColumn(0, 4, (tNum / tDenom) / eV); 195 analysisManager->FillNtupleDColumn(0, 5, (stpRadius * radius) / nm); 196 analysisManager->FillNtupleDColumn(0, 6, Einc / eV); 197 analysisManager->AddNtupleRow(); 198 199 // G4cout << "---radius/nm=" << radius/nm << G4endl; 200 // G4cout << "----end of radius--- " << radius/nm << G4endl; 201 202 radius *= stpRadius; 203 204 } // loop on radii 205 } 206 207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 208