Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MoleculeCounter Class Reference

#include <G4MoleculeCounter.hh>

+ Inheritance diagram for G4MoleculeCounter:

Classes

struct  Search
 

Public Types

using ReactantList = std::vector<Reactant*>
 
using CounterMapType = std::map<Reactant*, NbMoleculeAgainstTime>
 
using RecordedMolecules = std::unique_ptr<ReactantList>
 
- Public Types inherited from G4VMoleculeCounter
using Reactant = const G4MolecularConfiguration
 

Public Member Functions

void Initialize () override
 
void ResetCounter () override
 
void DontRegister (const G4MoleculeDefinition *) override
 
bool IsRegistered (const G4MoleculeDefinition *) override
 
void RegisterAll () override
 
int GetNMoleculesAtTime (Reactant *molecule, double time)
 
const NbMoleculeAgainstTimeGetNbMoleculeAgainstTime (Reactant *molecule)
 
RecordedMolecules GetRecordedMolecules ()
 
RecordedTimes GetRecordedTimes ()
 
void SetVerbose (G4int)
 
G4int GetVerbose ()
 
void Dump ()
 
G4bool IsTimeCheckedForConsistency () const
 
void CheckTimeForConsistency (G4bool flag)
 
- Public Member Functions inherited from G4VMoleculeCounter
void Use (G4bool flag=true)
 
G4bool InUse ()
 

Static Public Member Functions

static G4MoleculeCounterInstance ()
 
static void SetTimeSlice (double)
 
- Static Public Member Functions inherited from G4VMoleculeCounter
static void SetInstance (G4VMoleculeCounter *)
 
static void DeleteInstance ()
 
static G4VMoleculeCounterInstance ()
 
static void InitializeInstance ()
 

Protected Member Functions

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
 
G4bool SearchTimeMap (Reactant *molecule)
 
int SearchUpperBoundTime (double time, bool sameTypeOfMolecule)
 
 G4MoleculeCounter ()
 
 ~G4MoleculeCounter () override
 
- Protected Member Functions inherited from G4VMoleculeCounter
 G4VMoleculeCounter ()=default
 
virtual ~G4VMoleculeCounter ()=default
 

Protected Attributes

CounterMapType fCounterMap
 
std::map< const G4MoleculeDefinition *, G4boolfDontRegister
 
G4int fVerbose
 
G4bool fCheckTimeIsConsistentWithScheduler
 
std::unique_ptr< SearchfpLastSearch
 
- Protected Attributes inherited from G4VMoleculeCounter
G4bool fUse = false
 

Friends

class G4Molecule
 
class G4VMoleculeCounter
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VMoleculeCounter
static G4ThreadLocal G4VMoleculeCounterfpInstance = nullptr
 

Detailed Description

Definition at line 69 of file G4MoleculeCounter.hh.

Member Typedef Documentation

◆ CounterMapType

◆ ReactantList

Definition at line 73 of file G4MoleculeCounter.hh.

◆ RecordedMolecules

Definition at line 75 of file G4MoleculeCounter.hh.

Constructor & Destructor Documentation

◆ G4MoleculeCounter()

G4MoleculeCounter::G4MoleculeCounter ( )
protected

Definition at line 72 of file G4MoleculeCounter.cc.

Referenced by Instance().

◆ ~G4MoleculeCounter()

G4MoleculeCounter::~G4MoleculeCounter ( )
overrideprotecteddefault

Member Function Documentation

◆ AddAMoleculeAtTime()

void G4MoleculeCounter::AddAMoleculeAtTime ( Reactant * molecule,
G4double time,
const G4ThreeVector * position = nullptr,
int number = 1 )
overrideprotectedvirtual

Implements G4VMoleculeCounter.

Definition at line 209 of file G4MoleculeCounter.cc.

213{
214 if (fDontRegister[molecule->GetDefinition()])
215 {
216 return;
217 }
218
219 if (fVerbose != 0)
220 {
221 G4cout << "G4MoleculeCounter::AddAMoleculeAtTime : " << molecule->GetName()
222 << " at time : " << G4BestUnit(time, "Time") << G4endl;
223 }
224
225 auto counterMap_i = fCounterMap.find(molecule);
226
227 if (counterMap_i == fCounterMap.end())
228 {
229 fCounterMap[molecule][time] = number;
230 }
231 else if (counterMap_i->second.empty())
232 {
233 counterMap_i->second[time] = number;
234 }
235 else
236 {
237 auto end = counterMap_i->second.rbegin();
238
239 if (end->first <= time ||
240 fabs(end->first - time) <= G4::MoleculeCounter::TimePrecision::fPrecision)
241 // Case 1 = new time comes after last recorded data
242 // Case 2 = new time is about the same as the last recorded one
243 {
244 double newValue = end->second + number;
245 counterMap_i->second[time] = newValue;
246 }
247 else
248 {
249 // if(fabs(time - G4Scheduler::Instance()->GetGlobalTime()) >
250 // G4Scheduler::Instance()->GetTimeTolerance())
251 {
253 errMsg << "Time of species "
254 << molecule->GetName() << " is "
255 << G4BestUnit(time, "Time") << " while "
256 << " global time is "
257 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
258 << G4endl;
259 G4Exception("G4MoleculeCounter::AddAMoleculeAtTime",
260 "TIME_DONT_MATCH",
261 FatalException, errMsg);
262 }
263 }
264 }
265}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4BestUnit(a, b)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
CounterMapType fCounterMap
static G4Scheduler * Instance()
static G4ThreadLocal double fPrecision

◆ CheckTimeForConsistency()

void G4MoleculeCounter::CheckTimeForConsistency ( G4bool flag)

Definition at line 504 of file G4MoleculeCounter.cc.

505{
507}

◆ DontRegister()

void G4MoleculeCounter::DontRegister ( const G4MoleculeDefinition * molDef)
overridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 480 of file G4MoleculeCounter.cc.

481{
482 fDontRegister[molDef] = true;
483}

◆ Dump()

void G4MoleculeCounter::Dump ( )

Definition at line 429 of file G4MoleculeCounter.cc.

430{
431 for (const auto& it : fCounterMap)
432 {
433 auto pReactant = it.first;
434
435 G4cout << " --- > For " << pReactant->GetName() << G4endl;
436
437 for (const auto& it2 : it.second)
438 {
439 G4cout << " " << G4BestUnit(it2.first, "Time")
440 << " " << it2.second << G4endl;
441 }
442 }
443}

Referenced by RemoveAMoleculeAtTime().

◆ GetNbMoleculeAgainstTime()

const NbMoleculeAgainstTime & G4MoleculeCounter::GetNbMoleculeAgainstTime ( Reactant * molecule)

Definition at line 459 of file G4MoleculeCounter.cc.

460{
461 return fCounterMap[molecule];
462}

◆ GetNMoleculesAtTime()

int G4MoleculeCounter::GetNMoleculesAtTime ( Reactant * molecule,
double time )

Definition at line 200 of file G4MoleculeCounter.cc.

202{
203 G4bool sameTypeOfMolecule = SearchTimeMap(molecule);
204 return SearchUpperBoundTime(time, sameTypeOfMolecule);
205}
bool G4bool
Definition G4Types.hh:86
G4bool SearchTimeMap(Reactant *molecule)
int SearchUpperBoundTime(double time, bool sameTypeOfMolecule)

Referenced by G4DNAEventScheduler::ClearAndReChargeCounter().

◆ GetRecordedMolecules()

G4MoleculeCounter::RecordedMolecules G4MoleculeCounter::GetRecordedMolecules ( )

Definition at line 369 of file G4MoleculeCounter.cc.

370{
371 if (fVerbose > 1)
372 {
373 G4cout << "Entering in G4MoleculeCounter::RecordMolecules" << G4endl;
374 }
375
376 RecordedMolecules output(new ReactantList());
377
378 for (const auto & it : fCounterMap)
379 {
380 output->push_back(it.first);
381 }
382 return output;
383}
std::unique_ptr< ReactantList > RecordedMolecules
std::vector< Reactant * > ReactantList

Referenced by G4DNAEventScheduler::ClearAndReChargeCounter().

◆ GetRecordedTimes()

RecordedTimes G4MoleculeCounter::GetRecordedTimes ( )

Definition at line 387 of file G4MoleculeCounter.cc.

388{
389 RecordedTimes output(new std::set<G4double>);
390
391 for(const auto& it : fCounterMap)
392 {
393 for(const auto& it2 : it.second)
394 {
395 //time = it2->first;
396 output->insert(it2.first);
397 }
398 }
399
400 return output;
401}
std::unique_ptr< std::set< G4double > > RecordedTimes

◆ GetVerbose()

G4int G4MoleculeCounter::GetVerbose ( )

Definition at line 473 of file G4MoleculeCounter.cc.

474{
475 return fVerbose;
476}

◆ Initialize()

void G4MoleculeCounter::Initialize ( )
overridevirtual

Implements G4VMoleculeCounter.

Definition at line 84 of file G4MoleculeCounter.cc.

85{
87 while ((mol_iterator)())
88 {
89 if (!IsRegistered(mol_iterator.value()->GetDefinition()))
90 {
91 continue;
92 }
93
94 fCounterMap[mol_iterator.value()]; // initialize the second map
95 }
96}
bool IsRegistered(const G4MoleculeDefinition *) override
static G4MoleculeTable * Instance()
G4ConfigurationIterator GetConfigurationIterator()

◆ Instance()

G4MoleculeCounter * G4MoleculeCounter::Instance ( )
static

Definition at line 61 of file G4MoleculeCounter.cc.

62{
63 if (fpInstance == nullptr)
64 {
66 }
67 return dynamic_cast<G4MoleculeCounter*>(fpInstance);
68}
static G4ThreadLocal G4VMoleculeCounter * fpInstance

Referenced by G4DNAEventScheduler::ClearAndReChargeCounter().

◆ IsRegistered()

bool G4MoleculeCounter::IsRegistered ( const G4MoleculeDefinition * molDef)
overridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 487 of file G4MoleculeCounter.cc.

488{
489 return fDontRegister.find(molDef) == fDontRegister.end();
490}

Referenced by Initialize().

◆ IsTimeCheckedForConsistency()

G4bool G4MoleculeCounter::IsTimeCheckedForConsistency ( ) const

Definition at line 499 of file G4MoleculeCounter.cc.

500{
502}

◆ RegisterAll()

void G4MoleculeCounter::RegisterAll ( )
overridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 494 of file G4MoleculeCounter.cc.

495{
496 fDontRegister.clear();
497}

◆ RemoveAMoleculeAtTime()

void G4MoleculeCounter::RemoveAMoleculeAtTime ( Reactant * ,
G4double time,
const G4ThreeVector * position = nullptr,
int number = 1 )
overrideprotectedvirtual

Implements G4VMoleculeCounter.

Definition at line 269 of file G4MoleculeCounter.cc.

273{
274 if (fDontRegister[pMolecule->GetDefinition()])
275 {
276 return;
277 }
278
279 if (fVerbose != 0)
280 {
281 G4cout << "G4MoleculeCounter::RemoveAMoleculeAtTime : "
282 << pMolecule->GetName() << " at time : " << G4BestUnit(time, "Time")
283 << G4endl;
284 }
285
287 {
288 if (fabs(time - G4Scheduler::Instance()->GetGlobalTime()) >
289 G4Scheduler::Instance()->GetTimeTolerance())
290 {
292 errMsg << "Time of species "
293 << pMolecule->GetName() << " is "
294 << G4BestUnit(time, "Time") << " while "
295 << " global time is "
296 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
297 << G4endl;
298 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
299 "TIME_DONT_MATCH",
300 FatalException, errMsg);
301 }
302 }
303
304 NbMoleculeAgainstTime& nbMolPerTime = fCounterMap[pMolecule];
305
306 if (nbMolPerTime.empty())
307 {
308 pMolecule->PrintState();
309 Dump();
310 G4String errMsg =
311 "You are trying to remove molecule " + pMolecule->GetName() +
312 " from the counter while this kind of molecules has not been registered yet";
313 G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime", "",
314 FatalErrorInArgument, errMsg);
315
316 return;
317 }
318
319 auto 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}
@ FatalErrorInArgument
std::map< G4double, G4int, G4::MoleculeCounter::TimePrecision > NbMoleculeAgainstTime

◆ ResetCounter()

void G4MoleculeCounter::ResetCounter ( )
overridevirtual

Implements G4VMoleculeCounter.

Definition at line 447 of file G4MoleculeCounter.cc.

448{
449 if (fVerbose != 0)
450 {
451 G4cout << " ---> G4MoleculeCounter::ResetCounter" << G4endl;
452 }
453 fCounterMap.clear();
454 fpLastSearch.reset(nullptr);
455}
std::unique_ptr< Search > fpLastSearch

Referenced by G4DNAEventScheduler::ClearAndReChargeCounter().

◆ SearchTimeMap()

G4bool G4MoleculeCounter::SearchTimeMap ( Reactant * molecule)
protected

Definition at line 107 of file G4MoleculeCounter.cc.

108{
109 if (fpLastSearch == nullptr)
110 {
111 fpLastSearch = std::make_unique<Search>();
112 }
113 else
114 {
115 if (fpLastSearch->fLowerBoundSet &&
116 fpLastSearch->fLastMoleculeSearched->first == molecule)
117 {
118 return true;
119 }
120 }
121
122 auto mol_it = fCounterMap.find(molecule);
123 fpLastSearch->fLastMoleculeSearched = mol_it;
124
125 if (mol_it != fCounterMap.end())
126 {
127 fpLastSearch->fLowerBoundTime = fpLastSearch->fLastMoleculeSearched->second
128 .end();
129 fpLastSearch->fLowerBoundSet = true;
130 }
131 else
132 {
133 fpLastSearch->fLowerBoundSet = false;
134 }
135
136 return false;
137}

Referenced by GetNMoleculesAtTime().

◆ SearchUpperBoundTime()

int G4MoleculeCounter::SearchUpperBoundTime ( double time,
bool sameTypeOfMolecule )
protected

Definition at line 141 of file G4MoleculeCounter.cc.

143{
144 auto mol_it = fpLastSearch->fLastMoleculeSearched;
145 if (mol_it == fCounterMap.end())
146 {
147 return 0;
148 }
149
150 NbMoleculeAgainstTime& timeMap = mol_it->second;
151 if (timeMap.empty())
152 {
153 return 0;
154 }
155
156 if (sameTypeOfMolecule)
157 {
158 if (fpLastSearch->fLowerBoundSet && fpLastSearch->fLowerBoundTime != timeMap.end())
159 {
160 if (fpLastSearch->fLowerBoundTime->first < time)
161 {
162 auto upperToLast = fpLastSearch->fLowerBoundTime;
163 upperToLast++;
164
165 if (upperToLast == timeMap.end())
166 {
167 return fpLastSearch->fLowerBoundTime->second;
168 }
169
170 if (upperToLast->first > time)
171 {
172 return fpLastSearch->fLowerBoundTime->second;
173 }
174 }
175 }
176 }
177
178 auto up_time_it = timeMap.upper_bound(time);
179
180 if (up_time_it == timeMap.end())
181 {
182 auto last_time = timeMap.rbegin();
183 return last_time->second;
184 }
185 if (up_time_it == timeMap.begin())
186 {
187 return 0;
188 }
189
190 up_time_it--;
191
192 fpLastSearch->fLowerBoundTime = up_time_it;
193 fpLastSearch->fLowerBoundSet = true;
194
195 return fpLastSearch->fLowerBoundTime->second;
196}

Referenced by GetNMoleculesAtTime().

◆ SetTimeSlice()

void G4MoleculeCounter::SetTimeSlice ( double timeSlice)
static

Definition at line 100 of file G4MoleculeCounter.cc.

◆ SetVerbose()

void G4MoleculeCounter::SetVerbose ( G4int level)

Definition at line 466 of file G4MoleculeCounter.cc.

467{
468 fVerbose = level;
469}

Friends And Related Symbol Documentation

◆ G4Molecule

friend class G4Molecule
friend

Definition at line 148 of file G4MoleculeCounter.hh.

◆ G4VMoleculeCounter

friend class G4VMoleculeCounter
friend

Definition at line 149 of file G4MoleculeCounter.hh.

Member Data Documentation

◆ fCheckTimeIsConsistentWithScheduler

G4bool G4MoleculeCounter::fCheckTimeIsConsistentWithScheduler
protected

◆ fCounterMap

◆ fDontRegister

std::map<const G4MoleculeDefinition*, G4bool> G4MoleculeCounter::fDontRegister
protected

◆ fpLastSearch

std::unique_ptr<Search> G4MoleculeCounter::fpLastSearch
protected

Definition at line 146 of file G4MoleculeCounter.hh.

Referenced by ResetCounter(), SearchTimeMap(), and SearchUpperBoundTime().

◆ fVerbose

G4int G4MoleculeCounter::fVerbose
protected

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