Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAScavengerMaterial.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//
27#include <memory>
29#include "G4StateManager.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4DNABoundingBox.hh"
35#include "G4VChemistryWorld.hh"
36#include "G4Scheduler.hh"
37#include "G4UnitsTable.hh"
38#include "G4MoleculeTable.hh"
39
40using namespace std;
41
42//------------------------------------------------------------------------------
43
45 G4VChemistryWorld* pChemistryInfo)
47 , fpChemistryInfo(pChemistryInfo)
48 , fIsInitialized(false)
49 , fCounterAgainstTime(false)
50 , fVerbose(0)
51{
52 Initialize();
53}
54
55//------------------------------------------------------------------------------
56
58{
59 if(fIsInitialized)
60 {
61 return;
62 }
63
64 if(fpChemistryInfo->size() == 0)
65 {
66 G4cout << "G4DNAScavengerMaterial existed but empty" << G4endl;
67 }
68 Reset();
69 fIsInitialized = true;
70}
71
73 MolType matConf) const
74{
75 // no change these molecules
76 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf)
77 {
78 G4ExceptionDescription exceptionDescription;
79 exceptionDescription << "matConf : "<<matConf->GetName();
80 G4Exception("G4DNAScavengerMaterial::GetNumberMoleculePerVolumeUnitForMaterialConf", "G4DNAScavengerMaterial001",
81 FatalErrorInArgument, exceptionDescription);
82 }
83
84 auto iter = fScavengerTable.find(matConf);
85 if(iter == fScavengerTable.end())
86 {
87 return 0;
88 }
89 else
90 {
91 if(iter->second >= 1)
92 {
93 return (floor)(iter->second);
94 }
95 else
96 {
97 return 0;
98 }
99 }
100}
101
103 MolType matConf, G4double time)
104{
105 // no change these molecules
106 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf ||
107 G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)") ==
108 matConf || // suppose that pH is not changed during simu
109 G4MoleculeTable::Instance()->GetConfiguration("OHm(B)") == matConf)
110 {
111 // G4cout<<"moletype : "<<matConf->GetName()<<G4endl;
112 // kobs is already counted these molecule concentrations
113 return;
114 }
115 if(!find(matConf)) // matConf must greater than 0
116 {
117 return;
118 }
119 fScavengerTable[matConf]--;
120 if(fScavengerTable[matConf] < 0) // projection
121 {
122 assert(false);
123 }
124
125 if(fCounterAgainstTime)
126 {
127 RemoveAMoleculeAtTime(matConf, time);
128 }
129}
130
132 MolType matConf, G4double time)
133{
134 // no change these molecules
135
136 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf ||
137 G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)") ==
138 matConf || // pH has no change
139 G4MoleculeTable::Instance()->GetConfiguration("OHm(B)") == matConf)
140 {
141 // G4cout<<"moletype : "<<matConf->GetName()<<G4endl;
142 // kobs is already counted these molecule concentrations
143 return;
144 }
145
146 auto it = fScavengerTable.find(matConf);
147 if(it == fScavengerTable.end()) // matConf must be in fScavengerTable
148 {
149 return;
150 }
151 fScavengerTable[matConf]++;
152
153 if(fCounterAgainstTime)
154 {
155 AddAMoleculeAtTime(matConf, time);
156 }
157}
158
160{
161 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
162 auto iter = fpChemistryInfo->begin();
163 G4cout << "**************************************************************"
164 << G4endl;
165 for(; iter != fpChemistryInfo->end(); iter++)
166 {
167 auto containedConf = iter->first;
168 // auto concentration = iter->second;
169 auto concentration =
170 fScavengerTable[containedConf] / (Avogadro * pConfinedBox->Volume());
171 G4cout << "Scavenger:" << containedConf->GetName() << " : "
172 << concentration / 1.0e-6 /*mm3 to L*/ << " (M) with : "
173 << fScavengerTable[containedConf] << " (molecules)"
174 << "in: " << pConfinedBox->Volume() / (um * um * um) << " (um3)"
175 << G4endl;
176 if(fScavengerTable[containedConf] < 1)
177 {
178 G4cout << "!!!!!!!!!!!!! this molecule has less one molecule for "
179 "considered volume"
180 << G4endl;
181 // assert(false);
182 }
183 if(fVerbose != 0)
184 {
185 Dump();
186 }
187 }
188 G4cout << "**************************************************************"
189 << G4endl;
190}
191
193{
194 if(fpChemistryInfo == nullptr)
195 {
196 return;
197 }
198
199 if(fpChemistryInfo->size() == 0)
200 {
201 return;
202 }
203
204 fScavengerTable.clear();
205 fCounterMap.clear();
206 fpLastSearch.reset(nullptr);
207
208 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
209 auto iter = fpChemistryInfo->begin();
210 for(; iter != fpChemistryInfo->end(); iter++)
211 {
212 auto containedConf = iter->first;
213 auto concentration = iter->second;
214 fScavengerTable[containedConf] =
215 floor(Avogadro * concentration * pConfinedBox->Volume());
216 fCounterMap[containedConf][1 * picosecond] =
217 floor(Avogadro * concentration * pConfinedBox->Volume());
218 }
219 if(fVerbose != 0){PrintInfo();}
220}
221
222//------------------------------------------------------------------------------
223
225 MolType molecule, G4double time, const G4ThreeVector* /*position*/,
226 int number)
227{
228 if(fVerbose != 0)
229 {
230 G4cout << "G4DNAScavengerMaterial::AddAMoleculeAtTime : "
231 << molecule->GetName() << " at time : " << G4BestUnit(time, "Time")
232 << G4endl;
233 }
234
235 auto counterMap_i = fCounterMap.find(molecule);
236
237 if(counterMap_i == fCounterMap.end())
238 {
239 fCounterMap[molecule][time] = number;
240 }
241 else if(counterMap_i->second.empty())
242 {
243 counterMap_i->second[time] = number;
244 }
245 else
246 {
247 auto end = counterMap_i->second.rbegin();
248
249 if(end->first <= time || fabs(end->first - time) <=
251 {
252 G4double newValue = end->second + number;
253 counterMap_i->second[time] = newValue;
254 if(newValue != (floor)(fScavengerTable[molecule])) // protection
255 {
256 G4String errMsg = "You are trying to add wrong molecule ";
257 G4Exception("AddAMoleculeAtTime", "",
258 FatalErrorInArgument, errMsg);
259
260 }
261 }
262 }
263}
264
265//------------------------------------------------------------------------------
266
268 MolType pMolecule, G4double time, const G4ThreeVector* /*position*/,
269 int number)
270{
271 NbMoleculeInTime& nbMolPerTime = fCounterMap[pMolecule];
272
273 if(fVerbose != 0)
274 {
275 auto it_ = nbMolPerTime.rbegin();
276 G4cout << "G4DNAScavengerMaterial::RemoveAMoleculeAtTime : "
277 << pMolecule->GetName() << " at time : " << G4BestUnit(time, "Time")
278
279 << " form : " << it_->second << G4endl;
280 }
281
282 if(nbMolPerTime.empty())
283 {
284 Dump();
285 G4String errMsg = "You are trying to remove molecule " +
286 pMolecule->GetName() +
287 " from the counter while this kind of molecules has not "
288 "been registered yet";
289 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "",
290 FatalErrorInArgument, errMsg);
291
292 return;
293 }
294 else
295 {
296 auto it = nbMolPerTime.rbegin();
297
298 if(it == nbMolPerTime.rend())
299 {
300 it--;
301
302 G4String errMsg = "There was no " + pMolecule->GetName() +
303 " recorded at the time or even before the time asked";
304 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "",
305 FatalErrorInArgument, errMsg);
306 }
307
308 G4double finalN = it->second - number;
309 if(finalN < 0)
310 {
311 Dump();
312
313 G4cout << "fScavengerTable : " << pMolecule->GetName() << " : "
314 << (fScavengerTable[pMolecule]) << G4endl;
315
317 errMsg << "After removal of " << number << " species of "
318 << " " << it->second << " " << pMolecule->GetName()
319 << " the final number at time " << G4BestUnit(time, "Time")
320 << " is less than zero and so not valid." << G4endl;
321 G4cout << " Global time is "
322 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
323 << ". Previous selected time is " << G4BestUnit(it->first, "Time")
324 << G4endl;
325 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "N_INF_0",
326 FatalException, errMsg);
327 }
328 nbMolPerTime[time] = finalN;
329
330 if(finalN != (floor)(fScavengerTable[pMolecule])) // protection
331 {
332 assert(false);
333 }
334 }
335}
336
338{
339 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
340 auto V = pConfinedBox->Volume();
341 for(const auto& it : fCounterMap)
342 {
343 auto pReactant = it.first;
344
345 G4cout << " --- > For " << pReactant->GetName() << G4endl;
346
347 for(const auto& it2 : it.second)
348 {
349 G4cout << " " << G4BestUnit(it2.first, "Time") << " "
350 << it2.second / (Avogadro * V * 1.0e-6 /*mm3 to L*/) << G4endl;
351 }
352 }
353}
354
356{
357 if(!fCounterAgainstTime)
358 {
359 G4cout << "fCounterAgainstTime == false" << G4endl;
360 assert(false);
361 }
362
363 G4bool sameTypeOfMolecule = SearchTimeMap(molecule);
364 auto output = SearchUpperBoundTime(time, sameTypeOfMolecule);
365 if(output < 0)
366 {
368 errMsg << "N molecules not valid < 0 : "<<
369 molecule->GetName() <<" N : "<< output << G4endl;
370 G4Exception("G4DNAScavengerMaterial::GetNMoleculesAtTime", "",
371 FatalErrorInArgument, errMsg);
372 }
373 return output;
374}
375
377{
378 if(fpLastSearch == nullptr)
379 {
380 fpLastSearch = std::make_unique<Search>();
381 }
382 else
383 {
384 if(fpLastSearch->fLowerBoundSet &&
385 fpLastSearch->fLastMoleculeSearched->first == molecule)
386 {
387 return true;
388 }
389 }
390
391 auto mol_it = fCounterMap.find(molecule);
392 fpLastSearch->fLastMoleculeSearched = mol_it;
393
394 if(mol_it != fCounterMap.end())
395 {
396 fpLastSearch->fLowerBoundTime =
397 fpLastSearch->fLastMoleculeSearched->second.end();
398 fpLastSearch->fLowerBoundSet = true;
399 }
400 else
401 {
402 fpLastSearch->fLowerBoundSet = false;
403 }
404
405 return false;
406}
407
408//------------------------------------------------------------------------------
409
411 G4bool sameTypeOfMolecule)
412{
413 auto mol_it = fpLastSearch->fLastMoleculeSearched;
414 if(mol_it == fCounterMap.end())
415 {
416 return 0;
417 }
418
419 NbMoleculeInTime& timeMap = mol_it->second;
420 if(timeMap.empty())
421 {
422 return 0;
423 }
424
425 if(sameTypeOfMolecule)
426 {
427 if(fpLastSearch->fLowerBoundSet &&
428 fpLastSearch->fLowerBoundTime != timeMap.end())
429 {
430 if(fpLastSearch->fLowerBoundTime->first < time)
431 {
432 auto upperToLast = fpLastSearch->fLowerBoundTime;
433 upperToLast++;
434
435 if(upperToLast == timeMap.end())
436 {
437 return fpLastSearch->fLowerBoundTime->second;
438 }
439
440 if(upperToLast->first > time)
441 {
442 return fpLastSearch->fLowerBoundTime->second;
443 }
444 }
445 }
446 }
447
448 auto up_time_it = timeMap.upper_bound(time);
449
450 if(up_time_it == timeMap.end())
451 {
452 auto last_time = timeMap.rbegin();
453 return last_time->second;
454 }
455 if(up_time_it == timeMap.begin())
456 {
457 return 0;
458 }
459
460 up_time_it--;
461
462 fpLastSearch->fLowerBoundTime = up_time_it;
463 fpLastSearch->fLowerBoundSet = true;
464
465 return fpLastSearch->fLowerBoundTime->second;
466}
@ FatalException
@ 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
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double Volume() const
std::map< G4double, int64_t, G4::MoleculeCounter::TimePrecision > NbMoleculeInTime
G4double GetNumberMoleculePerVolumeUnitForMaterialConf(MolType) const
MaterialMap::iterator end()
int64_t SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule)
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
void AddAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
int64_t GetNMoleculesAtTime(MolType molecule, G4double time)
void AddNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
G4DNAScavengerMaterial()=default
G4bool SearchTimeMap(MolType molecule)
void RemoveAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
const G4String & GetName() const
static G4MoleculeTable * Instance()
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
std::map< MolType, double >::iterator begin()
G4DNABoundingBox * GetChemistryBoundary() const
std::map< MolType, double >::iterator end()
static G4ThreadLocal double fPrecision