Warning, file /geant4/examples/advanced/dna/moleculardna/src/TimeStepAction.cc was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027 #include "TimeStepAction.hh"
0028
0029 #include "ChromosomeMapper.hh"
0030 #include "DNAGeometry.hh"
0031 #include "DNAHit.hh"
0032 #include "DetectorConstruction.hh"
0033 #include "EventAction.hh"
0034 #include "OctreeNode.hh"
0035
0036 #include "G4DNAMolecularReactionTable.hh"
0037 #include "G4ITTrackHolder.hh"
0038 #include "G4ITTrackingManager.hh"
0039 #include "G4Molecule.hh"
0040 #include "G4RunManager.hh"
0041 #include "G4Scheduler.hh"
0042
0043
0044 TimeStepAction::TimeStepAction(EventAction* event)
0045 : G4UserTimeStepAction(),
0046 fEventAction(event),
0047 fRadicalKillDistance(4.5 * nm),
0048 fpChemistryTrackHolder(G4ITTrackHolder::Instance())
0049 {
0050
0051
0052 }
0053
0054
0055
0056 TimeStepAction::~TimeStepAction() = default;
0057
0058
0059 void TimeStepAction::StartProcessing()
0060 {
0061 auto det = dynamic_cast<const DetectorConstruction*>(
0062 G4RunManager::GetRunManager()->GetUserDetectorConstruction());
0063 fDNAGeometry = det->GetDNAGeometry();
0064 if (fDNAGeometry == nullptr) {
0065 G4ExceptionDescription exceptionDescription;
0066 exceptionDescription << "fDNAGeometry is null";
0067 G4Exception(
0068 "TimeStepAction"
0069 "StartProcessing",
0070 "no fDNAGeometry", FatalException, exceptionDescription);
0071 }
0072
0073 fRadicalKillDistance = fDNAGeometry->GetRadicalKillDistance();
0074 }
0075
0076
0077
0078 void TimeStepAction::UserPostTimeStepAction()
0079 {
0080 RadicalKillDistance();
0081 }
0082
0083
0084
0085 void TimeStepAction::UserPreTimeStepAction() {}
0086
0087
0088
0089 void TimeStepAction::UserReactionAction(const G4Track&, const G4Track&,
0090 const std::vector<G4Track*>*)
0091 {}
0092
0093
0094
0095 void TimeStepAction::RadicalKillDistance()
0096 {
0097 if (fpChemistryTrackHolder == nullptr) {
0098 G4ExceptionDescription exceptionDescription;
0099 exceptionDescription << "fpChemistryTrackHolder is null";
0100 G4Exception(
0101 "TimeStepAction"
0102 "RadicalKillDistance",
0103 "NO_fpChemistryTrackHolder", FatalException, exceptionDescription);
0104 }
0105 G4Track* trackToKill;
0106 G4TrackManyList::iterator it_begin = fpChemistryTrackHolder->GetMainList()->begin();
0107 G4TrackManyList::iterator it_end = fpChemistryTrackHolder->GetMainList()->end();
0108 while (it_begin != it_end) {
0109 trackToKill = nullptr;
0110
0111 const G4VTouchable* touchable = it_begin->GetTouchable();
0112 if (touchable == nullptr) {
0113 ++it_begin;
0114 continue;
0115 }
0116
0117 G4LogicalVolume* trackLV = touchable->GetVolume()->GetLogicalVolume();
0118 G4LogicalVolume* motherLV = touchable->GetVolume()->GetMotherLogical();
0119
0120 OctreeNode* octree_track = fDNAGeometry->GetTopOctreeNode(trackLV);
0121 OctreeNode* octree_mother = fDNAGeometry->GetTopOctreeNode(motherLV);
0122
0123 if ((octree_track == nullptr) && (octree_mother == nullptr)) {
0124 trackToKill = *it_begin;
0125 }
0126 else if (octree_track != nullptr) {
0127 const G4AffineTransform& trans = it_begin->GetTouchable()->GetHistory()->GetTopTransform();
0128 G4ThreeVector pos = trans.TransformPoint(it_begin->GetPosition());
0129 size_t n = octree_track->SearchOctree(pos, fRadicalKillDistance).size();
0130 if (n == 0) {
0131 trackToKill = *it_begin;
0132 }
0133 if (fDNAGeometry->IsInsideHistone(trackLV, pos)) {
0134 trackToKill = *it_begin;
0135 }
0136 }
0137 else {
0138 const G4AffineTransform& trans = it_begin->GetTouchable()->GetHistory()->GetTopTransform();
0139 G4ThreeVector pos = trans.TransformPoint(it_begin->GetPosition());
0140 if (fDNAGeometry->IsInsideHistone(trackLV, pos)) {
0141 trackToKill = *it_begin;
0142 }
0143 }
0144 ++it_begin;
0145 if (trackToKill != nullptr) {
0146 fpChemistryTrackHolder->PushToKill(trackToKill);
0147 G4Scheduler::Instance()->SetInteractionStep(true);
0148 }
0149 }
0150 }
0151
0152