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

#include <G4DNAExcitation.hh>

+ Inheritance diagram for G4DNAExcitation:

Public Member Functions

 G4DNAExcitation (const G4String &processName="DNAExcitation", G4ProcessType type=fElectromagnetic)
 
virtual ~G4DNAExcitation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
virtual G4VEmProcessGetEmProcess (const G4String &name)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
 
void SetLambdaBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4VEmModelEmModel (size_t index=0) const
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4int GetNumberOfModels () const
 
G4int GetNumberOfRegionModels (size_t couple_index) const
 
G4VEmModelGetRegionModel (G4int idx, size_t couple_index) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4VEmModelGetCurrentModel () 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 SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 &)
 
- 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 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 Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double kinEnergy, size_t index)
 
virtual 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
 
G4bool IsIntegral () const
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
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
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- 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
G4LossTableManagerlManager
 
G4EmParameterstheParameters
 
G4EmBiasingManagerbiasManager
 
const G4ParticleDefinitiontheGamma
 
const G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGamma fParticleChange
 
std::vector< G4DynamicParticle * > secParticles
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t currentCoupleIndex
 
size_t basedCoupleIndex
 
G4int mainSecondaries
 
G4int secID
 
G4int fluoID
 
G4int augerID
 
G4int biasID
 
G4bool isTheMaster
 
G4double mfpKinEnergy
 
G4double preStepKinEnergy
 
G4double preStepLogKinEnergy
 
G4double preStepLambda
 
- 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 G4DNAExcitation.hh.

Constructor & Destructor Documentation

◆ G4DNAExcitation()

G4DNAExcitation::G4DNAExcitation ( const G4String processName = "DNAExcitation",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 36 of file G4DNAExcitation.cc.

37 :
38 G4VEmProcess(processName, type), isInitialised(false)
39{
41}
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4DNAExcitation()

G4DNAExcitation::~G4DNAExcitation ( )
virtual

Definition at line 45 of file G4DNAExcitation.cc.

46{
47}

Member Function Documentation

◆ InitialiseProcess()

void G4DNAExcitation::InitialiseProcess ( const G4ParticleDefinition p)
protectedvirtual

Implements G4VEmProcess.

Definition at line 67 of file G4DNAExcitation.cc.

68{
69 if(!isInitialised)
70 {
71 isInitialised = true;
72 SetBuildTableFlag(false);
73
75
76 if(name == "e-")
77 {
78
79 // Emfietzoglou model
80 /*
81 if(!EmModel()) SetEmModel(new G4DNAEmfietzoglouExcitationModel);
82 EmModel()->SetLowEnergyLimit(8.23*eV);
83 EmModel()->SetHighEnergyLimit(10*MeV);
84 */
85 // Born model
86 if(!EmModel())
87 {
89 SetEmModel(born);
90 born->SetLowEnergyLimit(9 * eV);
91 born->SetHighEnergyLimit(1 * MeV);
92 }
93 AddEmModel(1, EmModel());
94 }
95
96 else if(name == "e+")
97 {
98 if(!EmModel())
99 {
101 SetEmModel(lepts);
102 lepts->SetLowEnergyLimit(1 * eV);
103 lepts->SetHighEnergyLimit(1 * MeV);
104 }
105 AddEmModel(1, EmModel());
106 }
107
108 else if(name == "proton")
109 {
110 if(!EmModel(0)) // MK: Is this a correct test ?
111 {
114 SetEmModel(miller);
115 miller->SetLowEnergyLimit(10 * eV);
116 miller->SetHighEnergyLimit(500 * keV);
117
119 SetEmModel(born);
120 born->SetLowEnergyLimit(500 * keV);
121 born->SetHighEnergyLimit(100 * MeV);
122 }
123
124 AddEmModel(1, EmModel());
125 if(EmModel(1)) AddEmModel(2, EmModel(1));
126 }
127
128 else if(name == "hydrogen")
129 {
130 if(!EmModel())
131 {
134 SetEmModel(miller);
135 miller->SetLowEnergyLimit(10 * eV);
136 miller->SetHighEnergyLimit(500 * keV);
137 }
138 AddEmModel(1, EmModel());
139 }
140
141 else if(name == "alpha" || name == "alpha+" || name == "helium")
142 {
143 if(!EmModel())
144 {
147 SetEmModel(miller);
148 miller->SetLowEnergyLimit(1 * keV);
149 miller->SetHighEnergyLimit(400 * MeV);
150 }
151 AddEmModel(1, EmModel());
152 }
153 }
154}
G4DNABornExcitationModel1 G4DNABornExcitationModel
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
G4VEmModel * EmModel(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 G4DNAExcitation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 51 of file G4DNAExcitation.cc.

52{
53
56
57 return (&p == G4Electron::Electron() || &p == G4Positron::Positron()
59 || &p == instance->GetIon("hydrogen")
60 || &p == instance->GetIon("alpha++")
61 || &p == instance->GetIon("alpha+")
62 || &p == instance->GetIon("helium"));
63}
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87

◆ PrintInfo()

void G4DNAExcitation::PrintInfo ( )
virtual

Reimplemented from G4VEmProcess.

Definition at line 158 of file G4DNAExcitation.cc.

159{
160 if(EmModel(2))
161 {
162 G4cout << " Total cross sections computed from " << EmModel(1)->GetName()
163 << " and " << EmModel(2)->GetName() << " models" << G4endl;
164 }
165 else
166 {
167 G4cout << " Total cross sections computed from "
168 << EmModel()->GetName()
169 << G4endl;
170 }
171}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
Definition: G4VEmModel.hh:827

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