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