Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAScavengerProcess Class Reference

#include <G4DNAScavengerProcess.hh>

Inheritance diagram for G4DNAScavengerProcess:

Classes

struct  G4DNAScavengerProcessState

Public Types

using MolType = const G4MolecularConfiguration*
using Data = const G4DNAMolecularReactionData

Public Member Functions

 G4DNAScavengerProcess (const G4String &aName, const G4DNABoundingBox &box, G4ProcessType type=fUserDefined)
 ~G4DNAScavengerProcess () override
 G4DNAScavengerProcess (const G4DNAScavengerProcess &)=delete
G4DNAScavengerProcessoperator= (const G4DNAScavengerProcess &)=delete
void StartTracking (G4Track *) override
void SetReaction (MolType, Data *pData)
void BuildPhysicsTable (const G4ParticleDefinition &) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 ~G4VITProcess () override
 G4VITProcess (const G4VITProcess &other)
G4VITProcessoperator= (const G4VITProcess &other)
G4bool operator== (const G4VITProcess &right) const
G4bool operator!= (const G4VITProcess &right) const
size_t GetProcessID () const
G4shared_ptr< G4ProcessState_LockGetProcessState ()
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
void ResetProcessState ()
void StartTracking (G4Track *) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
G4double GetInteractionTimeLeft ()
void ResetNumberOfInteractionLengthLeft () override
G4bool ProposesTimeStep () const
Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 G4VProcess (const G4VProcess &right)
virtual ~G4VProcess ()
G4VProcessoperator= (const G4VProcess &)=delete
G4bool operator== (const G4VProcess &right) const
G4bool operator!= (const G4VProcess &right) const
G4double GetCurrentInteractionLength () const
void SetPILfactor (G4double value)
G4double GetPILfactor () const
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4bool IsApplicable (const G4ParticleDefinition &)
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
const G4StringGetProcessName () const
G4ProcessType GetProcessType () const
void SetProcessType (G4ProcessType)
G4int GetProcessSubType () const
void SetProcessSubType (G4int)
virtual const G4VProcessGetCreatorProcess () const
virtual void EndTracking ()
virtual void SetProcessManager (const G4ProcessManager *)
virtual const G4ProcessManagerGetProcessManager ()
G4double GetNumberOfInteractionLengthLeft () const
G4double GetTotalNumberOfInteractionLengthTraversed () const
G4bool isAtRestDoItIsEnabled () const
G4bool isAlongStepDoItIsEnabled () const
G4bool isPostStepDoItIsEnabled () const
virtual void DumpInfo () const
virtual void ProcessDescription (std::ostream &outfile) const
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const
virtual void SetMasterProcess (G4VProcess *masterP)
const G4VProcessGetMasterProcess () const
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)

Protected Attributes

G4bool fIsInitialized
G4double fReturnedValue
G4ParticleChange fParticleChange
const G4MolecularConfigurationfpMolecularConfiguration
std::map< MolType, std::map< MolType, Data * > > fConfMap
std::vector< MolTypefpMaterialVector
MolType fpMaterialConf
const G4DNABoundingBoxfpBoundingBox
G4DNAScavengerMaterialfpScavengerMaterial {nullptr}
MolType fH3Op = G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)")
MolType fH2O = G4MoleculeTable::Instance()->GetConfiguration("H2O")
MolType fHOm = G4MoleculeTable::Instance()->GetConfiguration("OHm(B)")
Protected Attributes inherited from G4VITProcess
G4shared_ptr< G4ProcessStatefpState
G4bool fProposesTimeStep
Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
G4VParticleChangepParticleChange = nullptr
G4ParticleChange aParticleChange
G4double theNumberOfInteractionLengthLeft = -1.0
G4double currentInteractionLength = -1.0
G4double theInitialNumberOfInteractionLength = -1.0
G4String theProcessName
G4String thePhysicsTableFileName
G4ProcessType theProcessType = fNotDefined
G4int theProcessSubType = -1
G4double thePILfactor = 1.0
G4int verboseLevel = 0
G4bool enableAtRestDoIt = true
G4bool enableAlongStepDoIt = true
G4bool enablePostStepDoIt = true

Additional Inherited Members

Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
void CreateInfo ()
template<typename T>
T * GetState ()
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
virtual void ClearInteractionTimeLeft ()
virtual void ClearNumberOfInteractionLengthLeft ()
void SetInstantiateProcessState (G4bool flag)
G4bool InstantiateProcessState ()
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Detailed Description

Definition at line 37 of file G4DNAScavengerProcess.hh.

Member Typedef Documentation

◆ Data

◆ MolType

Constructor & Destructor Documentation

◆ G4DNAScavengerProcess() [1/2]

G4DNAScavengerProcess::G4DNAScavengerProcess ( const G4String & aName,
const G4DNABoundingBox & box,
G4ProcessType type = fUserDefined )
explicit

Definition at line 45 of file G4DNAScavengerProcess.cc.

48 : G4VITProcess(aName, type)
49 , fpBoundingBox(&box)
50{
52 enableAtRestDoIt = false;
53 enableAlongStepDoIt = false;
54 enablePostStepDoIt = true;
57 fIsInitialized = false;
59 fpMaterialConf = nullptr;
60 fProposesTimeStep = true;
62 verboseLevel = 0;
63}
const G4DNABoundingBox * fpBoundingBox
const G4MolecularConfiguration * fpMolecularConfiguration
void SetInstantiateProcessState(G4bool flag)
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
G4bool fProposesTimeStep
G4bool enableAtRestDoIt
G4int verboseLevel
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange
#define DBL_MAX
Definition templates.hh:62

Referenced by G4DNAScavengerProcess(), and operator=().

◆ ~G4DNAScavengerProcess()

G4DNAScavengerProcess::~G4DNAScavengerProcess ( )
override

Definition at line 65 of file G4DNAScavengerProcess.cc.

66{
67 for(auto& iter : fConfMap)
68 {
69 for(auto& iter2 : iter.second)
70 {
71 delete iter2.second;
72 }
73 }
74}
std::map< MolType, std::map< MolType, Data * > > fConfMap

◆ G4DNAScavengerProcess() [2/2]

G4DNAScavengerProcess::G4DNAScavengerProcess ( const G4DNAScavengerProcess & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4DNAScavengerProcess::AlongStepDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 80 of file G4DNAScavengerProcess.hh.

81 {
82 return nullptr;
83 }

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4double ,
G4double & ,
G4GPILSelection *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 72 of file G4DNAScavengerProcess.hh.

75 {
76 return -1.0;
77 }

◆ AtRestDoIt()

G4VParticleChange * G4DNAScavengerProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 66 of file G4DNAScavengerProcess.hh.

67 {
68 return nullptr;
69 }

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 60 of file G4DNAScavengerProcess.hh.

62 {
63 return -1.0;
64 }

◆ BuildPhysicsTable()

void G4DNAScavengerProcess::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4DNAScavengerProcess.cc.

82{
83 fpScavengerMaterial = dynamic_cast<G4DNAScavengerMaterial*>(
85 if(fpScavengerMaterial == nullptr)
86 {
87 return;
88 }
89 fIsInitialized = true;
90}
G4DNAScavengerMaterial * fpScavengerMaterial
static G4Scheduler * Instance()
G4VScavengerMaterial * GetScavengerMaterial() const

◆ operator=()

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

◆ PostStepDoIt()

G4VParticleChange * G4DNAScavengerProcess::PostStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 250 of file G4DNAScavengerProcess.cc.

252{
253 G4Molecule* molecule = GetMolecule(track);
254 auto molConf = molecule->GetMolecularConfiguration();
255 std::vector<G4Track*> products;
256#ifdef G4VERBOSE
257 if(verboseLevel > 1)
258 {
259 G4cout << ">>> Beginning of G4DNAScavengerProcess verbose>>> Returned value : " << G4BestUnit(fReturnedValue, "Time")
260 <<"molecule: "<<molConf->GetName()<<G4endl;
261 G4cout<<" selected Mat : "<<fpMaterialConf->GetName()<< G4endl;
262 }
263#endif
264
265 G4double reactionTime = track.GetGlobalTime();
266 auto data = fConfMap[molConf][fpMaterialConf];
267
268 if(data == nullptr)
269 {
270 G4ExceptionDescription exceptionDescription;
271 exceptionDescription
272 << "No reaction data for scavenger reaction between : "<<fpMaterialConf->GetName()
273 <<" + "<<molConf->GetName()<<G4endl;
274 G4Exception("G4DNAScavengerProcess::PostStepDoIt",
275 "G4DNAScavengerProcess0001111", FatalErrorInArgument,
276 exceptionDescription);
277 }
278
279 fpScavengerMaterial->SetEquilibrium(data, track.GetGlobalTime());
280
281 auto nbSecondaries = data->GetNbProducts();
282
283 for(G4int j = 0; j < nbSecondaries; ++j)
284 {
285 auto product = data->GetProduct(j);
286 auto isScavenger = fpScavengerMaterial->find(product);
287 if(isScavenger){
289 AddNumberMoleculePerVolumeUnitForMaterialConf(product,track.GetGlobalTime());
290 continue;
291 }
292 auto pProduct = new G4Molecule(data->GetProduct(j));
293 auto pProductTrack =
294 pProduct->BuildTrack(reactionTime, track.GetPosition());
295 pProductTrack->SetTrackStatus(fAlive);
296 G4ITTrackHolder::Instance()->Push(pProductTrack);
297 if(!G4ChemicalMoleculeFinder::Instance()->IsOctreeUsed()){
298 G4MoleculeFinder::Instance()->Push(pProductTrack);
299 }
300 products.push_back(pProductTrack);
301 }
302
303#ifdef G4VERBOSE
304 if(verboseLevel != 0)
305 {
306 G4cout << "At time : " << std::setw(7) << G4BestUnit(reactionTime, "Time")
307 << " Reaction : " << GetIT(track)->GetName() << " ("
308 << track.GetTrackID() << ") + " << fpMaterialConf->GetName() << " ("
309 << "B"
310 << ") -> ";
311 }
312#endif
313 if(nbSecondaries > 0)
314 {
315 for(G4int i = 0; i < nbSecondaries; ++i)
316 {
317#ifdef G4VERBOSE
318 if((verboseLevel != 0) && i != 0)
319 {
320 G4cout << " + ";
321 }
322
323 if(verboseLevel != 0)
324 {
325 auto product = data->GetProduct(i);
326 auto isScavenger = fpScavengerMaterial->find(product);
327 if(isScavenger)
328 {
329 G4cout<<product->GetName()<<" (B)";
330 }
331 else
332 {
333 auto trackSize = products.size();
334 if(trackSize > 0)
335 {
336 for(G4int it = 0; it < (G4int)trackSize; ++it)
337 {
338 if((verboseLevel != 0) && it != 0)
339 {
340 G4cout << " + ";
341 }
342 G4cout << GetIT(products.at(it))->GetName() << " ("
343 << products.at(it)->GetTrackID() << ")";
344 }
345 }
346 }
347 }
348#endif
349 }
350 }
351 else
352 {
353#ifdef G4VERBOSE
354 if(verboseLevel != 0)
355 {
356 G4cout << "No product";
357 }
358#endif
359 }
360#ifdef G4VERBOSE
361 if(verboseLevel != 0)
362 {
363 G4cout << G4endl;
364 }
365#endif
366
368 fParticleChange.Initialize(track);
369 fParticleChange.ProposeTrackStatus(fStopAndKill);
370 fpScavengerMaterial->ReduceNumberMoleculePerVolumeUnitForMaterialConf(
371 fpMaterialConf, reactionTime);
372 State(fPreviousTimeAtPreStepPoint) = -1;
373
375 || fpMaterialConf == fH2O
376 || fpMaterialConf == fHOm) { // these scavengers are not changed
378 }
379
380 return &fParticleChange;
381}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define State(X)
G4IT * GetIT(const G4Track *track)
Definition G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:75
#define G4BestUnit(a, b)
@ fAlive
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void Push(G4Track *track) override
static G4ITFinder * Instance()
void Push(G4Track *) override
static G4ITTrackHolder * Instance()
virtual const G4String & GetName() const =0
const G4MolecularConfiguration * GetMolecularConfiguration() const
static G4OctreeFinder * Instance()
void SetInteractionStep(G4bool InteractionStep)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 121 of file G4DNAScavengerProcess.cc.

124{
125 G4Molecule* molecule = GetMolecule(track);
126 auto molConf = molecule->GetMolecularConfiguration();
127 // this because process for moleculeDifinition not for configuration
128 // TODO: need change this
129 auto it = fConfMap.find(molConf);
130 if(it == fConfMap.end())
131 {
132 return DBL_MAX;
133 }
134 fpMolecularConfiguration = nullptr;
135 fpMaterialConf = nullptr;
136
137 fpMolecularConfiguration = molConf;
138 auto MaterialMap = it->second;
139
140 G4double r1 = G4UniformRand();
141 std::map<G4double /*Propensity*/, std::pair<MolType, G4double>>
142 ReactionDataMap;
143 G4double alpha0 = 0;
144
145 for(const auto& mat_it : MaterialMap)
146 {
147 auto matConf = mat_it.first;
148 G4double numMol =
149 fpScavengerMaterial->GetNumberMoleculePerVolumeUnitForMaterialConf(
150 matConf);
151 if(numMol == 0 && matConf != fH2O){
152 continue;}
153 auto data = mat_it.second;
154 auto reactionRate = data->GetObservedReactionRateConstant(); //_const
155 G4double propensity =
156 numMol * reactionRate / (fpBoundingBox->Volume() * Avogadro);
157
158 if(fH2O == matConf){
159 auto factor = reactionRate;
160 propensity = factor;
161 }
162
163 if(verboseLevel > 1)
164 {
165 G4cout << " Material of " << matConf->GetName() << " : " << propensity
166 << G4endl;
167 }
168
169 auto reactionData = std::make_pair(mat_it.first, propensity);
170 if(propensity == 0)
171 {
172 continue;
173 }
174 alpha0 += propensity;
175 ReactionDataMap[alpha0] = reactionData;
176 }
177 if(alpha0 == 0)
178 {
179 if(State(fIsInGoodMaterial))
180 {
182 State(fIsInGoodMaterial) = false;
183 }
184 return DBL_MAX;
185 }
186 auto rSelectedIter = ReactionDataMap.upper_bound(r1 * alpha0);
187
188 fpMaterialConf = rSelectedIter->second.first;
189
190
191 auto type = fConfMap[fpMolecularConfiguration][fpMaterialConf]->GetReactionType();
192 if(!fpScavengerMaterial->IsEquilibrium(type))
193 {
194 return DBL_MAX;
195 }
196
197 State(fIsInGoodMaterial) = true;
198 G4double previousTimeStep(-1.);
199
200 if(State(fPreviousTimeAtPreStepPoint) != -1)
201 {
202 previousTimeStep =
203 track.GetGlobalTime() - State(fPreviousTimeAtPreStepPoint);
204 }
205
206 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
207
208 // condition is set to "Not Forced"
209 *pForceCond = NotForced;
210
211 if((previousTimeStep <= 0.0) ||
212 (fpState->theNumberOfInteractionLengthLeft <= 0.0))
213 {
215 }
216 else if(previousTimeStep > 0.0)
217 {
219 }
220
221 fpState->currentInteractionLength = 1 / (rSelectedIter->second.second);
222 G4double value = DBL_MAX;
223 if(fpState->currentInteractionLength < DBL_MAX)
224 {
225 value = fpState->theNumberOfInteractionLengthLeft *
226 (fpState->currentInteractionLength);
227 }
228#ifdef G4VERBOSE
229 if(verboseLevel > 2)
230 {
231 G4cout << "Material : " << fpMaterialConf->GetName()
232 << " ID: " << track.GetTrackID()
233 << " Track Time : " << track.GetGlobalTime()
234 << " name : " << molConf->GetName()
235 << " Track Position : " << track.GetPosition()
236 << " Returned time : " << G4BestUnit(value, "Time") << G4endl;
237 }
238#endif
239
240 if(value < fReturnedValue)
241 {
242 fReturnedValue = value;
243 }
244
245 return value * -1;
246 // multiple by -1 to indicate to the tracking system that we are returning a
247 // time
248}
@ NotForced
#define G4UniformRand()
Definition Randomize.hh:52
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
void ResetNumberOfInteractionLengthLeft() override

◆ SetReaction()

void G4DNAScavengerProcess::SetReaction ( MolType molConf,
Data * pData )

Definition at line 99 of file G4DNAScavengerProcess.cc.

100{
102 {
103 G4ExceptionDescription exceptionDescription;
104 exceptionDescription
105 << "G4DNASecondOrderReaction was already initialised. ";
106 exceptionDescription << "You cannot set a reaction after initialisation.";
107 G4Exception("G4DNASecondOrderReaction::SetReaction",
108 "G4DNASecondOrderReaction001", FatalErrorInArgument,
109 exceptionDescription);
110 }
111 auto materialConf = pData->GetReactant1() == molConf ? pData->GetReactant2()
112 : pData->GetReactant1();
113 if(verboseLevel > 0)
114 {
115 G4cout << "G4DNAScavengerProcess::SetReaction : " << molConf->GetName()
116 << " materialConf : " << materialConf->GetName() << G4endl;
117 }
118 fConfMap[molConf][materialConf] = pData;
119}

◆ StartTracking()

void G4DNAScavengerProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 92 of file G4DNAScavengerProcess.cc.

93{
95 G4VITProcess::fpState = std::make_shared<G4DNAScavengerProcessState>();
97}
void StartTracking(G4Track *) override
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

Member Data Documentation

◆ fConfMap

std::map<MolType , std::map<MolType , Data*> > G4DNAScavengerProcess::fConfMap
protected

◆ fH2O

MolType G4DNAScavengerProcess::fH2O = G4MoleculeTable::Instance()->GetConfiguration("H2O")
protected

Definition at line 106 of file G4DNAScavengerProcess.hh.

Referenced by PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ fH3Op

MolType G4DNAScavengerProcess::fH3Op = G4MoleculeTable::Instance()->GetConfiguration("H3Op(B)")
protected

Definition at line 105 of file G4DNAScavengerProcess.hh.

Referenced by PostStepDoIt().

◆ fHOm

MolType G4DNAScavengerProcess::fHOm = G4MoleculeTable::Instance()->GetConfiguration("OHm(B)")
protected

Definition at line 107 of file G4DNAScavengerProcess.hh.

Referenced by PostStepDoIt().

◆ fIsInitialized

G4bool G4DNAScavengerProcess::fIsInitialized
protected

Definition at line 95 of file G4DNAScavengerProcess.hh.

Referenced by BuildPhysicsTable(), G4DNAScavengerProcess(), and SetReaction().

◆ fParticleChange

G4ParticleChange G4DNAScavengerProcess::fParticleChange
protected

Definition at line 97 of file G4DNAScavengerProcess.hh.

Referenced by G4DNAScavengerProcess(), and PostStepDoIt().

◆ fpBoundingBox

const G4DNABoundingBox* G4DNAScavengerProcess::fpBoundingBox
protected

◆ fpMaterialConf

MolType G4DNAScavengerProcess::fpMaterialConf
protected

◆ fpMaterialVector

std::vector<MolType> G4DNAScavengerProcess::fpMaterialVector
protected

Definition at line 101 of file G4DNAScavengerProcess.hh.

◆ fpMolecularConfiguration

const G4MolecularConfiguration* G4DNAScavengerProcess::fpMolecularConfiguration
protected

◆ fpScavengerMaterial

G4DNAScavengerMaterial* G4DNAScavengerProcess::fpScavengerMaterial {nullptr}
protected

Definition at line 104 of file G4DNAScavengerProcess.hh.

104{nullptr};

Referenced by BuildPhysicsTable(), PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ fReturnedValue

G4double G4DNAScavengerProcess::fReturnedValue
protected

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