Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMolecularReaction Class Reference

#include <G4DNAMolecularReaction.hh>

+ Inheritance diagram for G4DNAMolecularReaction:

Public Member Functions

 G4DNAMolecularReaction ()
 
 G4DNAMolecularReaction (G4VDNAReactionModel *)
 
 ~G4DNAMolecularReaction () override=default
 
 G4DNAMolecularReaction (const G4DNAMolecularReaction &other)=delete
 
G4DNAMolecularReactionoperator= (const G4DNAMolecularReaction &other)=delete
 
G4bool TestReactibility (const G4Track &, const G4Track &, double currentStepTime, bool userStepTimeLimit) override
 
std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction (G4ITReactionSet *, const double, const double, const bool) override
 
std::unique_ptr< G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &) override
 
void SetReactionModel (G4VDNAReactionModel *)
 
- Public Member Functions inherited from G4VITReactionProcess
 G4VITReactionProcess ()=default
 
virtual ~G4VITReactionProcess ()=default
 
 G4VITReactionProcess (const G4VITReactionProcess &other)=delete
 
G4VITReactionProcessoperator= (const G4VITReactionProcess &other)=delete
 
virtual void Initialize ()
 
virtual G4bool IsApplicable (const G4ITType &, const G4ITType &) const
 
virtual G4bool TestReactibility (const G4Track &, const G4Track &, double, bool)=0
 
virtual std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction (G4ITReactionSet *, const double, const double, const bool)=0
 
virtual std::unique_ptr< G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &)=0
 
virtual void SetReactionTable (const G4ITReactionTable *)
 

Protected Attributes

const G4DNAMolecularReactionTable *& fMolReactionTable
 
G4VDNAReactionModelfpReactionModel
 
- Protected Attributes inherited from G4VITReactionProcess
const G4ITReactionTablefpReactionTable = nullptr
 
G4String fName
 

Detailed Description

G4DNAMolecularReaction is the reaction process used in G4DNAMolecularStepByStepModel. It test whether molecules can react after testing. If so, the reaction is made.

Deprecated:
This class will be removed

Definition at line 63 of file G4DNAMolecularReaction.hh.

Constructor & Destructor Documentation

◆ G4DNAMolecularReaction() [1/3]

G4DNAMolecularReaction::G4DNAMolecularReaction ( )

Definition at line 49 of file G4DNAMolecularReaction.cc.

51 , fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
52 , fpReactionModel(nullptr)
53{
54}
const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModel * fpReactionModel
const G4ITReactionTable * fpReactionTable
G4VITReactionProcess()=default

◆ G4DNAMolecularReaction() [2/3]

G4DNAMolecularReaction::G4DNAMolecularReaction ( G4VDNAReactionModel pReactionModel)
explicit

Definition at line 56 of file G4DNAMolecularReaction.cc.

58{
59 fpReactionModel = pReactionModel;
60}

◆ ~G4DNAMolecularReaction()

G4DNAMolecularReaction::~G4DNAMolecularReaction ( )
overridedefault

◆ G4DNAMolecularReaction() [3/3]

G4DNAMolecularReaction::G4DNAMolecularReaction ( const G4DNAMolecularReaction other)
delete

Member Function Documentation

◆ FindReaction()

std::vector< std::unique_ptr< G4ITReactionChange > > G4DNAMolecularReaction::FindReaction ( G4ITReactionSet pReactionSet,
const double  currentStepTime,
const double  ,
const bool  reachedUserStepTimeLimit 
)
overridevirtual

Implements G4VITReactionProcess.

Definition at line 130 of file G4DNAMolecularReaction.cc.

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
@ fStopAndKill
std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &) override
G4bool TestReactibility(const G4Track &, const G4Track &, double currentStepTime, bool userStepTimeLimit) override
void SelectThisReaction(G4ITReactionPtr reaction)
void CleanAllReaction()
G4ITReactionPerTrackMap & GetReactionMap()
G4TrackStatus GetTrackStatus() const

◆ MakeReaction()

std::unique_ptr< G4ITReactionChange > G4DNAMolecularReaction::MakeReaction ( const G4Track trackA,
const G4Track trackB 
)
overridevirtual

Implements G4VITReactionProcess.

Definition at line 84 of file G4DNAMolecularReaction.cc.

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}
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
@ fAlive
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
Data * GetReactionData(Reactant *, Reactant *) const
virtual void Push(G4Track *track)
static G4ITFinder * Instance()
virtual void Push(G4Track *)
static G4ITTrackHolder * Instance()
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const

Referenced by FindReaction().

◆ operator=()

G4DNAMolecularReaction & G4DNAMolecularReaction::operator= ( const G4DNAMolecularReaction other)
delete

◆ SetReactionModel()

void G4DNAMolecularReaction::SetReactionModel ( G4VDNAReactionModel pReactionModel)

Definition at line 125 of file G4DNAMolecularReaction.cc.

126{
127 fpReactionModel = pReactionModel;
128}

◆ TestReactibility()

G4bool G4DNAMolecularReaction::TestReactibility ( const G4Track trackA,
const G4Track trackB,
double  currentStepTime,
bool  userStepTimeLimit 
)
overridevirtual

Implements G4VITReactionProcess.

Definition at line 62 of file G4DNAMolecularReaction.cc.

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}
bool G4bool
Definition: G4Types.hh:86
virtual G4bool FindReaction(const G4Track &, const G4Track &, G4double, G4double &, G4bool)=0
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0

Referenced by FindReaction().

Member Data Documentation

◆ fMolReactionTable

const G4DNAMolecularReactionTable*& G4DNAMolecularReaction::fMolReactionTable
protected

Definition at line 84 of file G4DNAMolecularReaction.hh.

Referenced by MakeReaction().

◆ fpReactionModel

G4VDNAReactionModel* G4DNAMolecularReaction::fpReactionModel
protected

The documentation for this class was generated from the following files: