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 publication: 29 // Med. Phys. 37 (2010) 4692-4708 30 // J. Comput. Phys. 274 (2014) 841-882 31 // The Geant4-DNA web site is available at http://geant4-dna.org 32 // 33 // Authors: J. Naoki D. Kondo (UCSF, US) 34 // J. Ramos-Mendez and B. Faddegon (UCSF, US) 35 // 36 /// \file ScoreStrandBreaks.cc 37 /// \brief Implementation of the ScoreStrandBreaks class 38 39 #include "ScoreStrandBreaks.hh" 40 41 #include "G4AnalysisManager.hh" 42 #include "G4DNAChemistryManager.hh" 43 #include "G4Event.hh" 44 #include "G4EventManager.hh" 45 #include "G4Scheduler.hh" 46 #include "G4UImessenger.hh" 47 #include "G4UnitsTable.hh" 48 49 #include <G4EventManager.hh> 50 #include <G4SystemOfUnits.hh> 51 #include <globals.hh> 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 54 55 ScoreStrandBreaks::ScoreStrandBreaks(G4String name, DetectorConstruction* detect, G4double* radius) 56 : G4VPrimitiveScorer(name), G4UImessenger(), fpDetector(detect) 57 { 58 fpOutputFileUI = new G4UIcmdWithAString("/scorer/StrandBreak/OutputFile", this); 59 fpOutputFileUI->SetGuidance("Set output file name"); 60 61 fpOutputTypeUI = new G4UIcmdWithAString("/scorer/StrandBreak/OutputFormat", this); 62 fpOutputTypeUI->SetGuidance("Set output file format: ASCII by default"); 63 64 fpBreakEnergyUI = new G4UIcmdWithADoubleAndUnit("/scorer/StrandBreak/BreakEnergy", this); 65 fpBreakEnergyUI->SetDefaultUnit("eV"); 66 fpBreakEnergyUI->SetGuidance("Direct SSB energy treshold"); 67 68 fRadius = radius; 69 fpGun = new MoleculeInserter(true); 70 // G4DNAChemistryManager::Instance()->SetGun(fpGun); 71 } 72 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 74 75 ScoreStrandBreaks::~ScoreStrandBreaks() 76 { 77 delete fpOutputFileUI; 78 delete fpBreakEnergyUI; 79 delete fpOutputTypeUI; 80 delete fpGun; 81 } 82 83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 84 85 void ScoreStrandBreaks::SetNewValue(G4UIcommand* command, G4String newValue) 86 { 87 if (command == fpOutputFileUI) { 88 fOutputName = newValue; 89 } 90 91 if (command == fpBreakEnergyUI) { 92 fBreakEnergy = fpBreakEnergyUI->GetNewDoubleValue(newValue); 93 } 94 95 if (command == fpOutputTypeUI) { 96 fOutputType = newValue; 97 G4StrUtil::to_lower(fOutputType); 98 } 99 } 100 101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 102 103 G4bool ScoreStrandBreaks::ProcessHits(G4Step* aStep, G4TouchableHistory*) 104 { 105 G4double edep = aStep->GetTotalEnergyDeposit(); 106 if (edep <= 0) { 107 return FALSE; 108 } 109 fEnergyDeposit += edep; 110 111 G4StepPoint* aStepPoint = aStep->GetTrack()->GetStep()->GetPreStepPoint(); 112 G4TouchableHandle aTouchable = aStepPoint->GetTouchableHandle(); 113 G4String volumeName = aTouchable->GetVolume()->GetName(); 114 G4StrUtil::to_lower(volumeName); 115 G4int strand = 0; 116 117 if (G4StrUtil::contains(volumeName, "deoxyribose") 118 || G4StrUtil::contains(volumeName, "phosphate")) 119 { 120 if (G4StrUtil::contains(volumeName, "1")) 121 strand = 1; 122 else if (G4StrUtil::contains(volumeName, "2")) 123 strand = 2; 124 125 G4int StrandID = strand; 126 G4int BaseID = aTouchable->GetReplicaNumber(0); 127 G4int PlasmidID = aTouchable->GetReplicaNumber(1); 128 G4int EventID = G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID(); 129 130 fEnergyDepositMap[EventID][PlasmidID][BaseID][StrandID] += edep; 131 return TRUE; 132 } 133 134 if (!fDNAInserted) { 135 // Insert DNA molecules 136 std::vector<std::vector<G4int>> DNADetails = fpDetector->GetDNADetails(); 137 std::vector<G4ThreeVector> DNAPositions = fpDetector->GetDNAPositions(); 138 std::vector<G4String> DNAMolecules = fpDetector->GetDNANames(); 139 140 fpGun->Clean(); 141 for (size_t i = 0; i < DNAPositions.size(); i++) { 142 fpGun->AddMolecule(DNAMolecules[i], DNAPositions[i], 1.0 * ps); 143 } 144 fpGun->DefineTracks(); 145 146 fDNAInserted = true; 147 } 148 return FALSE; 149 } 150 151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 152 153 void ScoreStrandBreaks::Initialize(G4HCofThisEvent*) {} 154 155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 156 157 void ScoreStrandBreaks::EndOfEvent(G4HCofThisEvent*) 158 { 159 // Get EventID 160 G4int EventID = G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID(); 161 162 // Absorbed Dose 163 G4double volume = 4.0 / 3 * 3.14159 * (*fRadius) * (*fRadius) * (*fRadius) / 8; 164 G4double density = 1.0 * g / cm3; 165 G4double mass = density * volume; 166 G4double Dose = fEnergyDeposit / mass; 167 fDoseArray[EventID] = Dose / gray; 168 169 // Direct Strand Breaks 170 for (auto EventAndElse : fEnergyDepositMap) { 171 G4int event = EventAndElse.first; 172 for (auto PlasmidAndElse : EventAndElse.second) { 173 G4int plasmid = PlasmidAndElse.first; 174 for (auto BaseAndElse : PlasmidAndElse.second) { 175 G4int base = BaseAndElse.first; 176 for (auto StrandAndEdep : BaseAndElse.second) { 177 G4int strand = StrandAndEdep.first; 178 G4double edep = StrandAndEdep.second; 179 if (edep > fBreakEnergy) fDirectDamageMap[event].push_back({plasmid, base, strand}); 180 } 181 } 182 } 183 } 184 fEnergyDepositMap.clear(); 185 186 // Indirect Strand Breaks 187 std::vector<std::vector<G4int>> DNADetails = fpDetector->GetDNADetails(); 188 std::vector<G4Track*> InsertedTracks = fpGun->GetInsertedTracks(); 189 std::vector<std::vector<G4int>> IndirectBreaks; 190 191 for (size_t i = 0; i < InsertedTracks.size(); i++) { 192 if (InsertedTracks[i]->GetTrackStatus() == fStopAndKill) 193 IndirectBreaks.push_back(DNADetails[i]); 194 } 195 196 for (size_t i = 0; i < IndirectBreaks.size(); i++) { 197 fIndirectDamageMap[EventID].push_back(IndirectBreaks[i]); 198 } 199 200 fEnergyDeposit = 0; 201 fDNAInserted = false; 202 fnbOfEvents++; 203 } 204 205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 206 207 void ScoreStrandBreaks::AbsorbResultsFromWorkerScorer(G4VPrimitiveScorer* workerScorer) 208 { 209 ScoreStrandBreaks* right = 210 dynamic_cast<ScoreStrandBreaks*>(dynamic_cast<G4VPrimitiveScorer*>(workerScorer)); 211 212 if (right == 0) return; 213 if (right == this) return; 214 215 DamageMap WorkerDirectDamageMap = right->fDirectDamageMap; 216 DamageMap WorkerIndirectDamageMap = right->fIndirectDamageMap; 217 std::map<G4int, G4double> WorkerDose = right->fDoseArray; 218 219 for (auto EventAndBreaks : WorkerDirectDamageMap) { 220 G4int EventID = EventAndBreaks.first; 221 std::vector<std::vector<G4int>> Breaks = EventAndBreaks.second; 222 for (size_t i = 0; i < Breaks.size(); i++) { 223 fDirectDamageMap[EventID].push_back(Breaks[i]); 224 } 225 } 226 227 for (auto EventAndBreaks : WorkerIndirectDamageMap) { 228 G4int EventID = EventAndBreaks.first; 229 std::vector<std::vector<G4int>> Breaks = EventAndBreaks.second; 230 for (size_t i = 0; i < Breaks.size(); i++) { 231 fIndirectDamageMap[EventID].push_back(Breaks[i]); 232 } 233 } 234 235 for (auto EventAndDose : WorkerDose) { 236 G4int EventID = EventAndDose.first; 237 G4double dose = EventAndDose.second; 238 fDoseArray[EventID] = dose; 239 } 240 241 fnbOfEvents += right->fnbOfEvents; 242 right->Clear(); 243 } 244 245 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 246 247 void ScoreStrandBreaks::Clear() 248 { 249 if (fDirectDamageMap.size() != 0) fDirectDamageMap.clear(); 250 251 if (fIndirectDamageMap.size() != 0) fIndirectDamageMap.clear(); 252 253 if (fEnergyDepositMap.size() != 0) fEnergyDepositMap.clear(); 254 255 if (fDoseArray.size() != 0) fDoseArray.clear(); 256 } 257 258 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 259 260 void ScoreStrandBreaks::DrawAll() {} 261 262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 263 264 void ScoreStrandBreaks::PrintAll() {} 265 266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 267 268 void ScoreStrandBreaks::OutputAndClear(G4double LET, G4double LET_STD) 269 { 270 if (G4Threading::IsWorkerThread()) return; 271 272 if (fOutputType != "ascii") { 273 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 274 275 if (fOutputType == "csv") 276 analysisManager->SetDefaultFileType("csv"); 277 278 else if (fOutputType == "root") 279 analysisManager->SetDefaultFileType("root"); 280 281 else if (fOutputType == "xml") 282 analysisManager->SetDefaultFileType("xml"); 283 284 WriteWithAnalysisManager(analysisManager); 285 } 286 287 else 288 ASCII(LET, LET_STD); 289 290 Clear(); 291 } 292 293 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 294 295 void ScoreStrandBreaks::WriteWithAnalysisManager(G4VAnalysisManager* analysisManager) 296 { 297 analysisManager->OpenFile(fOutputName); 298 int fNtupleID = analysisManager->CreateNtuple("SB", "Direct And Indirect SBs"); 299 analysisManager->CreateNtupleIColumn(fNtupleID, "EventID"); 300 analysisManager->CreateNtupleIColumn(fNtupleID, "PlasmidID"); 301 analysisManager->CreateNtupleIColumn(fNtupleID, "StrandID"); 302 analysisManager->CreateNtupleIColumn(fNtupleID, "BaseID"); 303 analysisManager->CreateNtupleIColumn(fNtupleID, "DamageID"); 304 analysisManager->CreateNtupleIColumn(fNtupleID, "BreakID"); 305 analysisManager->FinishNtuple(fNtupleID); 306 307 for (G4int i = 0; i < fnbOfEvents; i++) { 308 if (fDirectDamageMap.find(i) != fDirectDamageMap.end()) { 309 for (size_t j = 0; j < fDirectDamageMap[i].size(); j++) { 310 std::vector<G4int> Break = fDirectDamageMap[i][j]; 311 analysisManager->FillNtupleIColumn(fNtupleID, 0, i); 312 analysisManager->FillNtupleIColumn(fNtupleID, 1, Break[0]); 313 analysisManager->FillNtupleIColumn(fNtupleID, 2, Break[1]); 314 analysisManager->FillNtupleIColumn(fNtupleID, 3, Break[2]); 315 analysisManager->FillNtupleIColumn(fNtupleID, 4, 1); 316 analysisManager->AddNtupleRow(fNtupleID); 317 } 318 } 319 320 if (fIndirectDamageMap.find(i) != fIndirectDamageMap.end()) { 321 for (size_t j = 0; j < fIndirectDamageMap[i].size(); j++) { 322 std::vector<G4int> Break = fIndirectDamageMap[i][j]; 323 analysisManager->FillNtupleIColumn(fNtupleID, 0, i); 324 analysisManager->FillNtupleIColumn(fNtupleID, 1, Break[0]); 325 analysisManager->FillNtupleIColumn(fNtupleID, 2, Break[1]); 326 analysisManager->FillNtupleIColumn(fNtupleID, 3, Break[2]); 327 analysisManager->FillNtupleIColumn(fNtupleID, 4, 2); 328 analysisManager->AddNtupleRow(fNtupleID); 329 } 330 } 331 } 332 333 analysisManager->Write(); 334 analysisManager->CloseFile(); 335 } 336 337 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 338 339 void ScoreStrandBreaks::ASCII(G4double LET, G4double LET_STD) 340 { 341 std::ofstream dna(fOutputName + ".txt"); 342 343 dna << "# DNA SB Map File" << G4endl; 344 dna << "# LET = " << LET / (keV / um) << " +- " << LET_STD / (keV / um) << " keV / um" << G4endl; 345 dna << "#" << std::setw(11) << "Event-ID" << std::setw(12) << "Plasmid-ID" << std::setw(12) 346 << "BP-ID" << std::setw(12) << "Strand-ID" << std::setw(10) << "Break-ID" << std::setw(12) 347 << "Dose (Gy) " << G4endl; 348 349 for (G4int i = 0; i < fnbOfEvents; i++) { 350 if (fDirectDamageMap.find(i) != fDirectDamageMap.end()) { 351 for (size_t j = 0; j < fDirectDamageMap[i].size(); j++) { 352 std::vector<G4int> Break = fDirectDamageMap[i][j]; 353 dna << std::setw(12) << i << std::setw(12) << Break[0] << std::setw(12) << Break[1] 354 << std::setw(12) << Break[2] << std::setw(10) << 1 << std::setw(12) << fDoseArray[i] 355 << G4endl; 356 } 357 } 358 359 if (fIndirectDamageMap.find(i) != fIndirectDamageMap.end()) { 360 for (size_t j = 0; j < fIndirectDamageMap[i].size(); j++) { 361 std::vector<G4int> Break = fIndirectDamageMap[i][j]; 362 dna << std::setw(12) << i << std::setw(12) << Break[0] << std::setw(12) << Break[1] 363 << std::setw(12) << Break[2] << std::setw(10) << 2 << std::setw(12) << fDoseArray[i] 364 << G4endl; 365 } 366 } 367 } 368 369 dna.close(); 370 } 371 372 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......