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

#include <G4DNAIonisation.hh>

Inheritance diagram for G4DNAIonisation:

Public Member Functions

 G4DNAIonisation (const G4String &processName="DNAIonisation", G4ProcessType type=fElectromagnetic)
 ~G4DNAIonisation () override
G4bool IsApplicable (const G4ParticleDefinition &) override
virtual void PrintInfo ()
Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 ~G4VEmProcess () override
void ProcessDescription (std::ostream &outFile) const override
void PreparePhysicsTable (const G4ParticleDefinition &) override
void BuildPhysicsTable (const G4ParticleDefinition &) override
void StartTracking (G4Track *) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
virtual G4VEmProcessGetEmProcess (const G4String &name)
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
G4double GetCrossSection (const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4double MeanFreePath (const G4Track &track)
void SetLambdaBinning (G4int nbins)
void SetMinKinEnergy (G4double e)
void SetMinKinEnergyPrim (G4double e)
void SetMaxKinEnergy (G4double e)
G4PhysicsTableLambdaTable () const
G4PhysicsTableLambdaTablePrim () const
void SetLambdaTable (G4PhysicsTable *)
void SetLambdaTablePrim (G4PhysicsTable *)
std::vector< G4double > * EnergyOfCrossSectionMax () const
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
G4CrossSectionType CrossSectionType () const
void SetCrossSectionType (G4CrossSectionType val)
const G4ParticleDefinitionParticle () const
const G4ParticleDefinitionSecondaryParticle () const
G4VEmModelSelectModelForMaterial (G4double kinEnergy, std::size_t idxCouple) const
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel (G4VEmModel *, G4int index=0)
G4int NumberOfModels () const
G4VEmModelEmModel (std::size_t index=0) const
const G4VEmModelGetCurrentModel () const
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
const G4ElementGetCurrentElement () const
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
G4double CrossSectionBiasingFactor () const
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
void SetEmMasterProcess (const G4VEmProcess *)
void SetBuildTableFlag (G4bool val)
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
G4bool UseBaseMaterial () const
void BuildLambdaTable ()
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 G4VEmProcess (G4VEmProcess &)=delete
G4VEmProcessoperator= (const G4VEmProcess &right)=delete
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
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)
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 ()
virtual void ResetNumberOfInteractionLengthLeft ()
G4double GetNumberOfInteractionLengthLeft () const
G4double GetTotalNumberOfInteractionLengthTraversed () const
G4bool isAtRestDoItIsEnabled () const
G4bool isAlongStepDoItIsEnabled () const
G4bool isPostStepDoItIsEnabled () const
virtual void DumpInfo () 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 Member Functions

void InitialiseProcess (const G4ParticleDefinition *) override
Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
G4VEmModelSelectModel (G4double kinEnergy, std::size_t)
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
void DefineMaterial (const G4MaterialCutsCouple *couple)
G4int LambdaBinning () const
G4double MinKinEnergy () const
G4double MaxKinEnergy () const
G4double PolarAngleLimit () const
G4ParticleChangeForGammaGetParticleChange ()
void SetParticle (const G4ParticleDefinition *p)
void SetSecondaryParticle (const G4ParticleDefinition *p)
std::size_t CurrentMaterialCutsCoupleIndex () const
const G4MaterialCutsCoupleMaterialCutsCouple () const
G4bool ApplyCuts () const
G4double GetGammaEnergyCut ()
G4double GetElectronEnergyCut ()
void SetStartFromNullFlag (G4bool val)
void SetSplineFlag (G4bool val)
const G4ElementGetTargetElement () const
const G4IsotopeGetTargetIsotope () const
G4int DensityIndex (G4int idx) const
G4double DensityFactor (G4int idx) const
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Additional Inherited Members

Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
Protected Attributes inherited from G4VEmProcess
const G4MaterialCutsCouplecurrentCouple = nullptr
const G4MaterialcurrentMaterial = nullptr
G4EmBiasingManagerbiasManager = nullptr
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
G4double mfpKinEnergy = DBL_MAX
G4double preStepKinEnergy = 0.0
G4double preStepLambda = 0.0
G4int mainSecondaries = 1
G4int secID = _EM
G4int fluoID = _Fluorescence
G4int augerID = _AugerElectron
G4int biasID = _EM
G4int tripletID = _TripletElectron
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
std::size_t coupleIdxLambda = 0
std::size_t idxLambda = 0
G4bool isTheMaster = false
G4bool baseMat = false
std::vector< G4DynamicParticle * > secParticles
G4ParticleChangeForGamma fParticleChange
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

Detailed Description

Definition at line 44 of file G4DNAIonisation.hh.

Constructor & Destructor Documentation

◆ G4DNAIonisation()

G4DNAIonisation::G4DNAIonisation ( const G4String & processName = "DNAIonisation",
G4ProcessType type = fElectromagnetic )

Definition at line 40 of file G4DNAIonisation.cc.

41 :
42 G4VEmProcess(processName, type)
43{
45}
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
void SetProcessSubType(G4int)

◆ ~G4DNAIonisation()

G4DNAIonisation::~G4DNAIonisation ( )
overridedefault

Member Function Documentation

◆ InitialiseProcess()

void G4DNAIonisation::InitialiseProcess ( const G4ParticleDefinition * p)
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 74 of file G4DNAIonisation.cc.

75{
76 if(!isInitialised)
77 {
78 isInitialised = true;
79 SetBuildTableFlag(false);
80
81 G4String name = p->GetParticleName();
82
83 if(name == "e-")
84 {
85 if(EmModel() == nullptr)
86 {
87 auto born =
88 new G4DNABornIonisationModel();
89 SetEmModel(born);
90 born->SetLowEnergyLimit(11. * eV);
91 born->SetHighEnergyLimit(1. * MeV);
92 }
93 AddEmModel(1, EmModel());
94 }
95 else if(name == "e+")
96 {
97 if(EmModel() == nullptr)
98 {
99 auto lepts =
100 new G4LEPTSIonisationModel();
101 SetEmModel(lepts);
102 lepts->SetLowEnergyLimit(1. * eV);
103 lepts->SetHighEnergyLimit(1. * MeV);
104 }
105 AddEmModel(1, EmModel());
106 }
107
108 if(name == "proton")
109 {
110 if(EmModel(0) == nullptr) // MK : Is this a reliable test ? VI: it is useful
111 {
112 auto rudd =
113 new G4DNARuddIonisationModel();
114 rudd->SetLowEnergyLimit(0 * eV);
115 rudd->SetHighEnergyLimit(500 * keV);
116 SetEmModel(rudd);
117
118 auto born =
119 new G4DNABornIonisationModel();
120 born->SetLowEnergyLimit(500 * keV);
121 born->SetHighEnergyLimit(100 * MeV);
122 SetEmModel(born);
123 }
124
125 AddEmModel(1, EmModel());
126 if(EmModel(1) != nullptr) AddEmModel(2, EmModel(1));
127 }
128
129 if(name == "hydrogen")
130 {
131 if(EmModel() == nullptr)
132 {
133 auto rudd =
134 new G4DNARuddIonisationModel();
135 SetEmModel(rudd);
136 rudd->SetLowEnergyLimit(0 * eV);
137 rudd->SetHighEnergyLimit(100 * MeV);
138 }
139 AddEmModel(1, EmModel());
140 }
141
142 if(name == "alpha" || name == "alpha+" || name == "helium")
143 {
144 if(EmModel() == nullptr)
145 {
146 auto rudd =
147 new G4DNARuddIonisationModel();
148 SetEmModel(rudd);
149 rudd->SetLowEnergyLimit(0 * keV);
150 rudd->SetHighEnergyLimit(400 * MeV);
151 }
152 AddEmModel(1, EmModel());
153 }
154
155 // Extension to HZE proposed by Z. Francis
156
157 //SEB
158 if(/*name == "carbon" || name == "nitrogen" || name == "oxygen" || name == "iron" ||*/
159 name == "GenericIon")
160 //
161 {
162 if(EmModel() == nullptr)
163 {
164 auto ruddExt =
165 new G4DNARuddIonisationExtendedModel();
166 SetEmModel(ruddExt);
167 ruddExt->SetLowEnergyLimit(0 * keV);
168 //SEB: 1e6*MeV by default - updated in model class
169 //EmModel()->SetHighEnergyLimit(p->GetAtomicMass()*1e6*MeV);
170 ruddExt->SetHighEnergyLimit(1e6 * MeV);
171 }
172 AddEmModel(1, EmModel());
173 }
174 }
175}
const G4String & GetParticleName() const
G4VEmModel * EmModel(std::size_t index=0) const
void SetBuildTableFlag(G4bool val)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
const char * name(G4int ptype)

◆ IsApplicable()

G4bool G4DNAIonisation::IsApplicable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 54 of file G4DNAIonisation.cc.

55{
56 G4DNAGenericIonsManager *instance;
58
59 return (&p == G4Electron::Electron() || &p == G4Positron::Positron()
60 || &p == G4Proton::Proton() || &p == instance->GetIon("hydrogen")
61 || &p == instance->GetIon("alpha++")
62 || &p == instance->GetIon("alpha+")
63 || &p == instance->GetIon("helium")
64 //SEB
65 //|| &p == instance->GetIon("carbon")
66 //|| &p == instance->GetIon("nitrogen")
67 //|| &p == instance->GetIon("oxygen")
68 //|| &p == instance->GetIon("iron")
70}
G4TemplateRNGHelper< G4long > * G4TemplateRNGHelper< G4long >::instance
static G4DNAGenericIonsManager * Instance()
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4GenericIon * GenericIonDefinition()
static G4Positron * Positron()
Definition G4Positron.cc:90
static G4Proton * Proton()
Definition G4Proton.cc:90

◆ PrintInfo()

void G4DNAIonisation::PrintInfo ( )
virtual

Definition at line 179 of file G4DNAIonisation.cc.

180{
181 if(EmModel(1) != nullptr)
182 {
183 G4cout << " Total cross sections computed from " << EmModel(0)->GetName()
184 << " and " << EmModel(1)->GetName() << " models" << G4endl;
185 }
186 else
187 {
188 G4cout << " Total cross sections computed from "
189 << EmModel()->GetName()
190 << G4endl;
191 }
192}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const

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