Geant4-11
Public Types | Public Member Functions | Private Attributes
G4DNAUpdateSystemModel Class Reference

#include <G4DNAUpdateSystemModel.hh>

Inheritance diagram for G4DNAUpdateSystemModel:
G4VUpdateSystemModel

Public Types

using Index = G4Voxel::Index
 
using JumpingData = std::pair< MolType, Index >
 
using MolType = const G4MolecularConfiguration *
 
using ReactionData = const G4DNAMolecularReactionData
 

Public Member Functions

void CreateMolecule (const Index &index, MolType)
 
 G4DNAUpdateSystemModel ()
 
void JumpIn (const Index &index, MolType)
 
void JumpTo (const Index &index, MolType type)
 
void KillMolecule (const Index &index, MolType type)
 
void SetGlobalTime (const G4double &globalTime)
 
void SetMesh (G4DNAMesh *)
 
void SetVerbose (G4int verbose)
 
void UpdateSystem (const Index &index, const JumpingData &data)
 
void UpdateSystem (const Index &index, const ReactionData &data)
 
 ~G4DNAUpdateSystemModel () override=default
 

Private Attributes

G4double fGlobalTime
 
G4DNAMeshfpMesh
 
G4int fVerbose
 

Detailed Description

Definition at line 36 of file G4DNAUpdateSystemModel.hh.

Member Typedef Documentation

◆ Index

Definition at line 39 of file G4DNAUpdateSystemModel.hh.

◆ JumpingData

Definition at line 41 of file G4DNAUpdateSystemModel.hh.

◆ MolType

Definition at line 40 of file G4DNAUpdateSystemModel.hh.

◆ ReactionData

Definition at line 42 of file G4DNAUpdateSystemModel.hh.

Constructor & Destructor Documentation

◆ G4DNAUpdateSystemModel()

G4DNAUpdateSystemModel::G4DNAUpdateSystemModel ( )

Definition at line 34 of file G4DNAUpdateSystemModel.cc.

◆ ~G4DNAUpdateSystemModel()

G4DNAUpdateSystemModel::~G4DNAUpdateSystemModel ( )
overridedefault

Member Function Documentation

◆ CreateMolecule()

void G4DNAUpdateSystemModel::CreateMolecule ( const Index index,
MolType  type 
)

Definition at line 105 of file G4DNAUpdateSystemModel.cc.

106{
107 // for scavenger
108 auto pScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
110 if(pScavengerMaterial != nullptr && pScavengerMaterial->find(type))
111 {
113 type, fGlobalTime);
114 return;
115 }
116 // for molecule
117 auto& node = fpMesh->GetVoxelMapList(index);
118 auto iter = node.find(type);
119 if(iter != node.end())
120 {
121 iter->second++;
122 }
123 else
124 {
125 node[type] = 1;
126 }
127
129 {
131 }
132}
G4Voxel::MapList & GetVoxelMapList(Key key)
Definition: G4DNAMesh.cc:37
void AddNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194
static G4VMoleculeCounter * Instance()
virtual void AddAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1)=0

References G4VMoleculeCounter::AddAMoleculeAtTime(), G4DNAScavengerMaterial::AddNumberMoleculePerVolumeUnitForMaterialConf(), fGlobalTime, fpMesh, G4Scheduler::GetScavengerMaterial(), G4DNAMesh::GetVoxelMapList(), G4Scheduler::Instance(), G4VMoleculeCounter::Instance(), and G4VMoleculeCounter::InUse().

Referenced by UpdateSystem().

◆ JumpIn()

void G4DNAUpdateSystemModel::JumpIn ( const Index index,
MolType  type 
)

Definition at line 134 of file G4DNAUpdateSystemModel.cc.

135{
136 // for molecule
137 auto& node = fpMesh->GetVoxelMapList(index);
138 auto iter = node.find(type);
139 if(iter != node.end())
140 {
141 iter->second++;
142 }
143 else
144 {
145 node[type] = 1;
146 }
147}

References fpMesh, and G4DNAMesh::GetVoxelMapList().

Referenced by UpdateSystem().

◆ JumpTo()

void G4DNAUpdateSystemModel::JumpTo ( const Index index,
MolType  type 
)

Definition at line 81 of file G4DNAUpdateSystemModel.cc.

82{
83 auto& node = fpMesh->GetVoxelMapList(index);
84
85 auto iter = node.find(type);
86 if(iter != node.end())
87 {
88 if(iter->second <= 0)
89 {
90 G4cout << "G4DNAUpdateSystemModel::KillMolecule::molecule : "
91 << type->GetName() << " index : " << index
92 << " number : " << iter->second << G4endl;
93 assert(false);
94 }
95 iter->second--;
96 }
97 else
98 {
99 G4cout << "index : " << index << " " << type->GetName() << G4endl;
100 G4cout << "This molecule is not belong particle-base" << G4endl;
101 assert(false);
102 }
103}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

References fpMesh, G4cout, G4endl, G4MolecularConfiguration::GetName(), and G4DNAMesh::GetVoxelMapList().

Referenced by UpdateSystem().

◆ KillMolecule()

void G4DNAUpdateSystemModel::KillMolecule ( const Index index,
MolType  type 
)

Definition at line 42 of file G4DNAUpdateSystemModel.cc.

43{
44 // kill normal molecule
45 auto& node = fpMesh->GetVoxelMapList(index);
46 auto iter = node.find(type);
47 if(iter != node.end())
48 {
49 if(iter->second <= 0)
50 {
51 G4cout << "G4DNAUpdateSystemModel::KillMolecule::molecule : "
52 << type->GetName() << " index : " << index
53 << " number : " << iter->second << G4endl;
54 assert(false);
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 G4cout << "index : " << index << " " << type->GetName() << G4endl;
74 G4cout << "This molecule is not belong scavengers or particle-base"
75 << G4endl;
76 assert(false);
77 }
78 }
79}
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
virtual void RemoveAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1)=0

References fGlobalTime, fpMesh, G4cout, G4endl, G4MolecularConfiguration::GetName(), G4Scheduler::GetScavengerMaterial(), G4DNAMesh::GetVoxelMapList(), G4Scheduler::Instance(), G4VMoleculeCounter::Instance(), G4DNAScavengerMaterial::ReduceNumberMoleculePerVolumeUnitForMaterialConf(), and G4VMoleculeCounter::RemoveAMoleculeAtTime().

Referenced by UpdateSystem().

◆ SetGlobalTime()

void G4DNAUpdateSystemModel::SetGlobalTime ( const G4double globalTime)
inline

Definition at line 56 of file G4DNAUpdateSystemModel.hh.

56{ fGlobalTime = globalTime; }

References fGlobalTime.

◆ SetMesh()

void G4DNAUpdateSystemModel::SetMesh ( G4DNAMesh pMesh)

Definition at line 41 of file G4DNAUpdateSystemModel.cc.

41{ fpMesh = pMesh; }

References fpMesh.

◆ SetVerbose()

void G4DNAUpdateSystemModel::SetVerbose ( G4int  verbose)
inline

Definition at line 57 of file G4DNAUpdateSystemModel.hh.

57{ fVerbose = verbose; }

References fVerbose.

◆ UpdateSystem() [1/2]

void G4DNAUpdateSystemModel::UpdateSystem ( const Index index,
const JumpingData data 
)

Definition at line 260 of file G4DNAUpdateSystemModel.cc.

262{
263 auto reactant = std::get<0>(data);
264 auto JunpToIndex = std::get<1>(data);
265#ifdef G4VERBOSE
266 if(fVerbose > 1)
267 {
268 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
269 << " Jumping : " << reactant->GetName() << " from " << index
270 << " -> " << JunpToIndex << G4endl;
271 }
272#endif
273 JumpTo(index, reactant);
274 JumpIn(JunpToIndex, reactant);
275}
#define G4BestUnit(a, b)
void JumpIn(const Index &index, MolType)
void JumpTo(const Index &index, MolType type)

References fGlobalTime, fVerbose, G4BestUnit, G4cout, G4endl, JumpIn(), and JumpTo().

◆ UpdateSystem() [2/2]

void G4DNAUpdateSystemModel::UpdateSystem ( const Index index,
const ReactionData data 
)

Definition at line 149 of file G4DNAUpdateSystemModel.cc.

151{
152 auto reactant1 = data.GetReactant1();
153 auto reactant2 = data.GetReactant2();
154#ifdef G4VERBOSE
155 if(fVerbose != 0)
156 {
157 G4cout << "At time : " << std::setw(7) << G4BestUnit(fGlobalTime, "Time")
158 << " Reaction : " << reactant1->GetName() << " + "
159 << reactant2->GetName() << " -> ";
160 }
161#endif
162 const G4int nbProducts = data.GetNbProducts();
163 if(nbProducts != 0)
164 {
165 for(size_t j = 0; j < (size_t) nbProducts; ++j)
166 {
167#ifdef G4VERBOSE
168 if((fVerbose != 0) && j != 0)
169 {
170 G4cout << " + ";
171 }
172 if(fVerbose != 0)
173 {
174 G4cout << data.GetProduct(j)->GetName();
175 // for test
176 // G4cout<<" fGlobalTime : "<<fGlobalTime;
177 // end fortest
178 }
179#endif
180 CreateMolecule(index, data.GetProduct(j));
181 //#define DEBUG 1
182
183#ifdef DEBUG
184 if(G4MoleculeCounter::Instance()->InUse())
185 if(fpMesh->GetNumberOfType(data.GetProduct(j)) !=
186 G4MoleculeCounter::Instance()->GetCurrentNumberOf(
187 data.GetProduct(j)))
188 {
189 G4cout << "*********G4DNAUpdateSystemModel::DEBUG::GetNumberOfType("
190 << data.GetProduct(j)->GetName()
191 << ") : " << fpMesh->GetNumberOfType(data.GetProduct(j))
192 << G4endl;
193 G4cout << "G4MoleculeCounter::GetCurrentNumberOf ("
194 << data.GetProduct(j)->GetName() << ") : "
195 << G4MoleculeCounter::Instance()->GetCurrentNumberOf(
196 data.GetProduct(j))
197 << G4endl;
199 throw;
200 }
201
202#endif
203 }
204 }
205 else
206 {
207#ifdef G4VERBOSE
208 if(fVerbose != 0)
209 {
210 G4cout << "No product";
211 // for test
212 // G4cout<<" fGlobalTime : "<<fGlobalTime;
213 // end fortest
214 }
215#endif
216 }
217#ifdef G4VERBOSE
218 if(fVerbose != 0)
219 {
220 G4cout << G4endl;
221 }
222#endif
223 KillMolecule(index, reactant1);
224#ifdef DEBUG
225 if(G4MoleculeCounter::Instance()->InUse())
226 if(fpMesh->GetNumberOfType(reactant1) !=
227 G4MoleculeCounter::Instance()->GetCurrentNumberOf(reactant1))
228 {
229 G4cout << "*********G4DNAUpdateSystemModel::DEBUG::GetNumberOfType("
230 << reactant1->GetName()
231 << ") : " << fpMesh->GetNumberOfType(reactant1) << G4endl;
232 G4cout << "G4MoleculeCounter::GetCurrentNumberOf ("
233 << reactant1->GetName() << ") : "
234 << G4MoleculeCounter::Instance()->GetCurrentNumberOf(reactant1)
235 << G4endl;
237 throw;
238 }
239#endif
240 KillMolecule(index, reactant2);
241#ifdef DEBUG
242
243 if(G4MoleculeCounter::Instance()->InUse())
244 if(fpMesh->GetNumberOfType(reactant2) !=
245 G4MoleculeCounter::Instance()->GetCurrentNumberOf(reactant2))
246 {
247 G4cout << "*********G4DNAUpdateSystemModel::DEBUG::GetNumberOfType("
248 << reactant2->GetName()
249 << ") : " << fpMesh->GetNumberOfType(reactant2) << G4endl;
250 G4cout << "G4MoleculeCounter::GetCurrentNumberOf ("
251 << reactant2->GetName() << ") : "
252 << G4MoleculeCounter::Instance()->GetCurrentNumberOf(reactant2)
253 << G4endl;
255 throw;
256 }
257#endif
258}
int G4int
Definition: G4Types.hh:85
G4int GetNumberOfType(G4Voxel::MolType type) const
Definition: G4DNAMesh.cc:78
void CreateMolecule(const Index &index, MolType)
void KillMolecule(const Index &index, MolType type)
static G4MoleculeCounter * Instance()

References CreateMolecule(), G4MoleculeCounter::Dump(), fGlobalTime, fpMesh, fVerbose, G4BestUnit, G4cout, G4endl, G4MolecularConfiguration::GetName(), G4DNAMolecularReactionData::GetNbProducts(), G4DNAMesh::GetNumberOfType(), G4DNAMolecularReactionData::GetProduct(), G4DNAMolecularReactionData::GetReactant1(), G4DNAMolecularReactionData::GetReactant2(), G4MoleculeCounter::Instance(), and KillMolecule().

Field Documentation

◆ fGlobalTime

G4double G4DNAUpdateSystemModel::fGlobalTime
private

◆ fpMesh

G4DNAMesh* G4DNAUpdateSystemModel::fpMesh
private

◆ fVerbose

G4int G4DNAUpdateSystemModel::fVerbose
private

Definition at line 61 of file G4DNAUpdateSystemModel.hh.

Referenced by SetVerbose(), and UpdateSystem().


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