Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MoleculeCounter.hh
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// Author: Mathieu Karamitros
27
28// The code is developed in the framework of the ESA AO7146
29//
30// We would be very happy hearing from you, send us your feedback! :)
31//
32// In order for Geant4-DNA to be maintained and still open-source,
33// article citations are crucial.
34// If you use Geant4-DNA chemistry and you publish papers about your software,
35// in addition to the general paper on Geant4-DNA:
36//
37// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
38//
39// we would be very happy if you could please also cite the following
40// reference papers on chemistry:
41//
42// J. Comput. Phys. 274 (2014) 841-882
43// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
44
45#pragma once
46
47#include "G4VMoleculeCounter.hh"
48#include <map>
49#include <memory>
50#include <set>
51#include <vector>
52
53//------------------------------------------------------------------------------
54
55namespace G4 {
56namespace MoleculeCounter {
58{
59 bool operator()(const double& a, const double& b) const;
61};
62}
63}
64using NbMoleculeAgainstTime = std::map<G4double, G4int, G4::MoleculeCounter::TimePrecision>;
65using RecordedTimes = std::unique_ptr<std::set<G4double>>;
66
67//------------------------------------------------------------------------------
68
70{
71 //----------------------------------------------------------------------------
72public:
73 using ReactantList = std::vector<Reactant*>;
74 using CounterMapType = std::map<Reactant*, NbMoleculeAgainstTime>;
75 using RecordedMolecules = std::unique_ptr<ReactantList>;
76
78
79 void Initialize() override;
80 void ResetCounter() override;
81
82 /* The dynamics of the given molecule won't be saved into memory.*/
83 void DontRegister(const G4MoleculeDefinition*) override;
84 bool IsRegistered(const G4MoleculeDefinition*) override;
85 void RegisterAll() override;
86
87 //----------------------------------------------------------------------------
88
89 int GetNMoleculesAtTime(Reactant* molecule, double time);
91
94
95 void SetVerbose(G4int);
97
98 /* It sets the min time difference in between two time slices. */
99 static void SetTimeSlice(double);
100
101 void Dump();
102
105
106#ifdef MOLECULE_COUNTER_TESTING
107public:
108#else
109protected:
110#endif
112 G4double time,
113 const G4ThreeVector* position = nullptr,
114 int number = 1) override;
116 G4double time,
117 const G4ThreeVector* position = nullptr,
118 int number = 1) override;
119
120 //----------------------------------------------------------------------------
121protected:
122 G4bool SearchTimeMap(Reactant* molecule);
123 int SearchUpperBoundTime(double time, bool sameTypeOfMolecule);
124
125protected:
128
130 std::map<const G4MoleculeDefinition*, G4bool> fDontRegister;
131
134
135 struct Search
136 {
138 {
139 fLowerBoundSet = false;
140 }
141 CounterMapType::iterator fLastMoleculeSearched;
142 NbMoleculeAgainstTime::iterator fLowerBoundTime;
144 };
145
146 std::unique_ptr<Search> fpLastSearch;
147
148 friend class G4Molecule;
149 friend class G4VMoleculeCounter;
150};
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
std::unique_ptr< ReactantList > RecordedMolecules
std::map< Reactant *, NbMoleculeAgainstTime > CounterMapType
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
CounterMapType::iterator fLastMoleculeSearched
NbMoleculeAgainstTime::iterator fLowerBoundTime
static G4ThreadLocal double fPrecision
bool operator()(const double &a, const double &b) const
#define G4ThreadLocal
Definition: tls.hh:77