#include <G4DNASecondOrderReaction.hh>
Inheritance diagram for G4DNASecondOrderReaction:
Definition at line 46 of file G4DNASecondOrderReaction.hh.
G4DNASecondOrderReaction::G4DNASecondOrderReaction | ( | const G4String & | aName = "DNASecondOrderReaction" , |
|
G4ProcessType | type = fDecay | |||
) |
Definition at line 62 of file G4DNASecondOrderReaction.cc.
00062 : 00063 G4VITProcess(aName,type), 00064 InitProcessState(fpSecondOrderReactionState, fpState) 00065 { 00066 Create(); 00067 }
G4DNASecondOrderReaction::~G4DNASecondOrderReaction | ( | ) | [virtual] |
G4DNASecondOrderReaction::G4DNASecondOrderReaction | ( | const G4DNASecondOrderReaction & | ) |
Definition at line 69 of file G4DNASecondOrderReaction.cc.
00069 : 00070 G4VITProcess(rhs), 00071 InitProcessState(fpSecondOrderReactionState, fpState) 00072 { 00073 Create(); 00074 }
virtual G4VParticleChange* G4DNASecondOrderReaction::AlongStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
virtual G4double G4DNASecondOrderReaction::AlongStepGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4double | , | |||
G4double | , | |||
G4double & | , | |||
G4GPILSelection * | ||||
) | [inline, virtual] |
virtual G4VParticleChange* G4DNASecondOrderReaction::AtRestDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
virtual G4double G4DNASecondOrderReaction::AtRestGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [inline, virtual] |
void G4DNASecondOrderReaction::BuildPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
Reimplemented from G4VITProcess.
Definition at line 94 of file G4DNASecondOrderReaction.cc.
References fIsInitialized, fMolarMassOfMaterial, fpMaterial, fpMoleculeDensity, G4Material::GetMassOfMolecule(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), and G4DNAMolecularMaterial::Instance().
00095 { 00096 fpMoleculeDensity = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(fpMaterial); 00097 fMolarMassOfMaterial = fpMaterial->GetMassOfMolecule()*CLHEP::Avogadro*1e3; 00098 fIsInitialized = true; 00099 }
G4DNASecondOrderReaction & G4DNASecondOrderReaction::operator= | ( | const G4DNASecondOrderReaction & | ) |
Definition at line 80 of file G4DNASecondOrderReaction.cc.
00081 { 00082 if (this == &rhs) return *this; // handle self assignment 00083 00084 //assignment operator 00085 return *this; 00086 }
G4VParticleChange * G4DNASecondOrderReaction::PostStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 227 of file G4DNASecondOrderReaction.cc.
References G4DNADamages::AddIndirectDamage(), DBL_MAX, fParticleChange, fpMaterial, fReturnedValue, fStopAndKill, G4BestUnit, G4cout, G4endl, G4Track::GetGlobalTime(), GetMolecule(), G4Material::GetName(), G4Molecule::GetName(), G4Track::GetPosition(), G4ParticleChange::Initialize(), G4DNADamages::Instance(), G4ITTrackHolder::Instance(), G4VParticleChange::ProposeTrackStatus(), and G4VProcess::verboseLevel.
00228 { 00229 G4Molecule* molecule = GetMolecule(track); 00230 #ifdef G4VERBOSE 00231 if(verboseLevel > 1) 00232 { 00233 G4cout << "___________" << G4endl; 00234 G4cout << ">>> Beginning of G4DNASecondOrderReaction verbose" << G4endl; 00235 G4cout << ">>> Returned value : " << G4BestUnit(fReturnedValue,"Time") << G4endl; 00236 G4cout << ">>> Time Step : " << G4BestUnit(G4ITTrackHolder::Instance()->GetTimeStep(),"Time") << G4endl; 00237 G4cout << ">>> Reaction : " << molecule->GetName() << " + " << fpMaterial->GetName() << G4endl; 00238 G4cout << ">>> End of G4DNASecondOrderReaction verbose <<<" << G4endl; 00239 } 00240 #endif 00241 fReturnedValue = DBL_MAX; 00242 fParticleChange.Initialize(track); 00243 fParticleChange.ProposeTrackStatus(fStopAndKill); 00244 G4DNADamages::Instance()->AddIndirectDamage(fpMaterial->GetName(),molecule,track.GetPosition(),track.GetGlobalTime()); 00245 return &fParticleChange; 00246 }
G4double G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4ForceCondition * | condition | |||
) | [virtual] |
Implements G4VProcess.
Definition at line 126 of file G4DNASecondOrderReaction.cc.
References DBL_MAX, G4DynamicParticle::DumpInfo(), fConcentration, G4DNASecondOrderReaction::SecondOrderReactionState::fIsInGoodMaterial, fpMolecularConfiguration, G4DNASecondOrderReaction::SecondOrderReactionState::fPreviousTimeAtPreStepPoint, fpSecondOrderReactionState, G4VITProcess::fpState, fReactionRate, fReturnedValue, G4cout, G4endl, G4Track::GetCurrentStepNumber(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4Material::GetIndex(), G4Track::GetMaterial(), G4Molecule::GetMolecularConfiguration(), GetMolecule(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, G4VITProcess::ResetNumberOfInteractionLengthLeft(), G4VITProcess::SubtractNumberOfInteractionLengthLeft(), and G4VProcess::verboseLevel.
00129 { 00130 // G4cout << "G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength" << G4endl; 00131 // G4cout << "For reaction : " << fpMaterial->GetName() << " + " << fpMolecularConfiguration->GetName() << G4endl; 00132 00133 //_______________________________________________________________________ 00134 // Check whether the track is in the good material (maybe composite material) 00135 const G4Material* material = track.GetMaterial(); 00136 00137 G4Molecule* mol = GetMolecule(track); 00138 if(!mol) return DBL_MAX; 00139 if(mol->GetMolecularConfiguration() != fpMolecularConfiguration) 00140 { 00141 // G4cout <<"mol->GetMolecularConfiguration() != fpMolecularConfiguration" << G4endl; 00142 return DBL_MAX; 00143 } 00144 00145 G4double molDensity = (*fpMoleculeDensity)[material->GetIndex()]; 00146 00147 if(molDensity == 0.0) // ie : not found 00148 { 00149 if(fpSecondOrderReactionState->fIsInGoodMaterial) 00150 { 00151 ResetNumberOfInteractionLengthLeft(); 00152 fpSecondOrderReactionState->fIsInGoodMaterial = false; 00153 } 00154 00155 // G4cout << " Material " << fpMaterial->GetName() << " not found " 00156 // <<" | name of current material : " << material->GetName() 00157 // << G4endl; 00158 00159 return DBL_MAX; // Becareful return here !! 00160 } 00161 00162 // G4cout << " Va calculer le temps d'interaction " << G4endl; 00163 00164 fpSecondOrderReactionState->fIsInGoodMaterial = true; 00165 00166 // fConcentration = molDensity/fMolarMassOfMaterial; 00167 fConcentration = molDensity/CLHEP::Avogadro; 00168 00169 // G4cout << "Concentration : " << fConcentration / (g/mole)<< G4endl; 00170 00171 //_______________________________________________________________________ 00172 // Either initialize the lapse of time left 00173 // meaning => the track enters for the first time in the material 00174 // or substract the previous time step to the previously calculated lapse of time left 00175 // meaning => the track has not left this material since the previous call 00176 00177 G4double previousTimeStep(-1.); 00178 00179 if(track.GetCurrentStepNumber() > 0) 00180 previousTimeStep = track.GetGlobalTime() - fpSecondOrderReactionState->fPreviousTimeAtPreStepPoint ; 00181 00182 fpSecondOrderReactionState->fPreviousTimeAtPreStepPoint = track.GetGlobalTime(); 00183 00184 if ( (previousTimeStep < 0.0) || (fpState->theNumberOfInteractionLengthLeft<=0.0)) { 00185 // beggining of tracking (or just after DoIt of this process) 00186 ResetNumberOfInteractionLengthLeft(); 00187 } else if ( previousTimeStep > 0.0) { 00188 // subtract NumberOfInteractionLengthLeft 00189 SubtractNumberOfInteractionLengthLeft(previousTimeStep ); 00190 } else { 00191 // zero step 00192 // DO NOTHING 00193 } 00194 00195 // condition is set to "Not Forced" 00196 *pForceCond = NotForced; 00197 00198 // get mean free path 00199 fpState->currentInteractionLength = 1/(fReactionRate*fConcentration); 00200 00201 G4double value; 00202 if (fpState->currentInteractionLength <DBL_MAX) { 00203 value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength); 00204 } else { 00205 value = DBL_MAX; 00206 } 00207 #ifdef G4VERBOSE 00208 if (verboseLevel>2){ 00209 G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength "; 00210 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00211 track.GetDynamicParticle()->DumpInfo(); 00212 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00213 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl; 00214 } 00215 #endif 00216 00217 // G4cout << "currentInteractionLength : " << fpState->currentInteractionLength << G4endl; 00218 // G4cout << "Returned time : " << G4BestUnit(value,"Time") << G4endl; 00219 00220 if(value < fReturnedValue) 00221 fReturnedValue = value; 00222 00223 return value*-1; 00224 // multiple by -1 to indicate to the tracking system that we are returning a time 00225 }
void G4DNASecondOrderReaction::SetReaction | ( | const G4MolecularConfiguration * | , | |
const G4Material * | , | |||
double | ||||
) |
Definition at line 110 of file G4DNASecondOrderReaction.cc.
References FatalErrorInArgument, fIsInitialized, fpMaterial, fpMolecularConfiguration, fReactionRate, and G4Exception().
00112 { 00113 if(fIsInitialized) 00114 { 00115 G4ExceptionDescription exceptionDescription ; 00116 exceptionDescription << "G4DNASecondOrderReaction was already initialised. "; 00117 exceptionDescription << "You cannot set a reaction after initialisation."; 00118 G4Exception("G4DNASecondOrderReaction::SetReaction","G4DNASecondOrderReaction001", 00119 FatalErrorInArgument,exceptionDescription); 00120 } 00121 fpMolecularConfiguration = molConf; 00122 fpMaterial = mat; 00123 fReactionRate = reactionRate; 00124 }
void G4DNASecondOrderReaction::StartTracking | ( | G4Track * | ) | [virtual] |
Reimplemented from G4VITProcess.
Definition at line 102 of file G4DNASecondOrderReaction.cc.
References G4VITProcess::fpState, G4VITProcess::StartTracking(), and G4VProcess::StartTracking().
00103 { 00104 G4VProcess::StartTracking(track); 00105 G4VITProcess::fpState = new SecondOrderReactionState(); 00106 G4VITProcess::StartTracking(track); 00107 }
G4double G4DNASecondOrderReaction::fConcentration [protected] |
Definition at line 118 of file G4DNASecondOrderReaction.hh.
Referenced by PostStepGetPhysicalInteractionLength().
G4bool G4DNASecondOrderReaction::fIsInitialized [protected] |
Definition at line 112 of file G4DNASecondOrderReaction.hh.
Referenced by BuildPhysicsTable(), and SetReaction().
const G4Material* G4DNASecondOrderReaction::fpMaterial [protected] |
Definition at line 124 of file G4DNASecondOrderReaction.hh.
Referenced by BuildPhysicsTable(), PostStepDoIt(), and SetReaction().
const G4MolecularConfiguration* G4DNASecondOrderReaction::fpMolecularConfiguration [protected] |
Definition at line 123 of file G4DNASecondOrderReaction.hh.
Referenced by PostStepGetPhysicalInteractionLength(), and SetReaction().
const std::vector<double>* G4DNASecondOrderReaction::fpMoleculeDensity [protected] |
Definition at line 120 of file G4DNASecondOrderReaction.hh.
Referenced by PostStepGetPhysicalInteractionLength().
G4double G4DNASecondOrderReaction::fReactionRate [protected] |
Definition at line 117 of file G4DNASecondOrderReaction.hh.
Referenced by PostStepGetPhysicalInteractionLength(), and SetReaction().
G4double G4DNASecondOrderReaction::fReturnedValue [protected] |
Definition at line 114 of file G4DNASecondOrderReaction.hh.
Referenced by PostStepDoIt(), and PostStepGetPhysicalInteractionLength().