Geant4 11.1.1
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)
 
virtual ~G4VITProcess ()
 
 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 ()
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 BuildPhysicsTable (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 StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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
 
- 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 35 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 , fpScavengerMaterial(nullptr)
51{
53 enableAtRestDoIt = false;
54 enableAlongStepDoIt = false;
55 enablePostStepDoIt = true;
58 fIsInitialized = false;
60 fpMaterialConf = nullptr;
61 fProposesTimeStep = true;
63 verboseLevel = 0;
64}
const G4DNABoundingBox * fpBoundingBox
const G4MolecularConfiguration * fpMolecularConfiguration
G4DNAScavengerMaterial * fpScavengerMaterial
G4ParticleChange fParticleChange
void SetInstantiateProcessState(G4bool flag)
G4bool fProposesTimeStep
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4int verboseLevel
Definition: G4VProcess.hh:360
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4DNAScavengerProcess()

G4DNAScavengerProcess::~G4DNAScavengerProcess ( )
override

Definition at line 66 of file G4DNAScavengerProcess.cc.

67{
68 for(auto& iter : fConfMap)
69 {
70 for(auto& iter2 : iter.second)
71 {
72 delete iter2.second;
73 }
74 }
75}
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 78 of file G4DNAScavengerProcess.hh.

79 {
80 return nullptr;
81 }

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 70 of file G4DNAScavengerProcess.hh.

73 {
74 return -1.0;
75 }

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 64 of file G4DNAScavengerProcess.hh.

65 {
66 return nullptr;
67 }

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 58 of file G4DNAScavengerProcess.hh.

60 {
61 return -1.0;
62 }

◆ BuildPhysicsTable()

void G4DNAScavengerProcess::BuildPhysicsTable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VITProcess.

Definition at line 83 of file G4DNAScavengerProcess.cc.

84{
87 if(fpScavengerMaterial == nullptr)
88 {
89 return;
90 }
91 fIsInitialized = true;
92}
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 249 of file G4DNAScavengerProcess.cc.

251{
252 G4Molecule* molecule = GetMolecule(track);
253 auto molConf = molecule->GetMolecularConfiguration();
254 if(fpMolecularConfiguration != molConf)
255 {
258 State(fPreviousTimeAtPreStepPoint) = -1;
259 return &fParticleChange;
260 }
261 std::vector<G4Track*> products;
262#ifdef G4VERBOSE
263 if(verboseLevel > 1)
264 {
265 G4cout << "___________" << G4endl;
266 G4cout << ">>> Beginning of G4DNAScavengerProcess verbose" << G4endl;
267 G4cout << ">>> Returned value : " << G4BestUnit(fReturnedValue, "Time")
268 << G4endl;
269 G4cout << ">>> Time Step : "
270 << G4BestUnit(G4VScheduler::Instance()->GetTimeStep(), "Time")
271 << G4endl;
272 G4cout << ">>> Global Time : "
273 << G4BestUnit(G4VScheduler::Instance()->GetGlobalTime(), "Time")
274 << G4endl;
275 G4cout << ">>> Global Time Track : "
276 << G4BestUnit(track.GetGlobalTime(), "Time") << G4endl;
277 G4cout << ">>> Track Position : " << track.GetPosition() << G4endl;
278 G4cout << ">>> Reaction : " << molecule->GetName() << "("
279 << track.GetTrackID() << ") + " << fpMaterialConf->GetName()
280 << G4endl;
281 G4cout << ">>> End of G4DNAScavengerProcess verbose <<<" << G4endl;
282 }
283#endif
284
285 G4double reactionTime = track.GetGlobalTime();
287
288 auto nbSecondaries = data->GetNbProducts();
289
290 for(G4int j = 0; j < nbSecondaries; ++j)
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 G4MoleculeFinder::Instance()->Push(pProductTrack);
298 products.push_back(pProductTrack);
299 }
300
301#ifdef G4VERBOSE
302 if(verboseLevel != 0)
303 {
304 G4cout << "At time : " << std::setw(7) << G4BestUnit(reactionTime, "Time")
305 << " Reaction : " << GetIT(track)->GetName() << " ("
306 << track.GetTrackID() << ") + " << fpMaterialConf->GetName() << " ("
307 << "B"
308 << ") -> ";
309 }
310#endif
311 if(nbSecondaries > 0)
312 {
313 for(G4int i = 0; i < nbSecondaries; ++i)
314 {
315#ifdef G4VERBOSE
316 if((verboseLevel != 0) && i != 0)
317 {
318 G4cout << " + ";
319 }
320
321 if(verboseLevel != 0)
322 {
323 G4cout << GetIT(products.at(i))->GetName() << " ("
324 << products.at(i)->GetTrackID() << ")";
325 }
326#endif
327 }
328 }
329 else
330 {
331#ifdef G4VERBOSE
332 if(verboseLevel != 0)
333 {
334 G4cout << "No product";
335 }
336#endif
337 }
338#ifdef G4VERBOSE
339 if(verboseLevel != 0)
340 {
341 G4cout << G4endl;
342 }
343#endif
344
349 fpMaterialConf, reactionTime);
350 State(fPreviousTimeAtPreStepPoint) = -1;
351 return &fParticleChange;
352}
#define State(X)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
#define G4BestUnit(a, b)
@ fAlive
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
void Push(G4Track *track) override
static G4ITFinder * Instance()
virtual void Push(G4Track *)
static G4ITTrackHolder * Instance()
virtual const G4String & GetName() const =0
const G4String & GetName() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:530
const G4String & GetName() const
Definition: G4Molecule.cc:336
void Initialize(const G4Track &) override
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
void ProposeTrackStatus(G4TrackStatus status)
static G4VScheduler * Instance()
Definition: G4VScheduler.cc:48

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 123 of file G4DNAScavengerProcess.cc.

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

◆ SetReaction()

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

Definition at line 101 of file G4DNAScavengerProcess.cc.

102{
104 {
105 G4ExceptionDescription exceptionDescription;
106 exceptionDescription
107 << "G4DNASecondOrderReaction was already initialised. ";
108 exceptionDescription << "You cannot set a reaction after initialisation.";
109 G4Exception("G4DNASecondOrderReaction::SetReaction",
110 "G4DNASecondOrderReaction001", FatalErrorInArgument,
111 exceptionDescription);
112 }
113 auto materialConf = pData->GetReactant1() == molConf ? pData->GetReactant2()
114 : pData->GetReactant1();
115 if(verboseLevel > 0)
116 {
117 G4cout << "G4DNAScavengerProcess::SetReaction : " << molConf->GetName()
118 << " materialConf : " << materialConf->GetName() << G4endl;
119 }
120 fConfMap[molConf][materialConf] = pData;
121}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ StartTracking()

void G4DNAScavengerProcess::StartTracking ( G4Track track)
overridevirtual

Reimplemented from G4VITProcess.

Definition at line 94 of file G4DNAScavengerProcess.cc.

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

Member Data Documentation

◆ fConfMap

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

◆ fIsInitialized

G4bool G4DNAScavengerProcess::fIsInitialized
protected

Definition at line 93 of file G4DNAScavengerProcess.hh.

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

◆ fParticleChange

G4ParticleChange G4DNAScavengerProcess::fParticleChange
protected

Definition at line 95 of file G4DNAScavengerProcess.hh.

Referenced by G4DNAScavengerProcess(), and PostStepDoIt().

◆ fpBoundingBox

const G4DNABoundingBox* G4DNAScavengerProcess::fpBoundingBox
protected

Definition at line 101 of file G4DNAScavengerProcess.hh.

Referenced by PostStepGetPhysicalInteractionLength().

◆ fpMaterialConf

MolType G4DNAScavengerProcess::fpMaterialConf
protected

◆ fpMaterialVector

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

Definition at line 99 of file G4DNAScavengerProcess.hh.

◆ fpMolecularConfiguration

const G4MolecularConfiguration* G4DNAScavengerProcess::fpMolecularConfiguration
protected

◆ fpScavengerMaterial

G4DNAScavengerMaterial* G4DNAScavengerProcess::fpScavengerMaterial
protected

◆ fReturnedValue

G4double G4DNAScavengerProcess::fReturnedValue
protected

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