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