Geant4-11
Public Member Functions | Protected Attributes
G4DNAMakeReaction Class Reference

#include <G4DNAMakeReaction.hh>

Inheritance diagram for G4DNAMakeReaction:
G4VITReactionProcess

Public Member Functions

std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction (G4ITReactionSet *, const G4double, const G4double, const G4bool) override
 
 G4DNAMakeReaction ()
 
 G4DNAMakeReaction (const G4DNAMakeReaction &other)=delete
 
 G4DNAMakeReaction (G4VDNAReactionModel *)
 
virtual void Initialize ()
 
virtual G4bool IsApplicable (const G4ITType &, const G4ITType &) const
 
std::unique_ptr< G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &) override
 
G4DNAMakeReactionoperator= (const G4DNAMakeReaction &other)=delete
 
void SetReactionModel (G4VDNAReactionModel *)
 
virtual void SetReactionTable (const G4ITReactionTable *)
 
void SetTimeStepComputer (G4VITTimeStepComputer *)
 
G4bool TestReactibility (const G4Track &, const G4Track &, G4double currentStepTime, G4bool userStepTimeLimit) override
 
void UpdatePositionForReaction (G4Track &, G4Track &)
 
 ~G4DNAMakeReaction () override=default
 

Protected Attributes

const G4DNAMolecularReactionTable *& fMolReactionTable
 
G4String fName
 
G4VDNAReactionModelfpReactionModel
 
const G4ITReactionTablefpReactionTable = nullptr
 
G4VITTimeStepComputerfpTimeStepper
 
G4double fTimeStep
 

Detailed Description

Definition at line 37 of file G4DNAMakeReaction.hh.

Constructor & Destructor Documentation

◆ G4DNAMakeReaction() [1/3]

G4DNAMakeReaction::G4DNAMakeReaction ( )

Definition at line 41 of file G4DNAMakeReaction.cc.

43 , fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
44 , fpReactionModel(nullptr)
45 , fpTimeStepper(nullptr)
46 , fTimeStep(0)
47{
48}
G4VITTimeStepComputer * fpTimeStepper
G4VDNAReactionModel * fpReactionModel
const G4DNAMolecularReactionTable *& fMolReactionTable
const G4ITReactionTable * fpReactionTable
G4VITReactionProcess()=default

◆ G4DNAMakeReaction() [2/3]

G4DNAMakeReaction::G4DNAMakeReaction ( G4VDNAReactionModel pReactionModel)
explicit

Definition at line 50 of file G4DNAMakeReaction.cc.

52{
53 fpReactionModel = pReactionModel;
54}

References fpReactionModel.

◆ ~G4DNAMakeReaction()

G4DNAMakeReaction::~G4DNAMakeReaction ( )
overridedefault

◆ G4DNAMakeReaction() [3/3]

G4DNAMakeReaction::G4DNAMakeReaction ( const G4DNAMakeReaction other)
delete

Member Function Documentation

◆ FindReaction()

std::vector< std::unique_ptr< G4ITReactionChange > > G4DNAMakeReaction::FindReaction ( G4ITReactionSet pReactionSet,
const G4double  currentStepTime,
const  G4double,
const  G4bool 
)
overridevirtual

Implements G4VITReactionProcess.

Definition at line 185 of file G4DNAMakeReaction.cc.

189{
190 std::vector<std::unique_ptr<G4ITReactionChange>> ReactionInfo;
191 ReactionInfo.clear();
192 auto pReactionChange = dynamic_cast<G4DNAIndependentReactionTimeStepper*>(fpTimeStepper)->
193 FindReaction(pReactionSet,currentStepTime);
194
195 if (pReactionChange != nullptr)
196 {
197 ReactionInfo.push_back(std::move(pReactionChange));
198 }
199 return ReactionInfo;
200}
std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction(G4ITReactionSet *, const G4double, const G4double, const G4bool) override

References FindReaction(), and fpTimeStepper.

Referenced by FindReaction().

◆ Initialize()

virtual void G4VITReactionProcess::Initialize ( )
inlinevirtualinherited

First initialization (done once for all at the begin of the run) eg. check if the reaction table is given ...

Reimplemented in G4DNAIRT, and G4DNAIRT_geometries.

Definition at line 75 of file G4VITReactionProcess.hh.

75{;}

Referenced by G4ITModelProcessor::CalculateMinTimeStep().

◆ IsApplicable()

G4bool G4VITReactionProcess::IsApplicable ( const G4ITType ,
const G4ITType  
) const
virtualinherited

Definition at line 43 of file G4VITReactionProcess.cc.

44{
45 return true;
46}

◆ MakeReaction()

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

Implements G4VITReactionProcess.

Definition at line 71 of file G4DNAMakeReaction.cc.

73{
74 G4Track& tA = const_cast<G4Track&>(trackA);
75 G4Track& tB = const_cast<G4Track&>(trackB);
76 UpdatePositionForReaction( tA , tB );//TODO: should change it
77
78 std::unique_ptr<G4ITReactionChange> pChanges(new G4ITReactionChange());
79 pChanges->Initialize(trackA, trackB);
80
81 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
82 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
83
84 const auto pReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
85 const G4int nbProducts = pReactionData->GetNbProducts();
86 if (nbProducts)
87 {
88 const G4double D1 = pMoleculeA->GetDiffusionCoefficient();
89 const G4double D2 = pMoleculeB->GetDiffusionCoefficient();
90 const G4double sqrD1 = D1 == 0. ? 0. : std::sqrt(D1);
91 const G4double sqrD2 = D2 == 0. ? 0. : std::sqrt(D2);
92 const G4double inv_numerator = 1./(sqrD1 + sqrD2);
93 const G4ThreeVector reactionSite = sqrD2 * inv_numerator * tA.GetPosition()
94 + sqrD1 * inv_numerator * tB.GetPosition();
95
97 auto randP = (1-u) * tA.GetPosition() + u * tB.GetPosition();
98
99 for (G4int j = 0; j < nbProducts; ++j)
100 {
101 auto pProduct = new G4Molecule(pReactionData->GetProduct(j));
102 auto pProductTrack = pProduct->BuildTrack(trackA.GetGlobalTime(), (reactionSite + randP)/2);
103 pProductTrack->SetTrackStatus(fAlive);
104 G4ITTrackHolder::Instance()->Push(pProductTrack);
105 pChanges->AddSecondary(pProductTrack);
106 }
107 }
108 pChanges->KillParents(true);
109 return pChanges;
110}
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
@ fAlive
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
void UpdatePositionForReaction(G4Track &, G4Track &)
Data * GetReactionData(Reactant *, Reactant *) const
virtual void Push(G4Track *)
static G4ITTrackHolder * Instance()
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const

References fAlive, fMolReactionTable, G4UniformRand, G4Track::GetGlobalTime(), G4Molecule::GetMolecularConfiguration(), GetMolecule(), G4DNAMolecularReactionData::GetNbProducts(), G4Track::GetPosition(), G4DNAMolecularReactionTable::GetReactionData(), G4ITTrackHolder::Instance(), G4ITTrackHolder::Push(), and UpdatePositionForReaction().

◆ operator=()

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

◆ SetReactionModel()

void G4DNAMakeReaction::SetReactionModel ( G4VDNAReactionModel pReactionModel)

Definition at line 112 of file G4DNAMakeReaction.cc.

113{
114 fpReactionModel = pReactionModel;
115}

References fpReactionModel.

◆ SetReactionTable()

void G4VITReactionProcess::SetReactionTable ( const G4ITReactionTable pReactionTable)
virtualinherited

Definition at line 38 of file G4VITReactionProcess.cc.

39{
40 fpReactionTable = pReactionTable;
41}

References G4VITReactionProcess::fpReactionTable.

◆ SetTimeStepComputer()

void G4DNAMakeReaction::SetTimeStepComputer ( G4VITTimeStepComputer pStepper)

Definition at line 56 of file G4DNAMakeReaction.cc.

57{
58 fpTimeStepper = pStepper;
59}

References fpTimeStepper.

◆ TestReactibility()

G4bool G4DNAMakeReaction::TestReactibility ( const G4Track ,
const G4Track ,
G4double  currentStepTime,
G4bool  userStepTimeLimit 
)
overridevirtual

Implements G4VITReactionProcess.

Definition at line 61 of file G4DNAMakeReaction.cc.

65{
66 fTimeStep = currentStepTime;
67 return true;
68}

References fTimeStep.

◆ UpdatePositionForReaction()

void G4DNAMakeReaction::UpdatePositionForReaction ( G4Track trackA,
G4Track trackB 
)

Definition at line 117 of file G4DNAMakeReaction.cc.

119{
120 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
121 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
122 G4double D1 = pMoleculeA->GetDiffusionCoefficient();
123 G4double D2 = pMoleculeB->GetDiffusionCoefficient();
124
125 G4double reactionRadius = fpReactionModel->GetReactionRadius( pMoleculeA, pMoleculeB );
126 G4ThreeVector p1 = trackA.GetPosition();
127 G4ThreeVector p2 = trackB.GetPosition();
128
129 G4ThreeVector S1 = p1 - p2;
130 G4double distance = S1.mag();
131
132 if(D1 == 0)
133 {
134 trackB.SetPosition(p1);
135 return;
136 }
137 else if(D2 == 0)
138 {
139 trackA.SetPosition(p2);
140 return;
141 }
142
143 if(distance == 0)
144 {
145 G4ExceptionDescription exceptionDescription;
146 exceptionDescription << "Two particles are overlap: "
147 <<GetMolecule(trackA)->GetName()
148 <<" and "<<GetMolecule(trackB)->GetName()
149 <<" at "<<trackA.GetPosition();
150 G4Exception("G4DNAMakeReaction::PrepareForReaction()",
151 "G4DNAMakeReaction003",
152 FatalErrorInArgument,exceptionDescription);
153 }
154 S1.setMag(reactionRadius);
155
156 const G4double dt = fTimeStep;//irt - actualize molecule time
157
158 if(dt > 0)// irt > 0
159 {
160 G4double s12 = 2.0 * D1 * dt;
161 G4double s22 = 2.0 * D2 * dt;
162 G4double sigma = s12 + ( s12 * s12 ) / s22;
163 G4double alpha = reactionRadius * distance / (2 * (D1 + D2) * dt );
164
165 G4ThreeVector S2 = (p1 + ( s12 / s22 ) * p2) +
167 G4RandGauss::shoot(0.0, sigma),
168 G4RandGauss::shoot(0.0, sigma));
169
170 S1.setPhi(rad * G4UniformRand() * 2.0 * CLHEP::pi);
171
172 S1.setTheta(rad * std::acos( 1.0 + (1. / alpha) *
173 std::log(1.0 - G4UniformRand() *
174 (1.-std::exp(-2.0 * alpha)))));
175
176 const G4ThreeVector R1 = (D1 * S1 + D2 * S2) / (D1 + D2);
177 const G4ThreeVector R2 = D2 * (S2 - S1) / (D1 + D2);
178
179 trackA.SetPosition(R1);
180 trackB.SetPosition(R2);
181 }
182}
@ 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
static const G4double alpha
static constexpr double rad
Definition: G4SIunits.hh:129
CLHEP::Hep3Vector G4ThreeVector
void setTheta(double)
double mag() const
void setMag(double)
Definition: ThreeVector.cc:20
void setPhi(double)
const G4String & GetName() const
Definition: G4Molecule.cc:338
void SetPosition(const G4ThreeVector &aValue)
virtual G4double GetReactionRadius(const G4MolecularConfiguration *, const G4MolecularConfiguration *)=0
static constexpr double pi
Definition: SystemOfUnits.h:55
ThreeVector shoot(const G4int Ap, const G4int Af)

References alpha, FatalErrorInArgument, fpReactionModel, fTimeStep, G4Exception(), G4UniformRand, G4MolecularConfiguration::GetDiffusionCoefficient(), G4Molecule::GetMolecularConfiguration(), GetMolecule(), G4Molecule::GetName(), G4Track::GetPosition(), G4VDNAReactionModel::GetReactionRadius(), CLHEP::Hep3Vector::mag(), CLHEP::pi, rad, CLHEP::Hep3Vector::setMag(), CLHEP::Hep3Vector::setPhi(), G4Track::SetPosition(), CLHEP::Hep3Vector::setTheta(), and G4INCL::DeJongSpin::shoot().

Referenced by MakeReaction().

Field Documentation

◆ fMolReactionTable

const G4DNAMolecularReactionTable*& G4DNAMakeReaction::fMolReactionTable
protected

Definition at line 58 of file G4DNAMakeReaction.hh.

Referenced by MakeReaction().

◆ fName

G4String G4VITReactionProcess::fName
protectedinherited

Definition at line 91 of file G4VITReactionProcess.hh.

◆ fpReactionModel

G4VDNAReactionModel* G4DNAMakeReaction::fpReactionModel
protected

◆ fpReactionTable

const G4ITReactionTable* G4VITReactionProcess::fpReactionTable = nullptr
protectedinherited

Definition at line 90 of file G4VITReactionProcess.hh.

Referenced by G4VITReactionProcess::SetReactionTable().

◆ fpTimeStepper

G4VITTimeStepComputer* G4DNAMakeReaction::fpTimeStepper
protected

Definition at line 60 of file G4DNAMakeReaction.hh.

Referenced by FindReaction(), and SetTimeStepComputer().

◆ fTimeStep

G4double G4DNAMakeReaction::fTimeStep
protected

Definition at line 61 of file G4DNAMakeReaction.hh.

Referenced by TestReactibility(), and UpdatePositionForReaction().


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