Geant4 11.4.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 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 50 of file G4DNAMolecularReaction.cc.

51 :
53 , fpReactionModel(nullptr)
54{
55}
ReturnType & reference_cast(OriginalType &source)
const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModel * fpReactionModel
const G4ITReactionTable * fpReactionTable

Referenced by G4DNAMolecularReaction(), G4DNAMolecularReaction(), and operator=().

◆ G4DNAMolecularReaction() [2/3]

G4DNAMolecularReaction::G4DNAMolecularReaction ( G4VDNAReactionModel * pReactionModel)
explicit

Definition at line 57 of file G4DNAMolecularReaction.cc.

59{
60 fpReactionModel = pReactionModel;
61}

◆ ~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 134 of file G4DNAMolecularReaction.cc.

139{
140 std::vector<std::unique_ptr<G4ITReactionChange>> fReactionInfo;
141 fReactionInfo.clear();
142
143 if (pReactionSet == nullptr)
144 {
145 return fReactionInfo;
146 }
147
148 G4ITReactionPerTrackMap& reactionPerTrackMap = pReactionSet->GetReactionMap();
149 for (auto tracks_i = reactionPerTrackMap.cbegin();
150 tracks_i != reactionPerTrackMap.cend();
151 tracks_i = reactionPerTrackMap.cbegin())
152 {
153 G4Track* pTrackA = tracks_i->first;
154 if (pTrackA->GetTrackStatus() == fStopAndKill)
155 {
156 continue;
157 }
158
159 G4ITReactionPerTrackPtr reactionPerTrack = tracks_i->second;
160 G4ITReactionList& reactionList = reactionPerTrack->GetReactionList();
161
162 assert(reactionList.cbegin() != reactionList.cend());
163
164 for (auto it = reactionList.cbegin(); it != reactionList.cend(); it = reactionList.cbegin())
165 {
167 G4Track* pTrackB = reaction->GetReactant(pTrackA);
168 if (pTrackB->GetTrackStatus() == fStopAndKill)
169 {
170 continue;
171 }
172
173 if (pTrackB == pTrackA)
174 {
175 G4ExceptionDescription exceptionDescription;
176 exceptionDescription
177 << "The IT reaction process sent back a reaction between trackA and trackB. ";
178 exceptionDescription << "The problem is trackA == trackB";
179 G4Exception("G4ITModelProcessor::FindReaction",
180 "ITModelProcessor005",
182 exceptionDescription);
183 }
184
185 pReactionSet->SelectThisReaction(std::move(reaction));
186
187 if (TestReactibility(*pTrackA, *pTrackB, currentStepTime, reachedUserStepTimeLimit))
188 {
189 auto pReactionChange = MakeReaction(*pTrackA, *pTrackB);
190
191 if (pReactionChange)
192 {
193 fReactionInfo.push_back(std::move(pReactionChange));
194 break;
195 }
196 }
197 }
198 }
199
200 pReactionSet->CleanAllReaction();
201 return fReactionInfo;
202}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
std::map< G4Track *, G4ITReactionPerTrackPtr, compTrackPerID > G4ITReactionPerTrackMap
G4shared_ptr< G4ITReaction > G4ITReactionPtr
std::list< G4ITReactionPtr > G4ITReactionList
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
@ 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)
G4ITReactionPerTrackMap & GetReactionMap()
G4TrackStatus GetTrackStatus() const

◆ MakeReaction()

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

Implements G4VITReactionProcess.

Definition at line 85 of file G4DNAMolecularReaction.cc.

87{
88 std::unique_ptr<G4ITReactionChange> pChanges(new G4ITReactionChange());
89 pChanges->Initialize(trackA, trackB);
90
91 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
92 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
93
94 const auto pReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
95 // Notify molecule (reaction) counter
96 if (G4MoleculeCounterManager::Instance()->GetIsActive()) {
98 }
99 const G4int nbProducts = pReactionData->GetNbProducts();
100
101 if (nbProducts != 0)
102 {
103 const G4double D1 = pMoleculeA->GetDiffusionCoefficient();
104 const G4double D2 = pMoleculeB->GetDiffusionCoefficient();
105 const G4double sqrD1 = D1 == 0. ? 0. : std::sqrt(D1);
106 const G4double sqrD2 = D2 == 0. ? 0. : std::sqrt(D2);
107 const G4double inv_numerator = 1./(sqrD1 + sqrD2);
108 const G4ThreeVector reactionSite = sqrD2 * inv_numerator * trackA.GetPosition()
109 + sqrD1 * inv_numerator * trackB.GetPosition();
110
111 for (G4int j = 0; j < nbProducts; ++j)
112 {
113 auto pProduct = new G4Molecule(pReactionData->GetProduct(j));
114 auto pProductTrack = pProduct->BuildTrack(trackA.GetGlobalTime(), reactionSite);
115
116 pProductTrack->SetTrackStatus(fAlive);
117
118 G4ITTrackHolder::Instance()->Push(pProductTrack);
119
120 pChanges->AddSecondary(pProductTrack);
121 G4MoleculeFinder::Instance()->Push(pProductTrack);
122 }
123 }
124
125 pChanges->KillParents(true);
126 return pChanges;
127}
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:75
CLHEP::Hep3Vector G4ThreeVector
@ fAlive
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void Push(G4Track *track) override
static G4ITFinder * Instance()
void Push(G4Track *) override
static G4ITTrackHolder * Instance()
void RecordReaction(const G4DNAMolecularReactionData *, G4double, G4int=1)
static G4MoleculeCounterManager * Instance()
const G4MolecularConfiguration * GetMolecularConfiguration() const
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 129 of file G4DNAMolecularReaction.cc.

130{
131 fpReactionModel = pReactionModel;
132}

◆ TestReactibility()

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

Implements G4VITReactionProcess.

Definition at line 63 of file G4DNAMolecularReaction.cc.

67{
68 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
69 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
70
71 const G4double reactionRadius = fpReactionModel->GetReactionRadius(pMoleculeA, pMoleculeB);
72
73 G4double separationDistance = -1.;
74
75 if (currentStepTime == 0.)
76 {
77 userStepTimeLimit = false;
78 }
79
80 G4bool output = fpReactionModel->FindReaction(trackA, trackB, reactionRadius,
81 separationDistance, userStepTimeLimit);
82 return output;
83}
bool G4bool
Definition G4Types.hh:86

Referenced by FindReaction().

Member Data Documentation

◆ fMolReactionTable

const G4DNAMolecularReactionTable*& G4DNAMolecularReaction::fMolReactionTable
protected

Definition at line 84 of file G4DNAMolecularReaction.hh.

Referenced by G4DNAMolecularReaction(), and MakeReaction().

◆ fpReactionModel

G4VDNAReactionModel* G4DNAMolecularReaction::fpReactionModel
protected

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