Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAScavengerMaterial.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
27#ifndef G4DNASCAVENGERMATERIAL_HH
28#define G4DNASCAVENGERMATERIAL_HH
29#include "globals.hh"
30#include "G4ios.hh"
31#include <map>
32#include <vector>
33#include "G4MoleculeCounter.hh"
35#include "G4MoleculeTable.hh"
36
37class G4Material;
40
42{
43 public:
45 std::map<G4double, int64_t, G4::MoleculeCounter::TimePrecision>;
47 using MaterialMap = std::map<MolType, int64_t>;
48 using ReactantList = std::vector<MolType>;
49 using CounterMapType = std::map<MolType, NbMoleculeInTime>;
52 ~G4DNAScavengerMaterial() override = default;
55 void Initialize();
56
60
62 const G4ThreeVector* position = nullptr,
63 G4int number = 1);
65 const G4ThreeVector* position = nullptr,
66 G4int number = 1);
67
68 void Reset() override;
69
70 void PrintInfo();
71
72 MaterialMap::iterator end() { return fScavengerTable.end(); }
73 MaterialMap::iterator begin() { return fScavengerTable.begin(); }
74 size_t size() { return fScavengerTable.size(); }
75
77 {
78 auto it = fScavengerTable.find(type);
79 if(it != fScavengerTable.end())
80 {
81 return it->second > 0;
82 }
83
84 return false;
85
86 }
87
88 void SetCounterAgainstTime() { fCounterAgainstTime = true; }
89 void SetpH(const G4int& ph);
91
92 std::vector<MolType> GetScavengerList() const
93 {
94 std::vector<MolType> output;
95 for(const auto& it : fScavengerTable)
96 {
97 output.push_back(it.first);
98 }
99 return output;
100 }
101
102 void Dump();
103 int64_t GetNMoleculesAtTime(MolType molecule, G4double time);
104 G4bool SearchTimeMap(MolType molecule);
105 int64_t SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule);
106
107 private:
108 G4VChemistryWorld* fpChemistryInfo;
109 G4bool fIsInitialized;
110 MaterialMap fScavengerTable;
111 CounterMapType fCounterMap;
112 G4bool fCounterAgainstTime;
113 G4int fVerbose;
117 struct Search
118 {
119 Search() { fLowerBoundSet = false; }
120 CounterMapType::iterator fLastMoleculeSearched;
121 NbMoleculeInTime::iterator fLowerBoundTime;
122 G4bool fLowerBoundSet;
123 };
124
125 std::unique_ptr<Search> fpLastSearch;
126 void WaterEquilibrium();
127
128};
129#endif // G4DNASCAVENGERMATERIAL_HH
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double GetNumberMoleculePerVolumeUnitForMaterialConf(MolType) const
MaterialMap::iterator end()
int64_t SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule)
G4DNAScavengerMaterial & operator=(const G4DNAScavengerMaterial &)=delete
std::vector< MolType > GetScavengerList() const
std::map< MolType, NbMoleculeInTime > CounterMapType
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
~G4DNAScavengerMaterial() override=default
std::map< G4double, int64_t, G4::MoleculeCounter::TimePrecision > NbMoleculeInTime
void AddAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
MaterialMap::iterator begin()
std::map< MolType, int64_t > MaterialMap
int64_t GetNMoleculesAtTime(MolType molecule, G4double time)
std::vector< MolType > ReactantList
void AddNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
G4DNAScavengerMaterial()=default
G4bool SearchTimeMap(MolType molecule)
void RemoveAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
G4DNAScavengerMaterial(const G4DNAScavengerMaterial &right)=delete
const G4MolecularConfiguration * MolType
G4MolecularConfiguration * GetConfiguration(const G4String &, bool mustExist=true)
static G4MoleculeTable * Instance()