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

#include <G4DNABornExcitationModel1.hh>

+ Inheritance diagram for G4DNABornExcitationModel1:

Public Member Functions

 G4DNABornExcitationModel1 (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNABornExcitationModel")
 
 ~G4DNABornExcitationModel1 () override
 
G4DNABornExcitationModel1operator= (const G4DNABornExcitationModel1 &right)=delete
 
 G4DNABornExcitationModel1 (const G4DNABornExcitationModel1 &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 42 of file G4DNABornExcitationModel1.hh.

Constructor & Destructor Documentation

◆ G4DNABornExcitationModel1() [1/2]

G4DNABornExcitationModel1::G4DNABornExcitationModel1 ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNABornExcitationModel" )

Definition at line 40 of file G4DNABornExcitationModel1.cc.

41 :
42G4VEmModel(nam)
43{
44 fpMolWaterDensity = nullptr;
45 fHighEnergy = 0;
46 fLowEnergy = 0;
47 fParticleDefinition = nullptr;
48
49 verboseLevel = 0;
50 // Verbosity scale:
51 // 0 = nothing
52 // 1 = warning for energy non-conservation
53 // 2 = details of energy budget
54 // 3 = calculation of cross sections, file openings, sampling of atoms
55 // 4 = entering in methods
56
57 if (verboseLevel > 0)
58 {
59 G4cout << "Born excitation model is constructed " << G4endl;
60 }
62
63 // Selection of stationary mode
64
65 statCode = false;
66}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNABornExcitationModel1()

G4DNABornExcitationModel1::~G4DNABornExcitationModel1 ( )
override

Definition at line 70 of file G4DNABornExcitationModel1.cc.

71{
72 // Cross section
73
74 delete fTableData;
75}

◆ G4DNABornExcitationModel1() [2/2]

G4DNABornExcitationModel1::G4DNABornExcitationModel1 ( const G4DNABornExcitationModel1 & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNABornExcitationModel1::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 137 of file G4DNABornExcitationModel1.cc.

142{
143 if (verboseLevel > 3)
144 {
145 G4cout << "Calling CrossSectionPerVolume() of G4DNABornExcitationModel1"
146 << G4endl;
147 }
148
149 if(particleDefinition != fParticleDefinition) return 0;
150
151 // Calculate total cross section for model
152
153 G4double sigma=0;
154
155 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
156
157 if (ekin >= fLowEnergy && ekin <= fHighEnergy)
158 {
159 sigma = fTableData->FindValue(ekin);
160 }
161
162 if (verboseLevel > 2)
163 {
164 G4cout << "__________________________________" << G4endl;
165 G4cout << "G4DNABornExcitationModel1 - XS INFO START" << G4endl;
166 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
167 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
168 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
169 G4cout << "G4DNABornExcitationModel1 - XS INFO END" << G4endl;
170 }
171
172 return sigma*waterDensity;
173}
double G4double
Definition G4Types.hh:83
G4double FindValue(G4double e, G4int componentId=0) const override
std::size_t GetIndex() const

◆ GetPartialCrossSection()

G4double G4DNABornExcitationModel1::GetPartialCrossSection ( const G4Material * ,
G4int level,
const G4ParticleDefinition * particle,
G4double kineticEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 214 of file G4DNABornExcitationModel1.cc.

218{
219 if (fParticleDefinition != particle)
220 {
221 G4Exception("G4DNABornExcitationModel1::GetPartialCrossSection",
222 "bornParticleType",
224 "Model initialized for another particle type.");
225 }
226
227 return fTableData->GetComponent(level)->FindValue(kineticEnergy);
228}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
const G4VEMDataSet * GetComponent(G4int componentId) const override
virtual G4double FindValue(G4double x, G4int componentId=0) const =0

◆ Initialise()

void G4DNABornExcitationModel1::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & = *(new G4DataVector()) )
overridevirtual

Implements G4VEmModel.

Definition at line 79 of file G4DNABornExcitationModel1.cc.

81{
82
83 if (verboseLevel > 3)
84 {
85 G4cout << "Calling G4DNABornExcitationModel1::Initialise()" << G4endl;
86 }
87
88 if(fParticleDefinition != nullptr && fParticleDefinition != particle)
89 {
90 G4Exception("G4DNABornExcitationModel1::Initialise","em0001",
91 FatalException,"Model already initialized for another particle type.");
92 }
93
94 fParticleDefinition = particle;
95
96 if(particle->GetParticleName() == "e-")
97 {
98 fTableFile = "dna/sigma_excitation_e_born";
99 fLowEnergy = 9*eV;
100 fHighEnergy = 1*MeV;
101 }
102 else if(particle->GetParticleName() == "proton")
103 {
104 fTableFile = "dna/sigma_excitation_p_born";
105 fLowEnergy = 500. * keV;
106 fHighEnergy = 100. * MeV;
107 }
108
109 SetLowEnergyLimit(fLowEnergy);
110 SetHighEnergyLimit(fHighEnergy);
111
112 G4double scaleFactor = (1.e-22 / 3.343) * m*m;
113 fTableData = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
114 fTableData->LoadData(fTableFile);
115
116 if( verboseLevel>0 )
117 {
118 G4cout << "Born excitation model is initialized " << G4endl
119 << "Energy range: "
120 << LowEnergyLimit() / eV << " eV - "
121 << HighEnergyLimit() / keV << " keV for "
122 << particle->GetParticleName()
123 << G4endl;
124 }
125
126 // Initialize water density pointer
128
129 if (isInitialised)
130 { return;}
132 isInitialised = true;
133}
G4bool LoadData(const G4String &argFileName) override
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ operator=()

G4DNABornExcitationModel1 & G4DNABornExcitationModel1::operator= ( const G4DNABornExcitationModel1 & right)
delete

◆ SampleSecondaries()

void G4DNABornExcitationModel1::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 177 of file G4DNABornExcitationModel1.cc.

182{
183
184 if (verboseLevel > 3)
185 {
186 G4cout << "Calling SampleSecondaries() of G4DNABornExcitationModel1"
187 << G4endl;
188 }
189
190 G4double k = aDynamicParticle->GetKineticEnergy();
191
192 G4int level = RandomSelect(k);
193 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
194 G4double newEnergy = k - excitationEnergy;
195
196 if (newEnergy > 0)
197 {
199
200 if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
202
204 }
205
206 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
208 level,
209 theIncomingTrack);
210}
@ eExcitedMolecule
int G4int
Definition G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNABornExcitationModel1::SelectStationary ( G4bool input)
inline

Definition at line 104 of file G4DNABornExcitationModel1.hh.

105{
106 statCode = input;
107}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNABornExcitationModel1::fParticleChangeForGamma
protected

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