File indexing completed on 2025-12-16 09:29:05
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
0028
0029
0030 #include "TimeStepAction.hh"
0031 #include "G4UnitsTable.hh"
0032 #include "G4SystemOfUnits.hh"
0033
0034 #include "G4DNAChemistryManager.hh"
0035 #include "G4MoleculeCounter.hh"
0036 #include "G4MoleculeTable.hh"
0037 #include "G4Scheduler.hh"
0038 #include "UserMolecule.hh"
0039 #include "G4ChemTimeStepModel.hh"
0040 #include "G4EmParameters.hh"
0041 #include "Analysis.hh"
0042 #include "G4DNAMolecule.hh"
0043 #include "G4H3O.hh"
0044 #include "G4OH.hh"
0045 #include "G4H2O2.hh"
0046 #include "G4H2.hh"
0047 #include "G4Hydrogen.hh"
0048 #include "G4Electron_aq.hh"
0049 #include "G4ParticleTable.hh"
0050 #include "G4ParticleDefinition.hh"
0051
0052
0053 TimeStepAction::TimeStepAction()
0054 : G4UserTimeStepAction()
0055 {
0056 auto param = G4EmParameters::Instance();
0057 if (param->GetTimeStepModel() == G4ChemTimeStepModel::SBS) {
0058 AddTimeStep(1*picosecond,0.1*picosecond);
0059 AddTimeStep(10*picosecond,1*picosecond);
0060 AddTimeStep(100*picosecond,3*picosecond);
0061 AddTimeStep(1000*picosecond,10*picosecond);
0062 AddTimeStep(10000*picosecond,100*picosecond);
0063 }
0064 }
0065
0066
0067
0068
0069 TimeStepAction::TimeStepAction(const TimeStepAction& other)
0070 : G4UserTimeStepAction(other)
0071 {}
0072
0073
0074
0075 TimeStepAction& TimeStepAction::operator=(const TimeStepAction& rhs)
0076 {
0077 if (this == &rhs) return *this;
0078
0079 return *this;
0080 }
0081
0082
0083
0084 void TimeStepAction::StartProcessing()
0085 {
0086
0087 G4cout<<"Start chemistry "<<G4Scheduler::Instance()->GetGlobalTime()/s<<" s"<<G4endl;
0088 G4cout<<"Chemical endTime "<<G4Scheduler::Instance()->GetEndTime()/nanosecond<<" ns"<<G4endl;
0089 }
0090
0091
0092
0093 void TimeStepAction::UserReactionAction(const G4Track& a, const G4Track& b,
0094 const std::vector<G4Track*>* products)
0095 {
0096
0097
0098
0099
0100
0101
0102
0103 fReactif1 = 0;
0104 fReactif2 = 0;
0105 fProduct1 = 0;
0106 fProduct2 = 0;
0107
0108 fReactif1 = SetParticleFlag(a.GetDynamicParticle()->GetDefinition());
0109 fReactif2 = SetParticleFlag(b.GetDynamicParticle()->GetDefinition());
0110 if(products)
0111 {
0112 for(unsigned int i=0;i<products->size();++i)
0113 {
0114 switch(i)
0115 {
0116 case 0:
0117 fProduct1 = SetParticleFlag(((*products)[i])->GetDynamicParticle()->GetDefinition());
0118 break;
0119 case 1:
0120 fProduct2 = SetParticleFlag(((*products)[i])->GetDynamicParticle()->GetDefinition());;
0121 break;
0122 default:
0123 break;
0124 }
0125 }
0126 }
0127
0128
0129
0130
0131
0132 G4int base;
0133
0134 if( (fReactif1==3 || fReactif2 ==3)
0135 && (
0136 (fReactif1==8 || fReactif2 ==8)
0137 || (fReactif1==9 || fReactif2 ==9)
0138 || (fReactif1==10|| fReactif2 ==10)
0139 || (fReactif1==11|| fReactif2 ==11)
0140 || (fReactif1==12|| fReactif2 ==12)
0141 || (fReactif1==13|| fReactif2 ==13)
0142 ) )
0143 {
0144
0145 const G4Track* dnaMolecule = nullptr;
0146 const G4Track* radical = nullptr;
0147
0148 if (a.GetDynamicParticle()->GetDefinition() == G4OH::Definition()) {
0149 dnaMolecule = &b;
0150 radical = &a;
0151 }
0152 else {
0153 dnaMolecule = &a;
0154 radical = &b;
0155 }
0156
0157 auto _dnaMolecule = dynamic_cast<UserMolecule*> (dnaMolecule->GetUserInformation());
0158 if (_dnaMolecule) {
0159 G4int copyNo = _dnaMolecule->GetCopyNumber();
0160 G4int strand = _dnaMolecule->GetStrand();
0161
0162 if((fReactif1==8 || fReactif2 ==8)
0163 || (fReactif1==9 || fReactif2 ==9))
0164 {
0165 base=0;
0166 }
0167 else
0168 {
0169 base=1;
0170 }
0171
0172
0173 auto analysisManager = Analysis::GetAnalysis()->GetAnalysisManager();
0174 analysisManager->FillNtupleIColumn(1, 0, strand);
0175 analysisManager->FillNtupleIColumn(1, 1, copyNo);
0176 analysisManager->FillNtupleDColumn(1, 2, radical->GetPosition().getX()/nm);
0177 analysisManager->FillNtupleDColumn(1, 3, radical->GetPosition().getY()/nm);
0178 analysisManager->FillNtupleDColumn(1, 4, radical->GetPosition().getZ()/nm);
0179 analysisManager->FillNtupleDColumn(1, 5, G4Scheduler::Instance()->GetGlobalTime()/ns );
0180 analysisManager->FillNtupleIColumn(1, 6, base);
0181 analysisManager->AddNtupleRow(1);
0182 } else {
0183 G4String msg = "Error in Downcasting";
0184 G4Exception("TimeStepAction::UserReactionAction", "", FatalException, msg);
0185 }
0186
0187 }
0188
0189 }
0190
0191
0192
0193 G4int TimeStepAction::SetParticleFlag(const G4ParticleDefinition* partDef)
0194 {
0195 G4int flag (0);
0196
0197 if(partDef == G4H3O::Definition()) flag = 1;
0198 else if(partDef == G4OH::Definition()) flag = 3;
0199 else if(partDef == G4Electron_aq::Definition()) flag = 4;
0200 else if(partDef == G4Hydrogen::Definition()) flag = 5;
0201 else if(partDef == G4H2::Definition()) flag = 6;
0202 else if(partDef == G4H2O2::Definition()) flag = 7;
0203 else if(partDef == G4Deoxyribose::Definition()) flag = 8;
0204 else if(partDef == G4Phosphate::Definition()) flag = 9;
0205 else if(partDef == G4Adenine::Definition()) flag = 10;
0206 else if(partDef == G4Thymine::Definition()) flag = 11;
0207 else if(partDef == G4Guanine::Definition()) flag = 12;
0208 else if(partDef == G4Cytosine::Definition()) flag = 13;
0209 else if(partDef == G4Histone::Definition()) flag = 14;
0210 else if(partDef == G4DamagedDeoxyribose::Definition()) flag = 15;
0211 else if(partDef == G4DamagedAdenine::Definition()) flag = 16;
0212 else if(partDef == G4DamagedThymine::Definition()) flag = 17;
0213 else if(partDef == G4DamagedCytosine::Definition()) flag = 18;
0214 else if(partDef == G4DamagedGuanine::Definition()) flag = 19;
0215 else {
0216 G4ParticleDefinition* OHm =G4ParticleTable::GetParticleTable()->FindParticle("OHm");
0217 if (OHm && partDef == OHm) {
0218 flag = 2;
0219 }
0220 }
0221 return flag;
0222 }
0223
0224