Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MoleculeCounter.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
28#include "G4MoleculeCounter.hh"
29#include "G4MoleculeTable.hh"
30#include "G4UIcommand.hh"
31#include "G4UnitsTable.hh"
34#include "G4SystemOfUnits.hh"
35#include "G4Scheduler.hh" // TODO: remove this dependency
36#include <iomanip>
37
38using namespace std;
39
40namespace G4{
41namespace MoleculeCounter {
42
43bool TimePrecision::operator()(const double& a, const double& b) const
44{
45 if (std::fabs(a - b) < fPrecision)
46 {
47 return false;
48 }
49 else
50 {
51 return a < b;
52 }
53}
54
55G4ThreadLocal double TimePrecision::fPrecision = 0.5 * picosecond;
56}
57}
58
59//------------------------------------------------------------------------------
61{
62 if (!fpInstance)
63 {
65 }
66 return dynamic_cast<G4MoleculeCounter*>(fpInstance);
67}
68
69//------------------------------------------------------------------------------
70
72{
73 fVerbose = 0;
75}
76
77//------------------------------------------------------------------------------
78
80
81//------------------------------------------------------------------------------
82
84{
86 while ((mol_iterator)())
87 {
88 if (IsRegistered(mol_iterator.value()->GetDefinition()) == false)
89 {
90 continue;
91 }
92
93 fCounterMap[mol_iterator.value()]; // initialize the second map
94 }
95}
96
97//------------------------------------------------------------------------------
98
99void G4MoleculeCounter::SetTimeSlice(double timeSlice)
100{
102}
103
104//------------------------------------------------------------------------------
105
107{
108 if (fpLastSearch == nullptr)
109 {
110 fpLastSearch.reset(new Search());
111 }
112 else
113 {
114 if (fpLastSearch->fLowerBoundSet &&
115 fpLastSearch->fLastMoleculeSearched->first == molecule)
116 {
117 return true;
118 }
119 }
120
121 auto mol_it = fCounterMap.find(molecule);
122 fpLastSearch->fLastMoleculeSearched = mol_it;
123
124 if (mol_it != fCounterMap.end())
125 {
126 fpLastSearch->fLowerBoundTime = fpLastSearch->fLastMoleculeSearched->second
127 .end();
128 fpLastSearch->fLowerBoundSet = true;
129 }
130 else
131 {
132 fpLastSearch->fLowerBoundSet = false;
133 }
134
135 return false;
136}
137
138//------------------------------------------------------------------------------
139
141 bool sameTypeOfMolecule)
142{
143 auto mol_it = fpLastSearch->fLastMoleculeSearched;
144 if (mol_it == fCounterMap.end())
145 {
146 return 0;
147 }
148
149 NbMoleculeAgainstTime& timeMap = mol_it->second;
150 if (timeMap.empty())
151 {
152 return 0;
153 }
154
155 if (sameTypeOfMolecule == true)
156 {
157 if (fpLastSearch->fLowerBoundSet && fpLastSearch->fLowerBoundTime != timeMap.end())
158 {
159 if (fpLastSearch->fLowerBoundTime->first < time)
160 {
161 auto upperToLast = fpLastSearch->fLowerBoundTime;
162 upperToLast++;
163
164 if (upperToLast == timeMap.end())
165 {
166 return fpLastSearch->fLowerBoundTime->second;
167 }
168
169 if (upperToLast->first > time)
170 {
171 return fpLastSearch->fLowerBoundTime->second;
172 }
173 }
174 }
175 }
176
177 auto up_time_it = timeMap.upper_bound(time);
178
179 if (up_time_it == timeMap.end())
180 {
181 NbMoleculeAgainstTime::reverse_iterator last_time = timeMap.rbegin();
182 return last_time->second;
183 }
184 if (up_time_it == timeMap.begin())
185 {
186 return 0;
187 }
188
189 up_time_it--;
190
191 fpLastSearch->fLowerBoundTime = up_time_it;
192 fpLastSearch->fLowerBoundSet = true;
193
194 return fpLastSearch->fLowerBoundTime->second;
195}
196
197//------------------------------------------------------------------------------
198
200 double time)
201{
202 G4bool sameTypeOfMolecule = SearchTimeMap(molecule);
203 return SearchUpperBoundTime(time, sameTypeOfMolecule);
204}
205
206//------------------------------------------------------------------------------
207
209 G4double time,
210 const G4ThreeVector* /*position*/,
211 int number)
212{
213 if (fDontRegister[molecule->GetDefinition()])
214 {
215 return;
216 }
217
218 if (fVerbose)
219 {
220 G4cout << "G4MoleculeCounter::AddAMoleculeAtTime : " << molecule->GetName()
221 << " at time : " << G4BestUnit(time, "Time") << G4endl;
222 }
223
224 auto counterMap_i = fCounterMap.find(molecule);
225
226 if (counterMap_i == fCounterMap.end())
227 {
228 fCounterMap[molecule][time] = number;
229 }
230 else if (counterMap_i->second.empty())
231 {
232 counterMap_i->second[time] = number;
233 }
234 else
235 {
236 NbMoleculeAgainstTime::reverse_iterator end = counterMap_i->second.rbegin();
237
238 if (end->first <= time ||
239 fabs(end->first - time) <= G4::MoleculeCounter::TimePrecision::fPrecision)
240 // Case 1 = new time comes after last recorded data
241 // Case 2 = new time is about the same as the last recorded one
242 {
243 double newValue = end->second + number;
244 counterMap_i->second[time] = newValue;
245 }
246 else
247 {
248 // if(fabs(time - G4Scheduler::Instance()->GetGlobalTime()) >
249 // G4Scheduler::Instance()->GetTimeTolerance())
250 {
252 errMsg << "Time of species "
253 << molecule->GetName() << " is "
254 << G4BestUnit(time, "Time") << " while "
255 << " global time is "
256 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
257 << G4endl;
258 G4Exception("G4MoleculeCounter::AddAMoleculeAtTime",
259 "TIME_DONT_MATCH",
260 FatalException, errMsg);
261 }
262 }
263 }
264}
265
266//------------------------------------------------------------------------------
267
269 G4double time,
270 const G4ThreeVector* /*position*/,
271 int number)
272{
273 if (fDontRegister[pMolecule->GetDefinition()])
274 {
275 return;
276 }
277
278 if (fVerbose)
279 {
280 G4cout << "G4MoleculeCounter::RemoveAMoleculeAtTime : "
281 << pMolecule->GetName() << " at time : " << G4BestUnit(time, "Time")
282 << G4endl;
283 }
284
286 {
287 if (fabs(time - G4Scheduler::Instance()->GetGlobalTime()) >
288 G4Scheduler::Instance()->GetTimeTolerance())
289 {
291 errMsg << "Time of species "
292 << pMolecule->GetName() << " is "
293 << G4BestUnit(time, "Time") << " while "
294 << " global time is "
295 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
296 << G4endl;
297 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
298 "TIME_DONT_MATCH",
299 FatalException, errMsg);
300 }
301 }
302
303 NbMoleculeAgainstTime& nbMolPerTime = fCounterMap[pMolecule];
304
305 if (nbMolPerTime.empty())
306 {
307 pMolecule->PrintState();
308 Dump();
309 G4String errMsg =
310 "You are trying to remove molecule " + pMolecule->GetName() +
311 " from the counter while this kind of molecules has not been registered yet";
312 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime", "",
313 FatalErrorInArgument, errMsg);
314
315 return;
316 }
317 else
318 {
319 NbMoleculeAgainstTime::reverse_iterator it = nbMolPerTime.rbegin();
320
321 if (it == nbMolPerTime.rend())
322 {
323 it--;
324
325 G4String errMsg =
326 "There was no " + pMolecule->GetName() + " recorded at the time or even before the time asked";
327 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime", "",
328 FatalErrorInArgument, errMsg);
329 }
330
332 {
333 Dump();
335 errMsg << "Is time going back?? " << pMolecule->GetName()
336 << " is being removed at time " << G4BestUnit(time, "Time")
337 << " while last recorded time was "
338 << G4BestUnit(it->first, "Time") << ".";
339 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
340 "RETURN_TO_THE_FUTUR",
342 errMsg);
343 }
344
345 double finalN = it->second - number;
346
347 if (finalN < 0)
348 {
349 Dump();
351 errMsg << "After removal of " << number << " species of "
352 << pMolecule->GetName() << " the final number at time "
353 << G4BestUnit(time, "Time") << " is less than zero and so not valid."
354 << " Global time is "
355 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
356 << ". Previous selected time is "
357 << G4BestUnit(it->first, "Time")
358 << G4endl;
359 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
360 "N_INF_0",
361 FatalException, errMsg);
362 }
363
364 nbMolPerTime[time] = finalN;
365 }
366}
367
368//------------------------------------------------------------------------------
369
371{
372 if (fVerbose > 1)
373 {
374 G4cout << "Entering in G4MoleculeCounter::RecordMolecules" << G4endl;
375 }
376
377 RecordedMolecules output(new ReactantList());
378
379 for (const auto & it : fCounterMap)
380 {
381 output->push_back(it.first);
382 }
383 return output;
384}
385
386//------------------------------------------------------------------------------
387
389{
390 RecordedTimes output(new std::set<G4double>);
391
392 for(const auto& it : fCounterMap)
393 {
394 for(const auto& it2 : it.second)
395 {
396 //time = it2->first;
397 output->insert(it2.first);
398 }
399 }
400
401 return output;
402}
403
404//------------------------------------------------------------------------------
405
406// >>DEV<<
407//void G4MoleculeCounter::SignalReceiver(G4SpeciesInCM* /*speciesInCM*/,
408// size_t moleculeID,
409// int /*number*/,
410// G4SpeciesInCM::SpeciesChange speciesChange,
411// int diff)
412//{
413// switch(speciesChange)
414// {
415// case G4SpeciesInCM::eAdd:
416// AddAMoleculeAtTime(G4MoleculeTable::Instance()->GetConfiguration((int)moleculeID),
417// G4Scheduler::Instance()->GetGlobalTime(),
418// diff);
419// break;
420// case G4SpeciesInCM::eRemove:
421// RemoveAMoleculeAtTime(G4MoleculeTable::Instance()->GetConfiguration((int)moleculeID),
422// G4Scheduler::Instance()->GetGlobalTime(),
423// diff);
424// break;
425// }
426//}
427
428//------------------------------------------------------------------------------
429
431{
432 for (const auto& it : fCounterMap)
433 {
434 auto pReactant = it.first;
435
436 G4cout << " --- > For " << pReactant->GetName() << G4endl;
437
438 for (const auto& it2 : it.second)
439 {
440 G4cout << " " << G4BestUnit(it2.first, "Time")
441 << " " << it2.second << G4endl;
442 }
443 }
444}
445
446//------------------------------------------------------------------------------
447
449{
450 if (fVerbose)
451 {
452 G4cout << " ---> G4MoleculeCounter::ResetCounter" << G4endl;
453 }
454 fCounterMap.clear();
455 fpLastSearch.reset(0);
456}
457
458//------------------------------------------------------------------------------
459
461{
462 return fCounterMap[molecule];
463}
464
465//------------------------------------------------------------------------------
466
468{
469 fVerbose = level;
470}
471
472//------------------------------------------------------------------------------
473
475{
476 return fVerbose;
477}
478
479//------------------------------------------------------------------------------
480
482{
483 fDontRegister[molDef] = true;
484}
485
486//------------------------------------------------------------------------------
487
489{
490 if (fDontRegister.find(molDef) == fDontRegister.end())
491 {
492 return true;
493 }
494 return false;
495}
496
497//------------------------------------------------------------------------------
498
500{
501 fDontRegister.clear();
502}
503
505{
507}
508
510{
512}
513
@ 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
std::map< G4double, G4int, G4::MoleculeCounter::TimePrecision > NbMoleculeAgainstTime
std::unique_ptr< std::set< G4double > > RecordedTimes
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
const G4MoleculeDefinition * GetDefinition() const
std::unique_ptr< ReactantList > RecordedMolecules
void Initialize() override
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
const NbMoleculeAgainstTime & GetNbMoleculeAgainstTime(Reactant *molecule)
RecordedTimes GetRecordedTimes()
static void SetTimeSlice(double)
void RegisterAll() override
G4bool fCheckTimeIsConsistentWithScheduler
void DontRegister(const G4MoleculeDefinition *) override
G4bool SearchTimeMap(Reactant *molecule)
static G4MoleculeCounter * Instance()
void ResetCounter() override
int GetNMoleculesAtTime(Reactant *molecule, double time)
G4bool IsTimeCheckedForConsistency() const
std::unique_ptr< Search > fpLastSearch
CounterMapType fCounterMap
void AddAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1) override
void RemoveAMoleculeAtTime(Reactant *, G4double time, const G4ThreeVector *position=nullptr, int number=1) override
~G4MoleculeCounter() override
RecordedMolecules GetRecordedMolecules()
void CheckTimeForConsistency(G4bool flag)
int SearchUpperBoundTime(double time, bool sameTypeOfMolecule)
bool IsRegistered(const G4MoleculeDefinition *) override
std::vector< Reactant * > ReactantList
static G4MoleculeTable * Instance()
G4ConfigurationIterator GetConfigurationIterator()
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
static G4ThreadLocal G4VMoleculeCounter * fpInstance
static G4ThreadLocal double fPrecision
bool operator()(const double &a, const double &b) const
#define G4ThreadLocal
Definition: tls.hh:77