Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAUpdateSystemModel.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//
28#include "G4Molecule.hh"
30#include "G4UnitsTable.hh"
31#include "G4MoleculeCounter.hh"
33#include "G4Scheduler.hh"
36{}
37
38void G4DNAUpdateSystemModel::SetMesh(G4DNAMesh* pMesh) { fpMesh = pMesh; }
39void G4DNAUpdateSystemModel::KillMolecule(const Index& index, MolType type)
40{
41 // kill normal molecule
42 auto& node = fpMesh->GetVoxelMapList(index);
43 auto iter = node.find(type);
44 if(iter != node.end())
45 {
46 if(iter->second <= 0)
47 {
48 G4ExceptionDescription exceptionDescription;
49 exceptionDescription
50 << "G4DNAUpdateSystemModel::KillMolecule::molecule : "
51 << type->GetName() << " index : " << index
52 << " number : " << iter->second << G4endl;
53 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler002",
54 FatalErrorInArgument, exceptionDescription);
55 }
56 iter->second--;
57 if(G4VMoleculeCounter::Instance()->InUse())
58 {
60 }
61 }
62 else
63 {
64 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
66 if(pScavengerMaterial != nullptr)
67 {
69 type, fGlobalTime);
70 }
71 else
72 {
73 G4ExceptionDescription exceptionDescription;
74 exceptionDescription
75 << "index : " << index << " " << type->GetName()
76 << " This molecule is not belong scavengers or particle-base"
77 << G4endl;
78 G4Exception("G4DNAEventScheduler::Stepping", "G4DNAEventScheduler002",
79 FatalErrorInArgument, exceptionDescription);
80 }
81 }
82}
83
84void G4DNAUpdateSystemModel::JumpTo(const Index& index, MolType type)
85{
86 auto& node = fpMesh->GetVoxelMapList(index);
87 auto iter = node.find(type);
88 if(iter != node.end())
89 {
90 if(iter->second <= 0)
91 {
92 G4ExceptionDescription exceptionDescription;
93 exceptionDescription << "G4DNAUpdateSystemModel::JumpTo::molecule : "
94 << type->GetName() << " index : " << index
95 << " number : " << iter->second;
96 G4Exception("G4DNAUpdateSystemModel::JumpTo", "G4DNAUpdateSystemModel001",
97 FatalErrorInArgument, exceptionDescription);
98 }
99 iter->second--;
100 }
101 else
102 {
103 fpMesh->PrintVoxel(index);
104 G4ExceptionDescription exceptionDescription;
105 exceptionDescription << "index : " << index << " " << type->GetName()
106 << " There is no this type";
107 G4Exception("G4DNAUpdateSystemModel::JumpTo", "G4DNAUpdateSystemModel002",
108 FatalErrorInArgument, exceptionDescription);
109 }
110}
111
112void G4DNAUpdateSystemModel::CreateMolecule(const Index& index, MolType type)
113{
114 // for scavenger
115 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
117 if(pScavengerMaterial != nullptr && pScavengerMaterial->find(type))
118 {
120 type, fGlobalTime);
121 return;
122 }
123 // for molecule
124 auto& node = fpMesh->GetVoxelMapList(index);
125 auto iter = node.find(type);
126 if(iter != node.end())
127 {
128 iter->second++;
129 }
130 else
131 {
132 node[type] = 1;
133 }
134
136 {
138 }
139}
140
141void G4DNAUpdateSystemModel::JumpIn(const Index& index, MolType type)
142{
143 // for molecule
144 auto& node = fpMesh->GetVoxelMapList(index);
145 auto iter = node.find(type);
146 if(iter != node.end())
147 {
148 iter->second++;
149 }
150 else
151 {
152 node[type] = 1;
153 }
154}
155
157 const ReactionData& data)
158{
159 auto reactant1 = data.GetReactant1();
160 auto reactant2 = data.GetReactant2();
161#ifdef G4VERBOSE
162 if(fVerbose != 0)
163 {
164 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
165 << " Reaction : " << reactant1->GetName() << " + "
166 << reactant2->GetName() << " -> ";
167 }
168#endif
169 const G4int nbProducts = data.GetNbProducts();
170 if(nbProducts != 0)
171 {
172 for(G4int j = 0; j < nbProducts; ++j)
173 {
174#ifdef G4VERBOSE
175 if((fVerbose != 0) && j != 0)
176 {
177 G4cout << " + ";
178 }
179 if(fVerbose != 0)
180 {
181 G4cout << data.GetProduct(j)->GetName();
182 }
183#endif
184 CreateMolecule(index, data.GetProduct(j));
185 }
186 }
187 else
188 {
189#ifdef G4VERBOSE
190 if(fVerbose != 0)
191 {
192 G4cout << "No product";
193 }
194#endif
195 }
196#ifdef G4VERBOSE
197 if(fVerbose != 0)
198 {
199 G4cout << G4endl;
200 }
201#endif
202 KillMolecule(index, reactant1);
203 KillMolecule(index, reactant2);
204}
205
207 const JumpingData& data)
208{
209 auto reactant = std::get<0>(data);
210 auto JunpToIndex = std::get<1>(data);
211#ifdef G4VERBOSE
212 if(fVerbose > 1)
213 {
214 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
215 << " Jumping : " << reactant->GetName() << " from " << index
216 << " -> " << JunpToIndex << G4endl;
217 }
218#endif
219 JumpTo(index, reactant);
220 JumpIn(JunpToIndex, reactant);
221}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4BestUnit(a, b)
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void PrintVoxel(const Index &index)
Definition: G4DNAMesh.cc:102
Data & GetVoxelMapList(const Index &index)
Definition: G4DNAMesh.cc:62
Reactant * GetProduct(G4int i) const
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
void AddNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
std::pair< MolType, Index > JumpingData
void UpdateSystem(const Index &index, const ReactionData &data)
const G4String & GetName() const
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194
virtual void RemoveAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1)=0
static G4VMoleculeCounter * Instance()
virtual void AddAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1)=0