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 #include "TimeStepAction.hh" 28 29 #include "ChromosomeMapper.hh" 30 #include "DNAGeometry.hh" 31 #include "DNAHit.hh" 32 #include "DetectorConstruction.hh" 33 #include "EventAction.hh" 34 #include "OctreeNode.hh" 35 36 #include "G4DNAMolecularReactionTable.hh" 37 #include "G4ITTrackHolder.hh" 38 #include "G4ITTrackingManager.hh" 39 #include "G4Molecule.hh" 40 #include "G4RunManager.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 TimeStepAction::TimeStepAction(EventAction* event) 44 : G4UserTimeStepAction(), 45 fEventAction(event), 46 fRadicalKillDistance(4.5 * nm), 47 fpChemistryTrackHolder(G4ITTrackHolder::Instance()) 48 { 49 AddTimeStep(1 * picosecond, 0.5 * nanosecond); 50 // ctor 51 } 52 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 54 55 TimeStepAction::~TimeStepAction() = default; 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 57 58 void TimeStepAction::StartProcessing() 59 { 60 auto det = dynamic_cast<const DetectorConstruction*>( 61 G4RunManager::GetRunManager()->GetUserDetectorConstruction()); 62 fDNAGeometry = det->GetDNAGeometry(); 63 if (fDNAGeometry == nullptr) { 64 G4ExceptionDescription exceptionDescription; 65 exceptionDescription << "fDNAGeometry is null"; 66 G4Exception( 67 "TimeStepAction" 68 "StartProcessing", 69 "no fDNAGeometry", FatalException, exceptionDescription); 70 } 71 72 fRadicalKillDistance = fDNAGeometry->GetRadicalKillDistance(); 73 } 74 75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 76 77 void TimeStepAction::UserPostTimeStepAction() 78 { 79 RadicalKillDistance(); 80 } 81 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 83 84 void TimeStepAction::UserPreTimeStepAction() {} 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 88 void TimeStepAction::UserReactionAction(const G4Track&, const G4Track&, 89 const std::vector<G4Track*>*) 90 {} 91 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 93 94 void TimeStepAction::RadicalKillDistance() 95 { 96 if (fpChemistryTrackHolder == nullptr) { 97 G4ExceptionDescription exceptionDescription; 98 exceptionDescription << "fpChemistryTrackHolder is null"; 99 G4Exception( 100 "TimeStepAction" 101 "RadicalKillDistance", 102 "NO_fpChemistryTrackHolder", FatalException, exceptionDescription); 103 } 104 G4Track* trackToKill; 105 G4TrackManyList::iterator it_begin = fpChemistryTrackHolder->GetMainList()->begin(); 106 G4TrackManyList::iterator it_end = fpChemistryTrackHolder->GetMainList()->end(); 107 while (it_begin != it_end) { 108 trackToKill = nullptr; 109 110 const G4VTouchable* touchable = it_begin->GetTouchable(); 111 if (touchable == nullptr) { 112 ++it_begin; 113 continue; 114 } 115 116 G4LogicalVolume* trackLV = touchable->GetVolume()->GetLogicalVolume(); 117 G4LogicalVolume* motherLV = touchable->GetVolume()->GetMotherLogical(); 118 119 OctreeNode* octree_track = fDNAGeometry->GetTopOctreeNode(trackLV); 120 OctreeNode* octree_mother = fDNAGeometry->GetTopOctreeNode(motherLV); 121 122 if ((octree_track == nullptr) && (octree_mother == nullptr)) { 123 trackToKill = *it_begin; 124 } 125 else if (octree_track != nullptr) { 126 const G4AffineTransform& trans = it_begin->GetTouchable()->GetHistory()->GetTopTransform(); 127 G4ThreeVector pos = trans.TransformPoint(it_begin->GetPosition()); 128 size_t n = octree_track->SearchOctree(pos, fRadicalKillDistance).size(); 129 if (n == 0) { 130 trackToKill = *it_begin; 131 } 132 if (fDNAGeometry->IsInsideHistone(trackLV, pos)) { 133 trackToKill = *it_begin; 134 } 135 } 136 else { 137 const G4AffineTransform& trans = it_begin->GetTouchable()->GetHistory()->GetTopTransform(); 138 G4ThreeVector pos = trans.TransformPoint(it_begin->GetPosition()); 139 if (fDNAGeometry->IsInsideHistone(trackLV, pos)) { 140 trackToKill = *it_begin; 141 } 142 } 143 ++it_begin; 144 if (trackToKill != nullptr) { 145 fpChemistryTrackHolder->PushToKill(trackToKill); 146 } 147 } 148 } 149 150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 151