Geant4-11
Data Structures | Public Types | Public Member Functions | Private Attributes
G4DNAScavengerMaterial Class Reference

#include <G4DNAScavengerMaterial.hh>

Inheritance diagram for G4DNAScavengerMaterial:
G4VScavengerMaterial

Data Structures

struct  Search
 

Public Types

using CounterMapType = std::map< MolType, NbMoleculeAgainstTime >
 
using MaterialMap = std::map< MolType, G4double >
 
using MolType = const G4MolecularConfiguration *
 
using ReactantList = std::vector< MolType >
 

Public Member Functions

void AddAMoleculeAtTime (MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
 
void AddNumberMoleculePerVolumeUnitForMaterialConf (MolType, G4double)
 
MaterialMap::iterator begin ()
 
void Dump ()
 
MaterialMap::iterator end ()
 
G4bool find (MolType type)
 
 G4DNAScavengerMaterial ()=default
 
 G4DNAScavengerMaterial (const G4DNAScavengerMaterial &right)=delete
 
 G4DNAScavengerMaterial (G4VChemistryWorld *)
 
G4int GetNMoleculesAtTime (MolType molecule, G4double time)
 
G4double GetNumberMoleculePerVolumeUnitForMaterialConf (MolType) const
 
std::vector< MolTypeGetScavengerList () const
 
void Initialize ()
 
G4DNAScavengerMaterialoperator= (const G4DNAScavengerMaterial &)=delete
 
void PrintInfo ()
 
void ReduceNumberMoleculePerVolumeUnitForMaterialConf (MolType, G4double)
 
void RemoveAMoleculeAtTime (MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
 
void Reset () override
 
G4bool SearchTimeMap (MolType molecule)
 
G4int SearchUpperBoundTime (G4double time, G4bool sameTypeOfMolecule)
 
void SetCounterAgainstTime ()
 
size_t size ()
 
 ~G4DNAScavengerMaterial () override=default
 

Private Attributes

G4bool fCounterAgainstTime
 
CounterMapType fCounterMap
 
G4bool fIsInitialized
 
G4VChemistryWorldfpChemistryInfo
 
std::unique_ptr< SearchfpLastSearch
 
std::map< MolType, G4doublefScavengerTable
 
G4int fVerbose
 

Detailed Description

Definition at line 42 of file G4DNAScavengerMaterial.hh.

Member Typedef Documentation

◆ CounterMapType

Definition at line 48 of file G4DNAScavengerMaterial.hh.

◆ MaterialMap

Definition at line 46 of file G4DNAScavengerMaterial.hh.

◆ MolType

Definition at line 45 of file G4DNAScavengerMaterial.hh.

◆ ReactantList

Definition at line 47 of file G4DNAScavengerMaterial.hh.

Constructor & Destructor Documentation

◆ G4DNAScavengerMaterial() [1/3]

G4DNAScavengerMaterial::G4DNAScavengerMaterial ( )
default

◆ G4DNAScavengerMaterial() [2/3]

G4DNAScavengerMaterial::G4DNAScavengerMaterial ( G4VChemistryWorld pChemistryInfo)
explicit

Definition at line 44 of file G4DNAScavengerMaterial.cc.

References Initialize().

◆ ~G4DNAScavengerMaterial()

G4DNAScavengerMaterial::~G4DNAScavengerMaterial ( )
overridedefault

◆ G4DNAScavengerMaterial() [3/3]

G4DNAScavengerMaterial::G4DNAScavengerMaterial ( const G4DNAScavengerMaterial right)
delete

Member Function Documentation

◆ AddAMoleculeAtTime()

void G4DNAScavengerMaterial::AddAMoleculeAtTime ( MolType  molecule,
G4double  time,
const G4ThreeVector position = nullptr,
G4int  number = 1 
)

Definition at line 226 of file G4DNAScavengerMaterial.cc.

229{
230 if(fVerbose != 0)
231 {
232 G4cout << "G4DNAScavengerMaterial::AddAMoleculeAtTime : "
233 << molecule->GetName() << " at time : " << G4BestUnit(time, "Time")
234 << G4endl;
235 }
236
237 auto counterMap_i = fCounterMap.find(molecule);
238
239 if(counterMap_i == fCounterMap.end())
240 {
241 fCounterMap[molecule][time] = number;
242 }
243 else if(counterMap_i->second.empty())
244 {
245 counterMap_i->second[time] = number;
246 }
247 else
248 {
249 auto end = counterMap_i->second.rbegin();
250
251 if(end->first <= time || fabs(end->first - time) <=
253 {
254 G4double newValue = end->second + number;
255 counterMap_i->second[time] = newValue;
256 if(newValue != (floor)(fScavengerTable[molecule])) // protection
257 {
258 assert(false);
259 }
260 // G4cout<<" AddAMoleculeAtTime : "<<molecule->GetName()<< " :
261 // "<<newValue<<G4endl;
262 }
263 // else
264 // {
265 //
266 // G4ExceptionDescription errMsg;
267 // errMsg << "Time of species "
268 // << molecule->GetName() << " is "
269 // << G4BestUnit(time, "Time") << " while "
270 // << " global time is "
271 // <<" end->first : "<<end->first
272 // //<<
273 // G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(),
274 // "Time")
275 // << G4endl;
276 // G4Exception("G4DNAScavengerMaterial::AddAMoleculeAtTime",
277 // "TIME_DONT_MATCH",
278 // FatalException, errMsg);
279 // }
280 }
281}
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
MaterialMap::iterator end()
std::map< MolType, G4double > fScavengerTable
static G4ThreadLocal double fPrecision

References end(), fCounterMap, G4::MoleculeCounter::TimePrecision::fPrecision, fScavengerTable, fVerbose, G4BestUnit, G4cout, G4endl, and G4MolecularConfiguration::GetName().

Referenced by AddNumberMoleculePerVolumeUnitForMaterialConf().

◆ AddNumberMoleculePerVolumeUnitForMaterialConf()

void G4DNAScavengerMaterial::AddNumberMoleculePerVolumeUnitForMaterialConf ( MolType  matConf,
G4double  time 
)

Definition at line 133 of file G4DNAScavengerMaterial.cc.

135{
136 // no change these molecules
137
138 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf ||
139 G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)") ==
140 matConf || // pH has no change
141 G4MoleculeTable::Instance()->GetConfiguration("OHm(B)") == matConf)
142 {
143 // G4cout<<"moletype : "<<matConf->GetName()<<G4endl;
144 // kobs is already counted these molecule concentrations
145 return;
146 }
147
148 auto it = fScavengerTable.find(matConf);
149 if(it == fScavengerTable.end()) // matConf must be in fScavengerTable
150 {
151 return;
152 }
153 fScavengerTable[matConf]++;
154
156 {
157 AddAMoleculeAtTime(matConf, time);
158 }
159}
void AddAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)
static G4MoleculeTable * Instance()

References AddAMoleculeAtTime(), fCounterAgainstTime, fScavengerTable, and G4MoleculeTable::Instance().

Referenced by G4DNAUpdateSystemModel::CreateMolecule().

◆ begin()

MaterialMap::iterator G4DNAScavengerMaterial::begin ( )
inline

Definition at line 72 of file G4DNAScavengerMaterial.hh.

72{ return fScavengerTable.begin(); }

References fScavengerTable.

◆ Dump()

void G4DNAScavengerMaterial::Dump ( )

Definition at line 355 of file G4DNAScavengerMaterial.cc.

356{
357 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
358 auto V = pConfinedBox->Volume();
359 for(auto it : fCounterMap)
360 {
361 auto pReactant = it.first;
362
363 G4cout << " --- > For " << pReactant->GetName() << G4endl;
364
365 for(auto it2 : it.second)
366 {
367 G4cout << " " << G4BestUnit(it2.first, "Time") << " "
368 << it2.second / (Avogadro * V * 1.0e-6 /*mm3 to L*/) << G4endl;
369 }
370 }
371}
G4double Volume() const
G4DNABoundingBox * GetChemistryBoundary() const
float Avogadro
Definition: hepunit.py:252

References source.hepunit::Avogadro, fCounterMap, fpChemistryInfo, G4BestUnit, G4cout, G4endl, G4VChemistryWorld::GetChemistryBoundary(), and G4DNABoundingBox::Volume().

Referenced by PrintInfo(), and RemoveAMoleculeAtTime().

◆ end()

MaterialMap::iterator G4DNAScavengerMaterial::end ( )
inline

Definition at line 71 of file G4DNAScavengerMaterial.hh.

71{ return fScavengerTable.end(); }

References fScavengerTable.

Referenced by AddAMoleculeAtTime().

◆ find()

G4bool G4DNAScavengerMaterial::find ( MolType  type)
inline

Definition at line 75 of file G4DNAScavengerMaterial.hh.

76 {
77 auto it = fScavengerTable.find(type);
78 if(it != fScavengerTable.end())
79 {
80 return it->second > 0;
81 }
82 else
83 {
84 return false;
85 }
86 }

References fScavengerTable.

Referenced by ReduceNumberMoleculePerVolumeUnitForMaterialConf().

◆ GetNMoleculesAtTime()

int G4DNAScavengerMaterial::GetNMoleculesAtTime ( MolType  molecule,
G4double  time 
)

Definition at line 373 of file G4DNAScavengerMaterial.cc.

374{
376 {
377 G4cout << "fCounterAgainstTime == false" << G4endl;
378 assert(false);
379 }
380
381 G4bool sameTypeOfMolecule = SearchTimeMap(molecule);
382 return SearchUpperBoundTime(time, sameTypeOfMolecule);
383}
bool G4bool
Definition: G4Types.hh:86
G4int SearchUpperBoundTime(G4double time, G4bool sameTypeOfMolecule)
G4bool SearchTimeMap(MolType molecule)

References fCounterAgainstTime, G4cout, G4endl, SearchTimeMap(), and SearchUpperBoundTime().

◆ GetNumberMoleculePerVolumeUnitForMaterialConf()

G4double G4DNAScavengerMaterial::GetNumberMoleculePerVolumeUnitForMaterialConf ( MolType  matConf) const

Definition at line 72 of file G4DNAScavengerMaterial.cc.

74{
75 // no change these molecules
76 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf)
77 {
78 G4ExceptionDescription exceptionDescription;
79 exceptionDescription << "matConf : "<<matConf->GetName();
80 G4Exception("G4DNAScavengerMaterial::GetNumberMoleculePerVolumeUnitForMaterialConf", "G4DNAScavengerMaterial001",
81 FatalErrorInArgument, exceptionDescription);
82 }
83
84 auto iter = fScavengerTable.find(matConf);
85 if(iter == fScavengerTable.end())
86 {
87 // G4cout<<matConf->GetName()<<G4endl;
88 // throw std::runtime_error("this material is not existed");
89 return 0;
90 }
91 else
92 {
93 if(iter->second >= 1)
94 {
95 return (floor)(iter->second);
96 }
97 else
98 {
99 return 0;
100 }
101 }
102}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

References FatalErrorInArgument, fScavengerTable, G4Exception(), G4MolecularConfiguration::GetName(), and G4MoleculeTable::Instance().

Referenced by G4DNAGillespieDirectMethod::FindScavenging(), and G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength().

◆ GetScavengerList()

std::vector< MolType > G4DNAScavengerMaterial::GetScavengerList ( ) const
inline

Definition at line 90 of file G4DNAScavengerMaterial.hh.

91 {
92 std::vector<MolType> output;
93 for(const auto& it : fScavengerTable)
94 {
95 output.push_back(it.first);
96 }
97 return output;
98 }

References fScavengerTable.

◆ Initialize()

void G4DNAScavengerMaterial::Initialize ( )

Definition at line 57 of file G4DNAScavengerMaterial.cc.

58{
60 {
61 return;
62 }
63
64 if(fpChemistryInfo->size() == 0)
65 {
66 G4cout << "G4DNAScavengerMaterial existed but empty" << G4endl;
67 }
68 Reset();
69 fIsInitialized = true;
70}

References fIsInitialized, fpChemistryInfo, G4cout, G4endl, Reset(), and G4VChemistryWorld::size().

Referenced by G4DNAScavengerMaterial().

◆ operator=()

G4DNAScavengerMaterial & G4DNAScavengerMaterial::operator= ( const G4DNAScavengerMaterial )
delete

◆ PrintInfo()

void G4DNAScavengerMaterial::PrintInfo ( )

Definition at line 161 of file G4DNAScavengerMaterial.cc.

162{
163 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
164 auto iter = fpChemistryInfo->begin();
165 G4cout << "**************************************************************"
166 << G4endl;
167 for(; iter != fpChemistryInfo->end(); iter++)
168 {
169 auto containedConf = iter->first;
170 // auto concentration = iter->second;
171 auto concentration =
172 fScavengerTable[containedConf] / (Avogadro * pConfinedBox->Volume());
173 G4cout << "Scavenger:" << containedConf->GetName() << " : "
174 << concentration / 1.0e-6 /*mm3 to L*/ << " (M) with : "
175 << fScavengerTable[containedConf] << " (molecules)"
176 << "in: " << pConfinedBox->Volume() / (um * um * um) << " (um3)"
177 << G4endl;
178 if(fScavengerTable[containedConf] < 1)
179 {
180 G4cout << "!!!!!!!!!!!!! this molecule has less one molecule for "
181 "considered volume"
182 << G4endl;
183 // assert(false);
184 }
185 if(fVerbose != 0)
186 {
187 Dump();
188 }
189 }
190 G4cout << "**************************************************************"
191 << G4endl;
192}
static constexpr double um
Definition: G4SIunits.hh:93
std::map< MolType, double >::iterator begin()
std::map< MolType, double >::iterator end()

References source.hepunit::Avogadro, G4VChemistryWorld::begin(), Dump(), G4VChemistryWorld::end(), fpChemistryInfo, fScavengerTable, fVerbose, G4cout, G4endl, G4VChemistryWorld::GetChemistryBoundary(), and um.

Referenced by Reset().

◆ ReduceNumberMoleculePerVolumeUnitForMaterialConf()

void G4DNAScavengerMaterial::ReduceNumberMoleculePerVolumeUnitForMaterialConf ( MolType  matConf,
G4double  time 
)

Definition at line 104 of file G4DNAScavengerMaterial.cc.

106{
107 // no change these molecules
108 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == matConf ||
109 G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)") ==
110 matConf || // suppose that pH is not changed during simu
111 G4MoleculeTable::Instance()->GetConfiguration("OHm(B)") == matConf)
112 {
113 // G4cout<<"moletype : "<<matConf->GetName()<<G4endl;
114 // kobs is already counted these molecule concentrations
115 return;
116 }
117 if(!find(matConf)) // matConf must greater than 0
118 {
119 return;
120 }
121 fScavengerTable[matConf]--;
122 if(fScavengerTable[matConf] < 0) // projection
123 {
124 assert(false);
125 }
126
128 {
129 RemoveAMoleculeAtTime(matConf, time);
130 }
131}
void RemoveAMoleculeAtTime(MolType, G4double time, const G4ThreeVector *position=nullptr, G4int number=1)

References fCounterAgainstTime, find(), fScavengerTable, G4MoleculeTable::Instance(), and RemoveAMoleculeAtTime().

Referenced by G4DNAUpdateSystemModel::KillMolecule(), and G4DNAScavengerProcess::PostStepDoIt().

◆ RemoveAMoleculeAtTime()

void G4DNAScavengerMaterial::RemoveAMoleculeAtTime ( MolType  pMolecule,
G4double  time,
const G4ThreeVector position = nullptr,
G4int  number = 1 
)

Definition at line 285 of file G4DNAScavengerMaterial.cc.

288{
289 NbMoleculeAgainstTime& nbMolPerTime = fCounterMap[pMolecule];
290
291 if(fVerbose != 0)
292 {
293 auto it_ = nbMolPerTime.rbegin();
294 G4cout << "G4DNAScavengerMaterial::RemoveAMoleculeAtTime : "
295 << pMolecule->GetName() << " at time : " << G4BestUnit(time, "Time")
296
297 << " form : " << it_->second << G4endl;
298 }
299
300 if(nbMolPerTime.empty())
301 {
302 Dump();
303 G4String errMsg = "You are trying to remove molecule " +
304 pMolecule->GetName() +
305 " from the counter while this kind of molecules has not "
306 "been registered yet";
307 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "",
308 FatalErrorInArgument, errMsg);
309
310 return;
311 }
312 else
313 {
314 auto it = nbMolPerTime.rbegin();
315
316 if(it == nbMolPerTime.rend())
317 {
318 it--;
319
320 G4String errMsg = "There was no " + pMolecule->GetName() +
321 " recorded at the time or even before the time asked";
322 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "",
323 FatalErrorInArgument, errMsg);
324 }
325
326 G4double finalN = it->second - number;
327 if(finalN < 0)
328 {
329 Dump();
330
331 G4cout << "fScavengerTable : " << pMolecule->GetName() << " : "
332 << (fScavengerTable[pMolecule]) << G4endl;
333
335 errMsg << "After removal of " << number << " species of "
336 << " " << it->second << " " << pMolecule->GetName()
337 << " the final number at time " << G4BestUnit(time, "Time")
338 << " is less than zero and so not valid." << G4endl;
339 G4cout << " Global time is "
340 << G4BestUnit(G4Scheduler::Instance()->GetGlobalTime(), "Time")
341 << ". Previous selected time is " << G4BestUnit(it->first, "Time")
342 << G4endl;
343 G4Exception("G4DNAScavengerMaterial::RemoveAMoleculeAtTime", "N_INF_0",
344 FatalException, errMsg);
345 }
346 nbMolPerTime[time] = finalN;
347
348 if(finalN != (floor)(fScavengerTable[pMolecule])) // protection
349 {
350 assert(false);
351 }
352 }
353}
@ FatalException
std::map< G4double, G4int, G4::MoleculeCounter::TimePrecision > NbMoleculeAgainstTime
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101

References Dump(), FatalErrorInArgument, FatalException, fCounterMap, fScavengerTable, fVerbose, G4BestUnit, G4cout, G4endl, G4Exception(), G4MolecularConfiguration::GetName(), and G4Scheduler::Instance().

Referenced by ReduceNumberMoleculePerVolumeUnitForMaterialConf().

◆ Reset()

void G4DNAScavengerMaterial::Reset ( )
overridevirtual

Implements G4VScavengerMaterial.

Definition at line 194 of file G4DNAScavengerMaterial.cc.

195{
196 if(fpChemistryInfo == nullptr)
197 {
198 return;
199 }
200
201 if(fpChemistryInfo->size() == 0)
202 {
203 return;
204 }
205
206 fScavengerTable.clear();
207 fCounterMap.clear();
208 fpLastSearch.reset(nullptr);
209
210 auto pConfinedBox = fpChemistryInfo->GetChemistryBoundary();
211 auto iter = fpChemistryInfo->begin();
212 for(; iter != fpChemistryInfo->end(); iter++)
213 {
214 auto containedConf = iter->first;
215 auto concentration = iter->second;
216 fScavengerTable[containedConf] =
217 floor(Avogadro * concentration * pConfinedBox->Volume());
218 fCounterMap[containedConf][1 * picosecond] =
219 floor(Avogadro * concentration * pConfinedBox->Volume());
220 }
221 PrintInfo();
222}
static constexpr double picosecond
Definition: G4SIunits.hh:141
std::unique_ptr< Search > fpLastSearch

References source.hepunit::Avogadro, G4VChemistryWorld::begin(), G4VChemistryWorld::end(), fCounterMap, fpChemistryInfo, fpLastSearch, fScavengerTable, G4VChemistryWorld::GetChemistryBoundary(), picosecond, PrintInfo(), and G4VChemistryWorld::size().

Referenced by Initialize().

◆ SearchTimeMap()

G4bool G4DNAScavengerMaterial::SearchTimeMap ( MolType  molecule)

Definition at line 385 of file G4DNAScavengerMaterial.cc.

386{
387 if(fpLastSearch == nullptr)
388 {
389 fpLastSearch = std::make_unique<Search>();
390 }
391 else
392 {
393 if(fpLastSearch->fLowerBoundSet &&
394 fpLastSearch->fLastMoleculeSearched->first == molecule)
395 {
396 return true;
397 }
398 }
399
400 auto mol_it = fCounterMap.find(molecule);
401 fpLastSearch->fLastMoleculeSearched = mol_it;
402
403 if(mol_it != fCounterMap.end())
404 {
405 fpLastSearch->fLowerBoundTime =
406 fpLastSearch->fLastMoleculeSearched->second.end();
407 fpLastSearch->fLowerBoundSet = true;
408 }
409 else
410 {
411 fpLastSearch->fLowerBoundSet = false;
412 }
413
414 return false;
415}

References fCounterMap, and fpLastSearch.

Referenced by GetNMoleculesAtTime().

◆ SearchUpperBoundTime()

int G4DNAScavengerMaterial::SearchUpperBoundTime ( G4double  time,
G4bool  sameTypeOfMolecule 
)

Definition at line 419 of file G4DNAScavengerMaterial.cc.

421{
422 auto mol_it = fpLastSearch->fLastMoleculeSearched;
423 if(mol_it == fCounterMap.end())
424 {
425 return 0;
426 }
427
428 NbMoleculeAgainstTime& timeMap = mol_it->second;
429 if(timeMap.empty())
430 {
431 return 0;
432 }
433
434 if(sameTypeOfMolecule)
435 {
436 if(fpLastSearch->fLowerBoundSet &&
437 fpLastSearch->fLowerBoundTime != timeMap.end())
438 {
439 if(fpLastSearch->fLowerBoundTime->first < time)
440 {
441 auto upperToLast = fpLastSearch->fLowerBoundTime;
442 upperToLast++;
443
444 if(upperToLast == timeMap.end())
445 {
446 return fpLastSearch->fLowerBoundTime->second;
447 }
448
449 if(upperToLast->first > time)
450 {
451 return fpLastSearch->fLowerBoundTime->second;
452 }
453 }
454 }
455 }
456
457 auto up_time_it = timeMap.upper_bound(time);
458
459 if(up_time_it == timeMap.end())
460 {
461 auto last_time = timeMap.rbegin();
462 return last_time->second;
463 }
464 if(up_time_it == timeMap.begin())
465 {
466 return 0;
467 }
468
469 up_time_it--;
470
471 fpLastSearch->fLowerBoundTime = up_time_it;
472 fpLastSearch->fLowerBoundSet = true;
473
474 return fpLastSearch->fLowerBoundTime->second;
475}

References fCounterMap, and fpLastSearch.

Referenced by GetNMoleculesAtTime().

◆ SetCounterAgainstTime()

void G4DNAScavengerMaterial::SetCounterAgainstTime ( )
inline

Definition at line 88 of file G4DNAScavengerMaterial.hh.

88{ fCounterAgainstTime = true; }

References fCounterAgainstTime.

◆ size()

size_t G4DNAScavengerMaterial::size ( )
inline

Definition at line 73 of file G4DNAScavengerMaterial.hh.

73{ return fScavengerTable.size(); }

References fScavengerTable.

Field Documentation

◆ fCounterAgainstTime

G4bool G4DNAScavengerMaterial::fCounterAgainstTime
private

◆ fCounterMap

CounterMapType G4DNAScavengerMaterial::fCounterMap
private

◆ fIsInitialized

G4bool G4DNAScavengerMaterial::fIsInitialized
private

Definition at line 107 of file G4DNAScavengerMaterial.hh.

Referenced by Initialize().

◆ fpChemistryInfo

G4VChemistryWorld* G4DNAScavengerMaterial::fpChemistryInfo
private

Definition at line 106 of file G4DNAScavengerMaterial.hh.

Referenced by Dump(), Initialize(), PrintInfo(), and Reset().

◆ fpLastSearch

std::unique_ptr<Search> G4DNAScavengerMaterial::fpLastSearch
private

Definition at line 120 of file G4DNAScavengerMaterial.hh.

Referenced by Reset(), SearchTimeMap(), and SearchUpperBoundTime().

◆ fScavengerTable

std::map<MolType, G4double> G4DNAScavengerMaterial::fScavengerTable
private

◆ fVerbose

G4int G4DNAScavengerMaterial::fVerbose
private

Definition at line 111 of file G4DNAScavengerMaterial.hh.

Referenced by AddAMoleculeAtTime(), PrintInfo(), and RemoveAMoleculeAtTime().


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