![]() |
|
|||
File indexing completed on 2025-02-23 09:21:57
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // This example is provided by the Geant4-DNA collaboration 0027 // Any report or published results obtained using the Geant4-DNA software 0028 // shall cite the following Geant4-DNA collaboration publications: 0029 // Med. Phys. 45 (2018) e722-e739 0030 // Phys. Med. 31 (2015) 861-874 0031 // Med. Phys. 37 (2010) 4692-4708 0032 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 0033 // The Geant4-DNA web site is available at http://geant4-dna.org 0034 // 0035 /// \file TrackerSD.cc 0036 /// \brief Implementation of the TrackerSD class 0037 0038 #include "TrackerSD.hh" 0039 0040 #include "G4AnalysisManager.hh" 0041 #include "G4SDManager.hh" 0042 #include "G4SystemOfUnits.hh" 0043 #include "Randomize.hh" 0044 0045 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0046 0047 TrackerSD::TrackerSD(const G4String& name, const G4String& hitsCollectionName) 0048 : G4VSensitiveDetector(name), fHitsCollection(nullptr) 0049 { 0050 collectionName.insert(hitsCollectionName); 0051 } 0052 0053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0054 0055 TrackerSD::~TrackerSD() = default; 0056 0057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0058 0059 void TrackerSD::Initialize(G4HCofThisEvent* hce) 0060 { 0061 // Create hits collection 0062 fHitsCollection = new TrackerHitsCollection(SensitiveDetectorName, collectionName[0]); 0063 0064 // Add this collection in hce 0065 G4int hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); 0066 0067 hce->AddHitsCollection(hcID, fHitsCollection); 0068 } 0069 0070 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0071 0072 G4bool TrackerSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) 0073 { 0074 // Energy deposit 0075 G4double edep = aStep->GetTotalEnergyDeposit(); 0076 0077 if (edep == 0.) return false; 0078 0079 auto newHit = new TrackerHit(); 0080 0081 newHit->SetTrackID(aStep->GetTrack()->GetTrackID()); 0082 newHit->SetEdep(edep); 0083 newHit->SetPos(aStep->GetPostStepPoint()->GetPosition()); 0084 0085 if (aStep->GetTrack()->GetTrackID() == 1 && aStep->GetTrack()->GetParentID() == 0) { 0086 newHit->SetIncidentEnergy(aStep->GetTrack()->GetVertexKineticEnergy()); 0087 } 0088 fHitsCollection->insert(newHit); 0089 // newHit->Print(); 0090 0091 return true; 0092 } 0093 0094 void TrackerSD::SetRadius(const G4double& value) 0095 { 0096 fRadius = value; 0097 } 0098 0099 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0100 0101 void TrackerSD::EndOfEvent(G4HCofThisEvent*) 0102 { 0103 G4int nofHits = fHitsCollection->entries(); 0104 0105 G4double Einc = 0; 0106 0107 /* 0108 G4cout << G4endl 0109 << "-------->Hits Collection: in this event they are " 0110 << nofHits 0111 << " hits in the target volume " << G4endl; 0112 */ 0113 0114 // PROCESSING OF MICRODOSIMETRY Y & Z SPECTRA 0115 0116 // ************************************* 0117 // Please select herebelow : 0118 // the radius of the target sphere: 0119 // variable name = radius 0120 // it is set to 5 nm by default) 0121 0122 G4double radius = fRadius; 0123 0124 // 0125 0126 //*************** 0127 // y and z 0128 //*************** 0129 0130 // Select random hit 0131 G4int randHit = 0; // Runs from 0 to number of hits - 1 0132 randHit = static_cast<G4int>(G4UniformRand() * nofHits); 0133 0134 /* 0135 G4cout 0136 << "======> random selection of hit number randHit =" 0137 << randHit << G4endl; 0138 */ 0139 0140 // Get selected random hit position 0141 G4ThreeVector hitPos = (*fHitsCollection)[randHit]->GetPos(); 0142 // G4cout << "======> random hit position x/nm =" << hitPos.x()/nm << G4endl; 0143 // G4cout << "======> random hit position y/nm =" << hitPos.y()/nm << G4endl; 0144 // G4cout << "======> random hit position z/nm =" << hitPos.z()/nm << G4endl; 0145 0146 // Set random position of center of sphere within radius 0147 G4double chord = 4. * radius / 3; 0148 G4double density = 1 * g / cm3; 0149 G4double mass = (4. / 3) * CLHEP::pi * radius * radius * radius * density; 0150 0151 // Random placement of sphere: method 1 0152 /* 0153 G4ThreeVector randDir = G4RandomDirection(); 0154 G4double randRadius = G4UniformRand()*radius; 0155 G4ThreeVector randCenterPos = randRadius*randDir + hitPos; 0156 */ 0157 0158 // Random placement of sphere: method 2 0159 0160 G4double xRand = 1.01 * radius; 0161 G4double yRand = 1.01 * radius; 0162 G4double zRand = 1.01 * radius; 0163 G4double randRad = 1.01 * radius; 0164 do { 0165 xRand = (2 * G4UniformRand() - 1) * radius; 0166 yRand = (2 * G4UniformRand() - 1) * radius; 0167 zRand = (2 * G4UniformRand() - 1) * radius; 0168 randRad = std::sqrt(xRand * xRand + yRand * yRand + zRand * zRand); 0169 } while (randRad > radius); 0170 0171 G4ThreeVector randCenterPos(xRand + hitPos.x(), yRand + hitPos.y(), zRand + hitPos.z()); 0172 0173 // Search for neighbouring hits in the sphere and cumulate deposited energy 0174 // in epsilon 0175 G4double epsilon = 0; 0176 G4int nbEdep = 0; 0177 0178 for (G4int i = 0; i < nofHits; i++) { 0179 if ((*fHitsCollection)[i]->GetIncidentEnergy() > 0) 0180 Einc = (*fHitsCollection)[i]->GetIncidentEnergy(); 0181 0182 G4ThreeVector localPos = (*fHitsCollection)[i]->GetPos(); 0183 0184 // G4cout << i << " " << (*fHitsCollection)[i] << G4endl; 0185 // G4cout << i << " " << (*fHitsCollection)[i]->GetEdep()/eV << G4endl; 0186 0187 if ((localPos.x() - randCenterPos.x()) * (localPos.x() - randCenterPos.x()) 0188 + (localPos.y() - randCenterPos.y()) * (localPos.y() - randCenterPos.y()) 0189 + (localPos.z() - randCenterPos.z()) * (localPos.z() - randCenterPos.z()) 0190 <= radius * radius) 0191 0192 { 0193 epsilon = epsilon + (*fHitsCollection)[i]->GetEdep(); 0194 nbEdep = nbEdep + 1; 0195 } 0196 } 0197 0198 // For testing only 0199 /* 0200 G4cout << "======> for hit number #" << randHit << 0201 ", we collect " 0202 << nbEdep << " energy depositions in a sphere of radius " 0203 << radius/nm << " nm and mass " 0204 << mass/kg << " kg for a total of " 0205 << epsilon/eV << " eV or " 0206 << (epsilon/joule)/(mass/kg) << " Gy" << G4endl; 0207 G4cout << "-" << G4endl; 0208 */ 0209 0210 /* 0211 FILE* myFile; 0212 myFile=fopen("yz.txt","a"); 0213 fprintf(myFile,"%e %e %e\n",radius/nm,(epsilon/eV)/(chord/nm), 0214 (epsilon/joule)/(mass/kg)); 0215 fclose(myFile); 0216 */ 0217 0218 // Get analysis manager 0219 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 0220 0221 // Fill ntuple including weighting 0222 analysisManager->FillNtupleDColumn(0, radius / nm); 0223 analysisManager->FillNtupleDColumn(2, nofHits); 0224 analysisManager->FillNtupleDColumn(3, nbEdep); 0225 analysisManager->FillNtupleDColumn(4, (epsilon / eV) / (chord / nm)); 0226 analysisManager->FillNtupleDColumn(5, (epsilon / mass) / gray); 0227 analysisManager->FillNtupleDColumn(6, Einc / eV); 0228 analysisManager->AddNtupleRow(); 0229 } 0230 0231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |