Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VUserMoleculeReactionCounter< TIndex > Class Template Referenceabstract

#include <G4VUserMoleculeReactionCounter.hh>

Inheritance diagram for G4VUserMoleculeReactionCounter< TIndex >:

Classes

Public Member Functions

 G4VUserMoleculeReactionCounter ()
 G4VUserMoleculeReactionCounter (const G4String &, MoleculeReactionCounterType=MoleculeReactionCounterType::Basic)
 ~G4VUserMoleculeReactionCounter () override=default
void Initialize () final
void InitializeUser () override=0
void ResetCounter () override
void Dump () const override
void DumpCounterMapIndices () const override
void AbsorbCounter (const G4VMoleculeCounterInternalBase *) override
std::unique_ptr< G4VMoleculeReactionCounterIndexBuildSimpleIndex (const G4DNAMolecularReactionData *) const override=0
void RecordReaction (std::unique_ptr< G4VMoleculeReactionCounterIndex >, G4double, G4int=1) override
std::set< const G4DNAMolecularReactionData * > GetRecordedReactions () const override
std::set< G4doubleGetRecordedTimes () const override
const std::map< TIndex, InnerCounterMapType > & GetCounterMap () const
std::vector< TIndex > GetMapIndices () const
virtual G4int GetNbReactionsAtTime (const TIndex &, G4double) const
virtual G4int GetNbReactionsAtTime (Search &, const TIndex &, G4double) const
virtual std::vector< G4intGetNbReactionsAtTimes (const TIndex &, const std::vector< G4double > &) const
Public Member Functions inherited from G4VMoleculeReactionCounter
 ~G4VMoleculeReactionCounter () override=default
MoleculeReactionCounterType GetType () const
Public Member Functions inherited from G4VMoleculeCounterInternalBase
virtual ~G4VMoleculeCounterInternalBase ()=default
G4int GetId () const
void SetManagedId (G4int)
G4int GetManagedId () const
const G4StringGetName () const
G4int GetVerbose () const
void SetVerbose (G4int)
G4double GetActiveLowerBound () const
void SetActiveLowerBound (G4double, G4bool=true)
G4double GetActiveUpperBound () const
void SetActiveUpperBound (G4double, G4bool=true)
G4bool GetActiveLowerBoundInclusive () const
G4bool GetActiveUpperBoundInclusive () const
G4bool IsTimeBelowLowerBound (G4double) const
G4bool IsTimeAboveUpperBound (G4double) const
G4bool IsActiveAtGlobalTime (G4double) const
G4bool GetCheckTimeConsistencyWithScheduler () const
void SetCheckTimeConsistencyWithScheduler (G4bool=true)
G4bool GetCheckRecordedTimeConsistency () const
void SetCheckRecordedTimeConsistency (G4bool=true)
const G4MoleculeCounterTimeComparerGetTimeComparer () const
void SetTimeComparer (const G4MoleculeCounterTimeComparer &)

Protected Member Functions

G4bool SearchIndexUpdated (Search &, const TIndex &) const
G4int SearchUpperBoundTime (Search &, G4double, G4bool) const

Protected Attributes

std::map< TIndex, InnerCounterMapTypefCounterMap {}
Protected Attributes inherited from G4VMoleculeReactionCounter
MoleculeReactionCounterType fType {MoleculeReactionCounterType::Basic}
Protected Attributes inherited from G4VMoleculeCounterInternalBase
G4bool fIsInitialized {false}
G4int fId
G4int fManagedId {-1}
G4String fName {}
G4int fVerbose {0}
G4double fActiveLowerBound {0}
G4double fActiveUpperBound {std::numeric_limits<G4double>::max()}
G4bool fActiveLowerBoundInclusive {true}
G4bool fActiveUpperBoundInclusive {true}
G4bool fCheckTimeIsConsistentWithScheduler {true}
G4bool fCheckRecordedTimesAreConsistent {true}
G4MoleculeCounterTimeComparer fTimeComparer {}

Additional Inherited Members

Public Types inherited from G4VMoleculeReactionCounter
enum  MoleculeReactionCounterType { Other , Basic }
Static Public Member Functions inherited from G4VMoleculeCounterInternalBase
static void SetFixedTimePrecision (G4double)

Detailed Description

template<class TIndex>
class G4VUserMoleculeReactionCounter< TIndex >

Definition at line 40 of file G4VUserMoleculeReactionCounter.hh.

Constructor & Destructor Documentation

◆ G4VUserMoleculeReactionCounter() [1/2]

template<typename T>
G4VUserMoleculeReactionCounter< T >::G4VUserMoleculeReactionCounter ( )

◆ G4VUserMoleculeReactionCounter() [2/2]

◆ ~G4VUserMoleculeReactionCounter()

template<class TIndex>
G4VUserMoleculeReactionCounter< TIndex >::~G4VUserMoleculeReactionCounter ( )
overridedefault

Member Function Documentation

◆ AbsorbCounter()

template<typename TIndex>
void G4VUserMoleculeReactionCounter< TIndex >::AbsorbCounter ( const G4VMoleculeCounterInternalBase * pCounterBase)
overridevirtual

Implements G4VMoleculeCounterInternalBase.

Definition at line 379 of file G4VUserMoleculeReactionCounter.hh.

381{
382 if (pCounterBase == nullptr) {
384 errMsg << "Could not cast the pointer to type G4VUserMoleculeReactionCounter<"
386 << "Because the pointer is nullptr!" << G4endl;
387 G4Exception(G4String("G4VUserMoleculeReactionCounter<"
388 + G4::MoleculeCounter::GetTemplateTypeName<TIndex>() + ">::AbsorbCounter"),
389 "BAD_REFERENCE", FatalException, errMsg);
390 }
391
393
394 if (pCounter == nullptr) {
396 errMsg << "Could not cast the pointer to type G4VUserMoleculeReactionCounter<"
398 << "Because the objects aren't of the same type!" << G4endl;
399 G4Exception(G4String("G4VUserMoleculeReactionCounter<"
400 + G4::MoleculeCounter::GetTemplateTypeName<TIndex>() + ">::AbsorbCounter"),
401 "BAD_REFERENCE", FatalException, errMsg);
402 }
403
404 if (pCounter->GetType() != GetType()) {
406 errMsg << "You are trying to absorb a counter with different type!" << G4endl;
407 G4Exception(G4String("G4VUserMoleculeReactionCounter<"
408 + G4::MoleculeCounter::GetTemplateTypeName<TIndex>() + ">::AbsorbCounter"),
409 "TYPE_DIFF", JustWarning, errMsg);
410 }
411
412 for (auto const& worker_it : pCounter->GetCounterMap()) {
413 auto [master_it, indexIsNew] =
414 fCounterMap.emplace(worker_it.first, InnerCounterMapType{fTimeComparer});
415
417 for (auto const& [time, number] : worker_it.second) {
420
421 if (master_it->second.empty()) {
422 master_it->second.emplace(time, currentNumber);
423 }
424 else { // at least one element exists, so we can try to find the closest key
426 auto [it, _] = master_it->second.emplace(time, it_closest->second);
427 do {
428 it->second += currentNumber;
429 } while (++it != master_it->second.end());
430 }
431 }
432 }
433}
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
MoleculeReactionCounterType GetType() const
const std::map< TIndex, InnerCounterMapType > & GetCounterMap() const
std::map< TIndex, InnerCounterMapType > fCounterMap
std::map< TKey, TValue, TComp >::iterator FindClosestEntryForKey(std::map< TKey, TValue, TComp > &map, TKey key)

◆ BuildSimpleIndex()

template<class TIndex>
std::unique_ptr< G4VMoleculeReactionCounterIndex > G4VUserMoleculeReactionCounter< TIndex >::BuildSimpleIndex ( const G4DNAMolecularReactionData * ) const
overridepure virtual

Implements G4VMoleculeReactionCounter.

Implemented in G4MoleculeReactionCounter.

Referenced by BuildSimpleIndex().

◆ Dump()

template<typename T>
void G4VUserMoleculeReactionCounter< T >::Dump ( ) const
overridevirtual

Implements G4VMoleculeCounterInternalBase.

Definition at line 277 of file G4VUserMoleculeReactionCounter.hh.

278{
281}
void DumpCounterMapContents(const std::map< T, InnerCounterMapType > &map, G4bool includeEmpty=false)

◆ DumpCounterMapIndices()

template<typename T>
void G4VUserMoleculeReactionCounter< T >::DumpCounterMapIndices ( ) const
overridevirtual

Implements G4VMoleculeCounterInternalBase.

Definition at line 284 of file G4VUserMoleculeReactionCounter.hh.

285{
287}
void DumpCounterMapIndices(const std::map< T, InnerCounterMapType > &map, G4bool includeEmpty=false)

Referenced by Dump().

◆ GetCounterMap()

template<class TIndex>
const std::map< TIndex, InnerCounterMapType > & G4VUserMoleculeReactionCounter< TIndex >::GetCounterMap ( ) const
inline

Definition at line 75 of file G4VUserMoleculeReactionCounter.hh.

75{ return fCounterMap; }

◆ GetMapIndices()

template<typename TIndex>
std::vector< TIndex > G4VUserMoleculeReactionCounter< TIndex >::GetMapIndices ( ) const

Definition at line 243 of file G4VUserMoleculeReactionCounter.hh.

244{
245 if (fVerbose > 2) {
246 G4cout << "Entering in G4VUserMoleculeReactionCounter::GetMapIndices" << G4endl;
247 }
249}
const std::vector< T > GetMapIndices(const std::map< T, U > &_map)

◆ GetNbReactionsAtTime() [1/2]

template<typename TIndex>
G4int G4VUserMoleculeReactionCounter< TIndex >::GetNbReactionsAtTime ( const TIndex & index,
G4double time ) const
virtual

◆ GetNbReactionsAtTime() [2/2]

template<typename TIndex>
G4int G4VUserMoleculeReactionCounter< TIndex >::GetNbReactionsAtTime ( Search & search,
const TIndex & index,
G4double time ) const
virtual

Definition at line 135 of file G4VUserMoleculeReactionCounter.hh.

138{
141}
G4bool SearchIndexUpdated(Search &, const TIndex &) const
G4int SearchUpperBoundTime(Search &, G4double, G4bool) const

◆ GetNbReactionsAtTimes()

template<typename TIndex>
std::vector< G4int > G4VUserMoleculeReactionCounter< TIndex >::GetNbReactionsAtTimes ( const TIndex & index,
const std::vector< G4double > & times ) const
virtual

Definition at line 146 of file G4VUserMoleculeReactionCounter.hh.

148{
149 Search search = {};
151 for (auto time : times)
153 return counts;
154}

◆ GetRecordedReactions()

template<typename T>
std::set< const G4DNAMolecularReactionData * > G4VUserMoleculeReactionCounter< T >::GetRecordedReactions ( ) const
overridevirtual

Implements G4VMoleculeReactionCounter.

Definition at line 254 of file G4VUserMoleculeReactionCounter.hh.

255{
256 if (fVerbose > 2) {
257 G4cout << "Entering in G4VUserMoleculeReactionCounter::GetRecordedReactions" << G4endl;
258 }
260 for (const auto& it : fCounterMap) {
261 output.insert(it.first.GetReactionData());
262 }
263 return output;
264}

◆ GetRecordedTimes()

template<typename T>
std::set< G4double > G4VUserMoleculeReactionCounter< T >::GetRecordedTimes ( ) const
overridevirtual

Implements G4VMoleculeCounterInternalBase.

Definition at line 269 of file G4VUserMoleculeReactionCounter.hh.

270{
272}
std::set< G4double > GetRecordedTimes(const std::map< T, InnerCounterMapType > &map)

◆ Initialize()

template<typename TIndex>
void G4VUserMoleculeReactionCounter< TIndex >::Initialize ( )
finalvirtual

◆ InitializeUser()

template<class TIndex>
void G4VUserMoleculeReactionCounter< TIndex >::InitializeUser ( )
overridepure virtual

◆ RecordReaction()

template<class TIndex>
void G4VUserMoleculeReactionCounter< TIndex >::RecordReaction ( std::unique_ptr< G4VMoleculeReactionCounterIndex > ,
G4double ,
G4int = 1 )
overridevirtual

Implements G4VMoleculeReactionCounter.

Definition at line 159 of file G4VUserMoleculeReactionCounter.hh.

162{
163 const TIndex* mapIndex = dynamic_cast<TIndex*>(pIndex.get());
164
166 if (fVerbose > 3) {
167 G4cout << "G4VUserMoleculeReactionCounter<"
169 << ")::RecordReaction : " << mapIndex->GetReactionData()->GetReactionID()
170 << " at time : " << G4BestUnit(time, "Time") << G4endl;
171 G4cout << ":: [IsTimeAboveUpperBound] Skipping since IsTimeAboveUpperBound == true" << G4endl;
172 }
173 return;
174 }
175 else if (IsTimeBelowLowerBound(time)) {
176 if (fVerbose > 3) {
177 G4cout << "G4VUserMoleculeReactionCounter<"
179 << ")::RecordReaction : " << mapIndex->GetReactionData()->GetReactionID()
180 << " at time : " << G4BestUnit(time, "Time") << G4endl;
181 G4cout << ":: [IsTimeBelowLowerBound] Skipping since IsTimeBelowLowerBound == true" << G4endl;
182 }
183 return;
184 }
185
186 if (fVerbose > 2) {
187 G4cout << "G4VUserMoleculeReactionCounter<"
189 << ")::RecordReaction : " << mapIndex->GetReactionData()->GetReactionID()
190 << " at time : " << G4BestUnit(time, "Time") << G4endl;
191 }
192
194
195 if (indexIsNew || it->second.empty()) {
196 it->second.emplace(time, number);
197 // it->second[time] = number;
198 }
199 else {
201 // can only do consistency check if the counters are cleared before each event
202 {
203 auto end = it->second.rbegin();
204
205 if ((end->first <= time || std::fabs(end->first - time) <= fTimeComparer.GetPrecisionAtTime(time)))
206 // Case 1 = new time comes after last recorded data
207 // Case 2 = new time is about the same as the last recorded one
208 {
209 // it->second[time] = end->second + number;
210 auto [it_time, _] = it->second.emplace(time, end->second);
211 it_time->second += number;
212 }
213 else {
215 errMsg << "Time of reaction " << mapIndex->GetReactionData()->GetReactionID() << " is "
216 << G4BestUnit(time, "Time") << "while the global time is "
218 << "(last counter time: " << G4BestUnit(end->first, "Time") << ")" << G4endl;
219 G4Exception(G4String("G4VUserMoleculeReactionCounter<"
221 + ">::RecordReaction"),
222 "TIME_DONT_MATCH", FatalException, errMsg);
223 }
224 }
225 else {
226 // since counters are not cleared after chemical run (i.e., after event)
227 // there will already be numbers in the map, so...
228 // (1) find the closest time
229 // (2) emplace entry using closest value as init + number
230 // (3) add number to all "future" entries as well
232 auto [it_new, _] = it->second.emplace(time, it_closest->second);
233 do {
234 it_new->second += number;
235 } while (++it_new != it->second.end());
236 }
237 }
238}
#define G4BestUnit(a, b)

◆ ResetCounter()

template<typename T>
void G4VUserMoleculeReactionCounter< T >::ResetCounter ( )
overridevirtual

Implements G4VMoleculeCounterInternalBase.

Definition at line 292 of file G4VUserMoleculeReactionCounter.hh.

293{
294 if (fVerbose > 1) {
295 G4cout << "G4VUserMoleculeReactionCounter<" << G4::MoleculeCounter::GetTemplateTypeName<T>()
296 << ">(" << GetName() << ")::ResetCounter" << G4endl;
297 }
298 fCounterMap.clear();
299}

◆ SearchIndexUpdated()

template<typename TIndex>
G4bool G4VUserMoleculeReactionCounter< TIndex >::SearchIndexUpdated ( Search & search,
const TIndex & index ) const
protected

Definition at line 304 of file G4VUserMoleculeReactionCounter.hh.

305{
306 if (search.fLowerBoundSet && !(search.fLastIndexSearched->first < index)
308 {
309 return true;
310 }
311
312 auto mol_it = fCounterMap.find(index);
313 search.fLastIndexSearched = mol_it;
314
315 if (mol_it != fCounterMap.end()) {
316 search.fLowerBoundTime = search.fLastIndexSearched->second.end();
317 search.fLowerBoundSet = true;
318 }
319 else {
320 search.fLowerBoundSet = false;
321 }
322
323 return false;
324}

Referenced by GetNbReactionsAtTime().

◆ SearchUpperBoundTime()

template<typename T>
G4int G4VUserMoleculeReactionCounter< T >::SearchUpperBoundTime ( Search & search,
G4double time,
G4bool sameIndex ) const
protected

Definition at line 329 of file G4VUserMoleculeReactionCounter.hh.

330{
331 auto mol_it = search.fLastIndexSearched;
332 if (mol_it == fCounterMap.end()) {
333 return 0;
334 }
335
336 InnerCounterMapType const& timeMap = mol_it->second;
337 if (timeMap.empty()) {
338 return 0;
339 }
340
341 if (sameIndex) {
342 if (search.fLowerBoundSet && search.fLowerBoundTime != timeMap.end()) {
343 if (search.fLowerBoundTime->first < time) {
344 auto upperToLast = search.fLowerBoundTime;
345 upperToLast++;
346
347 if (upperToLast == timeMap.end()) {
348 return search.fLowerBoundTime->second;
349 }
350
351 if (upperToLast->first > time) {
352 return search.fLowerBoundTime->second;
353 }
354 }
355 }
356 }
357
358 auto up_time_it = timeMap.upper_bound(time);
359
360 if (up_time_it == timeMap.end()) {
361 auto last_time = timeMap.rbegin();
362 return last_time->second;
363 }
364 if (up_time_it == timeMap.begin()) {
365 return 0;
366 }
367
368 up_time_it--;
369
370 search.fLowerBoundTime = up_time_it;
371 search.fLowerBoundSet = true;
372
373 return search.fLowerBoundTime->second;
374}

Referenced by GetNbReactionsAtTime().

Member Data Documentation

◆ fCounterMap


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