Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4MoleculeCounter Class Reference

#include <G4MoleculeCounter.hh>

Inheritance diagram for G4MoleculeCounter:
Collaboration diagram for G4MoleculeCounter:

Classes

struct  Search
 

Public Types

typedef std::map
< G4MolecularConfiguration
*, NbMoleculeAgainstTime
CounterMapType
 
typedef std::unique_ptr
< std::vector
< G4MolecularConfiguration * > > 
RecordedMolecules
 

Public Member Functions

void Initialize () override
 
void ResetCounter () override
 
void AddAMoleculeAtTime (G4MolecularConfiguration *, G4double time, const G4ThreeVector *position=nullptr, int number=1) override
 
void RemoveAMoleculeAtTime (G4MolecularConfiguration *, G4double time, const G4ThreeVector *position=nullptr, int number=1) override
 
void DontRegister (const G4MoleculeDefinition *) override
 
bool IsRegistered (const G4MoleculeDefinition *) override
 
void RegisterAll () override
 
int GetNMoleculesAtTime (G4MolecularConfiguration *molecule, double time)
 
const NbMoleculeAgainstTimeGetNbMoleculeAgainstTime (G4MolecularConfiguration *molecule)
 
RecordedMolecules GetRecordedMolecules ()
 
RecordedTimes GetRecordedTimes ()
 
void SetVerbose (G4int)
 
G4int GetVerbose ()
 
void SetTimeSlice (double)
 
void Dump ()
 
G4bool IsTimeCheckedForConsistency () const
 
void CheckTimeForConsistency (G4bool flag)
 

Static Public Member Functions

static G4MoleculeCounterInstance ()
 
- Static Public Member Functions inherited from G4VMoleculeCounter
static void SetInstance (G4VMoleculeCounter *)
 
static void DeleteInstance ()
 
static G4VMoleculeCounterInstance ()
 
static void InitializeInstance ()
 
static void Use (G4bool flag=true)
 
static G4bool InUse ()
 

Protected Member Functions

 G4MoleculeCounter ()
 
virtual ~G4MoleculeCounter ()
 
G4bool SearchTimeMap (G4MolecularConfiguration *molecule)
 
int SearchUpperBoundTime (double time, bool sameTypeOfMolecule)
 
- Protected Member Functions inherited from G4VMoleculeCounter
 G4VMoleculeCounter ()
 
virtual ~G4VMoleculeCounter ()
 

Protected Attributes

CounterMapType fCounterMap
 
std::map< const
G4MoleculeDefinition *, G4bool
fDontRegister
 
G4int fVerbose
 
G4bool fCheckTimeIsConsistentWithScheduler
 
std::unique_ptr< SearchfpLastSearch
 

Friends

class G4Molecule
 
class G4VMoleculeCounter
 

Additional Inherited Members

- Static Protected Attributes inherited from G4VMoleculeCounter
static G4ThreadLocal
G4VMoleculeCounter
fpInstance = 0
 
static G4bool fUse = false
 

Detailed Description

Definition at line 79 of file G4MoleculeCounter.hh.

Member Typedef Documentation

typedef std::unique_ptr<std::vector<G4MolecularConfiguration*> > G4MoleculeCounter::RecordedMolecules

Definition at line 86 of file G4MoleculeCounter.hh.

Constructor & Destructor Documentation

G4MoleculeCounter::G4MoleculeCounter ( )
protected

Definition at line 51 of file G4MoleculeCounter.cc.

52 {
53  fVerbose = 0;
56  {
58  }
59 }
static constexpr double picosecond
Definition: G4SIunits.hh:161
static G4ThreadLocal double fPrecision
G4bool fCheckTimeIsConsistentWithScheduler
G4MoleculeCounter::~G4MoleculeCounter ( )
protectedvirtual

Definition at line 63 of file G4MoleculeCounter.cc.

64 {
65 }

Member Function Documentation

void G4MoleculeCounter::AddAMoleculeAtTime ( G4MolecularConfiguration molecule,
G4double  time,
const G4ThreeVector position = nullptr,
int  number = 1 
)
overridevirtual

Implements G4VMoleculeCounter.

Definition at line 217 of file G4MoleculeCounter.cc.

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

Here is the call graph for this function:

void G4MoleculeCounter::CheckTimeForConsistency ( G4bool  flag)
inline

Definition at line 170 of file G4MoleculeCounter.hh.

171  {
173  }
G4bool fCheckTimeIsConsistentWithScheduler
void G4MoleculeCounter::DontRegister ( const G4MoleculeDefinition molDef)
inlineoverridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 218 of file G4MoleculeCounter.hh.

219 {
220  fDontRegister[molDef] = true;
221 }
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
void G4MoleculeCounter::Dump ( )

Definition at line 440 of file G4MoleculeCounter.cc.

441 {
442  CounterMapType::iterator it = fCounterMap.begin();
443  CounterMapType::iterator end = fCounterMap.end();
444 
445  for(;it!=end;++it)
446  {
447  G4MolecularConfiguration* molConf = it->first;
448 
449  G4cout << " --- > For " << molConf->GetName() << G4endl;
450  NbMoleculeAgainstTime::iterator it2 = it->second.begin();
451  NbMoleculeAgainstTime::iterator end2 = it->second.end();
452 
453  for(;it2!=end2;++it2)
454  {
455  G4cout << " " << G4BestUnit(it2->first, "Time")
456  << " " << it2->second << G4endl;
457  }
458  }
459 }
const G4String & GetName() const
CounterMapType fCounterMap
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

const NbMoleculeAgainstTime & G4MoleculeCounter::GetNbMoleculeAgainstTime ( G4MolecularConfiguration molecule)
inline

Definition at line 196 of file G4MoleculeCounter.hh.

197 {
198  return fCounterMap[molecule];
199 }
CounterMapType fCounterMap
int G4MoleculeCounter::GetNMoleculesAtTime ( G4MolecularConfiguration molecule,
double  time 
)

Definition at line 208 of file G4MoleculeCounter.cc.

210 {
211  G4bool sameTypeOfMolecule = SearchTimeMap(molecule);
212  return SearchUpperBoundTime(time, sameTypeOfMolecule);
213 }
int SearchUpperBoundTime(double time, bool sameTypeOfMolecule)
G4bool SearchTimeMap(G4MolecularConfiguration *molecule)
bool G4bool
Definition: G4Types.hh:79
G4MoleculeCounter::RecordedMolecules G4MoleculeCounter::GetRecordedMolecules ( )

Definition at line 368 of file G4MoleculeCounter.cc.

369 {
370  if (fVerbose > 1)
371  {
372  G4cout << "Entering in G4MoleculeCounter::RecordMolecules" << G4endl;
373  }
374 
375  CounterMapType::iterator it;
376  RecordedMolecules output (new vector<G4MolecularConfiguration*>);
377 
378  for(it = fCounterMap.begin(); it != fCounterMap.end(); it++)
379  {
380  output->push_back(it->first);
381  }
382  return output;
383 }
CounterMapType fCounterMap
std::unique_ptr< std::vector< G4MolecularConfiguration * > > RecordedMolecules
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
RecordedTimes G4MoleculeCounter::GetRecordedTimes ( )

Definition at line 387 of file G4MoleculeCounter.cc.

388 {
389  RecordedTimes output(new std::set<G4double>);
390 
391  //G4double time;
392 
393  CounterMapType::iterator it;
394  CounterMapType::const_iterator ite;
395 
396  NbMoleculeAgainstTime::iterator it2;
397  NbMoleculeAgainstTime::const_iterator ite2;
398 
399  // iterate on each molecule
400  for (it = fCounterMap.begin(), ite = fCounterMap.end(); it != ite; ++it)
401  {
402  // iterate on each time
403  for (it2 = (it->second).begin(), ite2 = (it->second).end(); it2 != ite2;
404  ++it2)
405  {
406  //time = it2->first;
407  output->insert(it2->first);
408  }
409  }
410 
411  return output;
412 }
CounterMapType fCounterMap
std::unique_ptr< std::set< G4double > > RecordedTimes
G4int G4MoleculeCounter::GetVerbose ( )
inline

Definition at line 210 of file G4MoleculeCounter.hh.

211 {
212  return fVerbose;
213 }
void G4MoleculeCounter::Initialize ( )
overridevirtual

Implements G4VMoleculeCounter.

Definition at line 69 of file G4MoleculeCounter.cc.

70 {
71 // G4cout << "G4MoleculeCounter::Initialize" << G4endl;
72 
75  while ((mol_iterator)())
76  {
77  if(IsRegistered(mol_iterator.value()->GetDefinition()) == false)
78  {
79  continue;
80  }
81 
82  // G4cout << "G4MoleculeCounter::Initialize" << G4endl;
83  // G4cout << mol_iterator->value()->GetName() << G4endl;
84  fCounterMap[mol_iterator.value()]; // initialize the second map
85  }
86 }
CounterMapType fCounterMap
bool IsRegistered(const G4MoleculeDefinition *) override
static G4MoleculeTable * Instance()
G4ConfigurationIterator GetConfigurationIterator()

Here is the call graph for this function:

G4MoleculeCounter * G4MoleculeCounter::Instance ( )
static

Definition at line 44 of file G4MoleculeCounter.cc.

44  {
46  return dynamic_cast<G4MoleculeCounter*>(fpInstance);
47 }
static G4ThreadLocal G4VMoleculeCounter * fpInstance
bool G4MoleculeCounter::IsRegistered ( const G4MoleculeDefinition molDef)
inlineoverridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 226 of file G4MoleculeCounter.hh.

227 {
228  if(fDontRegister.find(molDef) == fDontRegister.end()) return true;
229  return false;
230 }
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
G4bool G4MoleculeCounter::IsTimeCheckedForConsistency ( ) const
inline

Definition at line 165 of file G4MoleculeCounter.hh.

166  {
168  }
G4bool fCheckTimeIsConsistentWithScheduler
void G4MoleculeCounter::RegisterAll ( )
inlineoverridevirtual

Reimplemented from G4VMoleculeCounter.

Definition at line 234 of file G4MoleculeCounter.hh.

235 {
236  fDontRegister.clear();
237 }
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
void G4MoleculeCounter::RemoveAMoleculeAtTime ( G4MolecularConfiguration molecule,
G4double  time,
const G4ThreeVector position = nullptr,
int  number = 1 
)
overridevirtual

Implements G4VMoleculeCounter.

Definition at line 272 of file G4MoleculeCounter.cc.

276 {
277  if (fDontRegister[molecule->GetDefinition()]) return;
278 
279  if (fVerbose)
280  {
281  G4cout << "G4MoleculeCounter::RemoveAMoleculeAtTime : "
282  << molecule->GetName() << " at time : " << G4BestUnit(time, "Time")
283  << G4endl;
284  }
285 
287  {
288  if(fabs(time - G4Scheduler::Instance()->GetGlobalTime()) >
289  G4Scheduler::Instance()->GetTimeTolerance())
290  {
291  G4ExceptionDescription errMsg;
292  errMsg << "Time of species "
293  << molecule->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[molecule];
305 
306  if (nbMolPerTime.empty())
307  {
308  molecule->PrintState();
309  Dump();
310  G4String errMsg =
311  "You are trying to remove molecule " + molecule->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  else
319  {
320  NbMoleculeAgainstTime::reverse_iterator it = nbMolPerTime.rbegin();
321 
322  if (it == nbMolPerTime.rend()){
323  it--;
324 
325  G4String errMsg = "There was no " + molecule->GetName()
326  + " recorded at the time or even before the time asked";
327  G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime", "",
328  FatalErrorInArgument, errMsg);
329  }
330 
331  if (time - it->first < -compDoubleWithPrecision::fPrecision){
332  Dump();
333  G4ExceptionDescription errMsg;
334  errMsg << "Is time going back?? " << molecule->GetName()
335  << " is being removed at time " << G4BestUnit(time, "Time")
336  << " while last recorded time was "
337  << G4BestUnit(it->first, "Time") << ".";
338  G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
339  "RETURN_TO_THE_FUTUR",
341  errMsg);
342  }
343 
344  double finalN = it->second - number;
345 
346  if(finalN < 0){
347  Dump();
348  G4ExceptionDescription errMsg;
349  errMsg << "After removal of " << number << " species of "
350  << molecule->GetName() << " the final number at time "
351  << G4BestUnit(time, "Time") << " is less than zero and so not valid."
352  << " Global time is "
353  << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
354  << ". Previous selected time is "
355  << G4BestUnit(it->first, "Time")
356  << G4endl;
357  G4Exception("G4MoleculeCounter::RemoveAMoleculeAtTime",
358  "N_INF_0",
359  FatalException, errMsg);
360  }
361 
362  nbMolPerTime[time] = finalN;
363  }
364 }
const G4String & GetName() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CounterMapType fCounterMap
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:102
std::map< const G4MoleculeDefinition *, G4bool > fDontRegister
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ThreadLocal double fPrecision
const G4MoleculeDefinition * GetDefinition() const
std::map< G4double, G4int, compDoubleWithPrecision > NbMoleculeAgainstTime
#define G4endl
Definition: G4ios.hh:61
G4bool fCheckTimeIsConsistentWithScheduler

Here is the call graph for this function:

void G4MoleculeCounter::ResetCounter ( )
inlineoverridevirtual

Implements G4VMoleculeCounter.

Definition at line 183 of file G4MoleculeCounter.hh.

184 {
185  if(fVerbose)
186  {
187  G4cout << " ---> G4MoleculeCounter::ResetCounter" << G4endl;
188  }
189  fCounterMap.clear();
190  fpLastSearch.reset(0);
191 }
CounterMapType fCounterMap
std::unique_ptr< Search > fpLastSearch
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4bool G4MoleculeCounter::SearchTimeMap ( G4MolecularConfiguration molecule)
protected

Definition at line 97 of file G4MoleculeCounter.cc.

98 {
99  if (fpLastSearch.get() == 0)
100  {
101  fpLastSearch.reset(new Search());
102  }
103  else
104  {
105  if (fpLastSearch->fLowerBoundSet &&
106  fpLastSearch->fLastMoleculeSearched->first == molecule)
107  return true;
108  }
109 
110  CounterMapType::iterator mol_it = fCounterMap.find(molecule);
111  fpLastSearch->fLastMoleculeSearched = mol_it;
112 
113  if (mol_it != fCounterMap.end()) // TODO
114  {
115  fpLastSearch->fLowerBoundTime = fpLastSearch->fLastMoleculeSearched->second
116  .end();
117  fpLastSearch->fLowerBoundSet = true;
118  }
119  else
120  {
121  fpLastSearch->fLowerBoundSet = false;
122  }
123 
124  return false;
125 }
CounterMapType fCounterMap
std::unique_ptr< Search > fpLastSearch
int G4MoleculeCounter::SearchUpperBoundTime ( double  time,
bool  sameTypeOfMolecule 
)
protected

Definition at line 129 of file G4MoleculeCounter.cc.

131 {
132  CounterMapType::iterator mol_it = fpLastSearch->fLastMoleculeSearched;
133  if (mol_it == fCounterMap.end()) return 0; // RETURN
134 
135  NbMoleculeAgainstTime& timeMap = mol_it->second;
136  if (timeMap.empty()) return 0;
137 
138  NbMoleculeAgainstTime::iterator end_time = timeMap.end();
139 
140  if (sameTypeOfMolecule == true)
141  {
142  //G4cout << "SAME MOLECULE" << G4endl;
143  if (fpLastSearch->fLowerBoundSet && fpLastSearch->fLowerBoundTime
144  != end_time)
145  {
146  if (fpLastSearch->fLowerBoundTime->first < time)
147  {
148  NbMoleculeAgainstTime::iterator upperToLast = fpLastSearch
149  ->fLowerBoundTime;
150  upperToLast++;
151 
152  if (upperToLast == end_time)
153  {
154  return fpLastSearch->fLowerBoundTime->second;
155  }
156 
157  if (upperToLast->first > time)
158  {
159  return fpLastSearch->fLowerBoundTime->second;
160  }
161  }
162  }
163  }
164  /*
165  else
166  {
167  G4cout << "--> Molecule has changed" << G4endl;
168  }
169  */
170  //G4cout << "Searching" << G4endl;
171  // With upper bound
172  NbMoleculeAgainstTime::iterator up_time_it = timeMap.upper_bound(time);
173 
174  if (up_time_it == end_time)
175  {
176  NbMoleculeAgainstTime::reverse_iterator last_time = timeMap.rbegin();
177 
178 // {
179  //G4cout << "RETURN LAST : " << G4BestUnit(time, "Time") << G4endl;
180  return last_time->second;
181 // }
182 
183 // G4cout << "RETURN 0 (1)" << G4endl;
184 // return 0; // RETURN
185  }
186  if (up_time_it == timeMap.begin())
187  {
188 // G4cout << "RETURN 0 (2)" << G4endl;
189  return 0; // RETURN
190  }
191 
192  //G4cout << "Going back : " << up_time_it->first << "-->";
193 
194  up_time_it--;
195 
196 // G4cout << up_time_it->first << G4endl;
197 
198  fpLastSearch->fLowerBoundTime = up_time_it;
199  fpLastSearch->fLowerBoundSet = true;
200 
201 // G4cout << "returning : " << fpLastSearch->fLowerBoundTime->second << G4endl;
202 
203  return fpLastSearch->fLowerBoundTime->second;
204 }
CounterMapType fCounterMap
std::unique_ptr< Search > fpLastSearch
std::map< G4double, G4int, compDoubleWithPrecision > NbMoleculeAgainstTime
void G4MoleculeCounter::SetTimeSlice ( double  timeSlice)

Definition at line 90 of file G4MoleculeCounter.cc.

91 {
93 }
static G4ThreadLocal double fPrecision
void G4MoleculeCounter::SetVerbose ( G4int  level)
inline

Definition at line 203 of file G4MoleculeCounter.hh.

204 {
205  fVerbose = level;
206 }

Friends And Related Function Documentation

friend class G4Molecule
friend

Definition at line 117 of file G4MoleculeCounter.hh.

friend class G4VMoleculeCounter
friend

Definition at line 118 of file G4MoleculeCounter.hh.

Member Data Documentation

G4bool G4MoleculeCounter::fCheckTimeIsConsistentWithScheduler
protected

Definition at line 97 of file G4MoleculeCounter.hh.

CounterMapType G4MoleculeCounter::fCounterMap
protected

Definition at line 93 of file G4MoleculeCounter.hh.

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

Definition at line 94 of file G4MoleculeCounter.hh.

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

Definition at line 109 of file G4MoleculeCounter.hh.

G4int G4MoleculeCounter::fVerbose
protected

Definition at line 96 of file G4MoleculeCounter.hh.


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