Geant4 11.4.0
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>
35#include "G4MoleculeTable.hh"
36#include "G4ChemEquilibrium.hh"
37
38class G4Material;
41
43{
44 public:
46 std::map<G4double, int64_t, G4::MoleculeCounter::FixedTimeComparer>;
48 using MaterialMap = std::map<MolType, int64_t>;
49 using ReactantList = std::vector<MolType>;
50 using CounterMapType = std::map<MolType, NbMoleculeInTime>;
53 ~G4DNAScavengerMaterial() override = default;
56 void Initialize();
57
61
63 const G4ThreeVector* position = nullptr,
64 G4int number = 1);
66 const G4ThreeVector* position = nullptr,
67 G4int number = 1);
68
69 void Reset() override;
70
71 void PrintInfo();
72
73 MaterialMap::iterator end() { return fScavengerTable.end(); }
74 MaterialMap::iterator begin() { return fScavengerTable.begin(); }
75 size_t size() { return fScavengerTable.size(); }
76
78 {
79 auto it = fScavengerTable.find(type);
80 if(it != fScavengerTable.end())
81 {
82 return it->second > 0;
83 }
84 return false;
85 }
86
87 void SetCounterAgainstTime() { fCounterAgainstTime = true; }
88 void SetpH(const G4int& ph);
90
91 std::vector<MolType> GetScavengerList() const
92 {
93 std::vector<MolType> output;
94 for(const auto& it : fScavengerTable)
95 {
96 output.push_back(it.first);
97 }
98 return output;
99 }
100
101 void Dump();
102 int64_t GetNMoleculesAtTime(MolType molecule, G4double time);
103 G4bool SearchTimeMap(MolType molecule);
104 int64_t SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule);
105 void ResetEquilibrium();
107 G4double time);
108 G4bool IsEquilibrium(const G4int& reactionType) const;
109
110private:
111 G4VChemistryWorld* fpChemistryInfo = nullptr;
112 G4bool fIsInitialized;
113 MaterialMap fScavengerTable;
114 CounterMapType fCounterMap;
115 G4bool fCounterAgainstTime;
116 G4int fVerbose;
120 struct Search
121 {
122 Search() { fLowerBoundSet = false; }
123 CounterMapType::iterator fLastMoleculeSearched;
124 NbMoleculeInTime::iterator fLowerBoundTime;
125 G4bool fLowerBoundSet;
126 };
127
128 std::unique_ptr<Search> fpLastSearch;
129 void WaterEquilibrium();
130
131 std::map<G4int,std::unique_ptr<G4ChemEquilibrium>> fEquilibriumProcesses;
132};
133#endif // G4DNASCAVENGERMATERIAL_HH
CLHEP::Hep3Vector G4ThreeVector
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()
G4bool IsEquilibrium(const G4int &reactionType) const
int64_t SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule)
G4DNAScavengerMaterial & operator=(const G4DNAScavengerMaterial &)=delete
G4bool SetEquilibrium(const G4DNAMolecularReactionData *pReaction, G4double time)
std::vector< MolType > GetScavengerList() const
std::map< G4double, int64_t, G4::MoleculeCounter::FixedTimeComparer > NbMoleculeInTime
std::map< MolType, NbMoleculeInTime > CounterMapType
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
~G4DNAScavengerMaterial() override=default
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()
G4VScavengerMaterial()=default