Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMolecularReaction.cc
Go to the documentation of this file.
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// Author: Mathieu Karamitros ([email protected])
28//
29// WARNING : This class is released as a prototype.
30// It might strongly evolve or even disapear in the next releases.
31//
32// History:
33// -----------
34// 10 Oct 2011 M.Karamitros created
35//
36// -------------------------------------------------------------------
37
42#include "G4Molecule.hh"
43#include "G4MoleculeFinder.hh"
44#include "G4ITReactionChange.hh"
45#include "G4ITReaction.hh"
46
47#include "G4ITTrackHolder.hh"
48
51 , fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
52 , fpReactionModel(nullptr)
53{
54}
55
58{
59 fpReactionModel = pReactionModel;
60}
61
63 const G4Track &trackB,
64 double currentStepTime,
65 bool userStepTimeLimit) /*const*/
66{
67 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
68 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
69
70 const G4double reactionRadius = fpReactionModel->GetReactionRadius(pMoleculeA, pMoleculeB);
71
72 G4double separationDistance = -1.;
73
74 if (currentStepTime == 0.)
75 {
76 userStepTimeLimit = false;
77 }
78
79 G4bool output = fpReactionModel->FindReaction(trackA, trackB, reactionRadius,
80 separationDistance, userStepTimeLimit);
81 return output;
82}
83
84std::unique_ptr<G4ITReactionChange> G4DNAMolecularReaction::MakeReaction(const G4Track &trackA,
85 const G4Track &trackB)
86{
87 std::unique_ptr<G4ITReactionChange> pChanges(new G4ITReactionChange());
88 pChanges->Initialize(trackA, trackB);
89
90 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
91 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
92
93 const auto pReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
94
95 const G4int nbProducts = pReactionData->GetNbProducts();
96
97 if (nbProducts)
98 {
99 const G4double D1 = pMoleculeA->GetDiffusionCoefficient();
100 const G4double D2 = pMoleculeB->GetDiffusionCoefficient();
101 const G4double sqrD1 = D1 == 0. ? 0. : std::sqrt(D1);
102 const G4double sqrD2 = D2 == 0. ? 0. : std::sqrt(D2);
103 const G4double inv_numerator = 1./(sqrD1 + sqrD2);
104 const G4ThreeVector reactionSite = sqrD2 * inv_numerator * trackA.GetPosition()
105 + sqrD1 * inv_numerator * trackB.GetPosition();
106
107 for (G4int j = 0; j < nbProducts; ++j)
108 {
109 auto pProduct = new G4Molecule(pReactionData->GetProduct(j));
110 auto pProductTrack = pProduct->BuildTrack(trackA.GetGlobalTime(), reactionSite);
111
112 pProductTrack->SetTrackStatus(fAlive);
113
114 G4ITTrackHolder::Instance()->Push(pProductTrack);
115
116 pChanges->AddSecondary(pProductTrack);
117 G4MoleculeFinder::Instance()->Push(pProductTrack);
118 }
119 }
120
121 pChanges->KillParents(true);
122 return pChanges;
123}
124
126{
127 fpReactionModel = pReactionModel;
128}
129
130std::vector<std::unique_ptr<G4ITReactionChange>> G4DNAMolecularReaction::FindReaction(
131 G4ITReactionSet* pReactionSet,
132 const double currentStepTime,
133 const double /*fGlobalTime*/,
134 const bool reachedUserStepTimeLimit)
135{
136 std::vector<std::unique_ptr<G4ITReactionChange>> fReactionInfo;
137 fReactionInfo.clear();
138
139 if (pReactionSet == nullptr)
140 {
141 return fReactionInfo;
142 }
143
144 G4ITReactionPerTrackMap& reactionPerTrackMap = pReactionSet->GetReactionMap();
145 for (auto tracks_i = reactionPerTrackMap.begin();
146 tracks_i != reactionPerTrackMap.end();
147 tracks_i = reactionPerTrackMap.begin())
148 {
149 G4Track* pTrackA = tracks_i->first;
150 if (pTrackA->GetTrackStatus() == fStopAndKill)
151 {
152 continue;
153 }
154
155 G4ITReactionPerTrackPtr reactionPerTrack = tracks_i->second;
156 G4ITReactionList& reactionList = reactionPerTrack->GetReactionList();
157
158 assert(reactionList.begin() != reactionList.end());
159
160 for (auto it = reactionList.begin(); it != reactionList.end(); it = reactionList.begin())
161 {
162 G4ITReactionPtr reaction(*it);
163 G4Track* pTrackB = reaction->GetReactant(pTrackA);
164 if (pTrackB->GetTrackStatus() == fStopAndKill)
165 {
166 continue;
167 }
168
169 if (pTrackB == pTrackA)
170 {
171 G4ExceptionDescription exceptionDescription;
172 exceptionDescription
173 << "The IT reaction process sent back a reaction between trackA and trackB. ";
174 exceptionDescription << "The problem is trackA == trackB";
175 G4Exception("G4ITModelProcessor::FindReaction",
176 "ITModelProcessor005",
178 exceptionDescription);
179 }
180
181 pReactionSet->SelectThisReaction(reaction);
182
183 if (TestReactibility(*pTrackA, *pTrackB, currentStepTime, reachedUserStepTimeLimit))
184 {
185 auto pReactionChange = MakeReaction(*pTrackA, *pTrackB);
186
187 if (pReactionChange)
188 {
189 fReactionInfo.push_back(std::move(pReactionChange));
190 break;
191 }
192 }
193 }
194 }
195
196 pReactionSet->CleanAllReaction();
197 return fReactionInfo;
198}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
Definition: G4ITReaction.hh:61
std::list< G4ITReactionPtr > G4ITReactionList
Definition: G4ITReaction.hh:63
G4shared_ptr< G4ITReaction > G4ITReactionPtr
Definition: G4ITReaction.hh:60
std::map< G4Track *, G4ITReactionPerTrackPtr, compTrackPerID > G4ITReactionPerTrackMap
Definition: G4ITReaction.hh:66
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
ReturnType & reference_cast(OriginalType &source)
@ fAlive
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Data * GetReactionData(Reactant *, Reactant *) const
void SetReactionModel(G4VDNAReactionModel *)
const G4DNAMolecularReactionTable *& fMolReactionTable
std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &) override
std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction(G4ITReactionSet *, const double, const double, const bool) override
G4VDNAReactionModel * fpReactionModel
G4bool TestReactibility(const G4Track &, const G4Track &, double currentStepTime, bool userStepTimeLimit) override
virtual void Push(G4Track *track)
static G4ITFinder * Instance()
void SelectThisReaction(G4ITReactionPtr reaction)
void CleanAllReaction()
G4ITReactionPerTrackMap & GetReactionMap()
virtual void Push(G4Track *)
static G4ITTrackHolder * Instance()
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
G4TrackStatus GetTrackStatus() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
virtual G4bool FindReaction(const G4Track &, const G4Track &, G4double, G4double &, G4bool)=0
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0