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

#include <G4LivermoreIonisationModel.hh>

+ Inheritance diagram for G4LivermoreIonisationModel:

Public Member Functions

 G4LivermoreIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &processName="LowEnergyIoni")
 
virtual ~G4LivermoreIonisationModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel ()
 
G4LivermoreIonisationModeloperator= (const G4LivermoreIonisationModel &right)=delete
 
 G4LivermoreIonisationModel (const G4LivermoreIonisationModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForLossfParticleChange
 
- 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
 
size_t currentCoupleIndex = 0
 
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 57 of file G4LivermoreIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreIonisationModel() [1/2]

G4LivermoreIonisationModel::G4LivermoreIonisationModel ( const G4ParticleDefinition p = nullptr,
const G4String processName = "LowEnergyIoni" 
)

Definition at line 76 of file G4LivermoreIonisationModel.cc.

77 :
78 G4VEmModel(nam), fParticleChange(nullptr),
79 crossSectionHandler(nullptr), energySpectrum(nullptr),
80 isInitialised(false)
81{
82 fIntrinsicLowEnergyLimit = 12.*eV;
83 fIntrinsicHighEnergyLimit = 100.0*GeV;
84
85 verboseLevel = 0;
87
88 transitionManager = G4AtomicTransitionManager::Instance();
89}
static G4AtomicTransitionManager * Instance()
G4ParticleChangeForLoss * fParticleChange
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4LivermoreIonisationModel()

G4LivermoreIonisationModel::~G4LivermoreIonisationModel ( )
virtual

Definition at line 93 of file G4LivermoreIonisationModel.cc.

94{
95 delete energySpectrum;
96 delete crossSectionHandler;
97}

◆ G4LivermoreIonisationModel() [2/2]

G4LivermoreIonisationModel::G4LivermoreIonisationModel ( const G4LivermoreIonisationModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LivermoreIonisationModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 174 of file G4LivermoreIonisationModel.cc.

180{
181 G4int iZ = G4int(Z);
182 if (!crossSectionHandler)
183 {
184 G4Exception("G4LivermoreIonisationModel::ComputeCrossSectionPerAtom",
185 "em1007",FatalException,
186 "The cross section handler is not correctly initialized");
187 return 0;
188 }
189
190 //The cut is already included in the crossSectionHandler
191 G4double cs =
192 crossSectionHandler->GetCrossSectionAboveThresholdForElement(energy,
193 cutEnergy,
194 iZ);
195
196 if (verboseLevel > 1)
197 {
198 G4cout << "G4LivermoreIonisationModel " << G4endl;
199 G4cout << "Cross section for delta emission > "
200 << cutEnergy/keV << " keV at "
201 << energy/keV << " keV and Z = " << iZ << " --> "
202 << cs/barn << " barn" << G4endl;
203 }
204 return cs;
205}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetCrossSectionAboveThresholdForElement(G4double energy, G4double cutEnergy, G4int Z)
G4double energy(const ThreeVector &p, const G4double m)

◆ ComputeDEDXPerVolume()

G4double G4LivermoreIonisationModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 211 of file G4LivermoreIonisationModel.cc.

215{
216 G4double sPower = 0.0;
217
218 const G4ElementVector* theElementVector = material->GetElementVector();
219 size_t NumberOfElements = material->GetNumberOfElements() ;
220 const G4double* theAtomicNumDensityVector =
221 material->GetAtomicNumDensityVector();
222
223 // loop for elements in the material
224 for (size_t iel=0; iel<NumberOfElements; iel++ )
225 {
226 G4int iZ = (G4int)((*theElementVector)[iel]->GetZ());
227 G4int nShells = transitionManager->NumberOfShells(iZ);
228 for (G4int n=0; n<nShells; n++)
229 {
230 G4double e = energySpectrum->AverageEnergy(iZ, 0.0,cutEnergy,
231 kineticEnergy, n);
232 G4double cs= crossSectionHandler->FindValue(iZ,kineticEnergy, n);
233 sPower += e * cs * theAtomicNumDensityVector[iel];
234 }
235 G4double esp = energySpectrum->Excitation(iZ,kineticEnergy);
236 sPower += esp * theAtomicNumDensityVector[iel];
237 }
238
239 if (verboseLevel > 2)
240 {
241 G4cout << "G4LivermoreIonisationModel " << G4endl;
242 G4cout << "Stopping power < " << cutEnergy/keV
243 << " keV at " << kineticEnergy/keV << " keV = "
244 << sPower/(keV/mm) << " keV/mm" << G4endl;
245 }
246
247 return sPower;
248}
std::vector< const G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
G4double FindValue(G4int Z, G4double e) const
virtual G4double AverageEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=nullptr) const =0
virtual G4double Excitation(G4int Z, G4double kineticEnergy) const =0

◆ GetVerboseLevel()

G4int G4LivermoreIonisationModel::GetVerboseLevel ( )
inline

Definition at line 86 of file G4LivermoreIonisationModel.hh.

86{return verboseLevel;};

◆ Initialise()

void G4LivermoreIonisationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 101 of file G4LivermoreIonisationModel.cc.

103{
104 //Check that the Livermore Ionisation is NOT attached to e+
105 if (particle != G4Electron::Electron())
106 {
107 G4Exception("G4LivermoreIonisationModel::Initialise",
108 "em0002",FatalException,
109 "Livermore Ionisation Model is applicable only to electrons");
110 }
111 transitionManager->Initialise();
112
113 //Read energy spectrum
114 if (energySpectrum)
115 {
116 delete energySpectrum;
117 energySpectrum = nullptr;
118 }
119 energySpectrum = new G4eIonisationSpectrum();
120 if (verboseLevel > 3)
121 G4cout << "G4VEnergySpectrum is initialized" << G4endl;
122
123 //Initialize cross section handler
124 if (crossSectionHandler)
125 {
126 delete crossSectionHandler;
127 crossSectionHandler = nullptr;
128 }
129
130 const size_t nbins = 20;
131 G4double emin = LowEnergyLimit();
132 G4double emax = HighEnergyLimit();
133 G4int ndec = G4int(std::log10(emax/emin) + 0.5);
134 if(ndec <= 0) { ndec = 1; }
135
136 G4VDataSetAlgorithm* interpolation = new G4SemiLogInterpolation();
137 crossSectionHandler =
138 new G4eIonisationCrossSectionHandler(energySpectrum,interpolation,
139 emin,emax,nbins*ndec);
140 crossSectionHandler->Clear();
141 crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
142 //This is used to retrieve cross section values later on
143 G4VEMDataSet* emdata =
144 crossSectionHandler->BuildMeanFreePathForMaterials(&cuts);
145 //The method BuildMeanFreePathForMaterials() is required here only to force
146 //the building of an internal table: the output pointer can be deleted
147 delete emdata;
148
149 if (verboseLevel > 0)
150 {
151 G4cout << "Livermore Ionisation model is initialized " << G4endl
152 << "Energy range: "
153 << LowEnergyLimit() / keV << " keV - "
154 << HighEnergyLimit() / GeV << " GeV"
155 << G4endl;
156 }
157
158 if (verboseLevel > 3)
159 {
160 G4cout << "Cross section data: " << G4endl;
161 crossSectionHandler->PrintData();
162 G4cout << "Parameters: " << G4endl;
163 energySpectrum->PrintData();
164 }
165
166 if(isInitialised) { return; }
168 isInitialised = true;
169}
void Initialise()
needs to be called once from other code before start of run
static G4Electron * Electron()
Definition: G4Electron.cc:93
void LoadShellData(const G4String &dataFile)
G4VEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=nullptr)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109
virtual void PrintData() const =0

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 252 of file G4LivermoreIonisationModel.cc.

258{
259
260 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
261
262 if (kineticEnergy <= fIntrinsicLowEnergyLimit)
263 {
266 return;
267 }
268
269 // Select atom and shell
270 G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
271 G4int shellIndex = crossSectionHandler->SelectRandomShell(Z, kineticEnergy);
272 const G4AtomicShell* shell = transitionManager->Shell(Z,shellIndex);
274
275 // Sample delta energy using energy interval for delta-electrons
276 G4double energyMax =
277 std::min(maxE,energySpectrum->MaxEnergyOfSecondaries(kineticEnergy));
278 G4double energyDelta = energySpectrum->SampleEnergy(Z, cutE, energyMax,
279 kineticEnergy, shellIndex);
280
281 if (energyDelta == 0.) //nothing happens
282 { return; }
283
285 G4DynamicParticle* delta = new G4DynamicParticle(electron,
286 GetAngularDistribution()->SampleDirectionForShell(aDynamicParticle, energyDelta,
287 Z, shellIndex,
288 couple->GetMaterial()),
289 energyDelta);
290
291 fvect->push_back(delta);
292
293 // Change kinematics of primary particle
294 G4ThreeVector direction = aDynamicParticle->GetMomentumDirection();
295 G4double totalMomentum = std::sqrt(kineticEnergy*(kineticEnergy + 2*electron_mass_c2));
296
297 G4ThreeVector finalP = totalMomentum*direction - delta->GetMomentum();
298 finalP = finalP.unit();
299
300 //This is the amount of energy available for fluorescence
301 G4double theEnergyDeposit = bindingEnergy;
302
303 // fill ParticleChange
304 // changed energy and momentum of the actual particle
305 G4double finalKinEnergy = kineticEnergy - energyDelta - theEnergyDeposit;
306 if(finalKinEnergy < 0.0)
307 {
308 theEnergyDeposit += finalKinEnergy;
309 finalKinEnergy = 0.0;
310 }
311 else
312 {
314 }
316
317 if (theEnergyDeposit < 0)
318 {
319 G4cout << "G4LivermoreIonisationModel: Negative energy deposit: "
320 << theEnergyDeposit/eV << " eV" << G4endl;
321 theEnergyDeposit = 0.0;
322 }
323
324 //Assign local energy deposit
326
327 if (verboseLevel > 1)
328 {
329 G4cout << "-----------------------------------------------------------" << G4endl;
330 G4cout << "Energy balance from G4LivermoreIonisation" << G4endl;
331 G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
332 G4cout << "-----------------------------------------------------------" << G4endl;
333 G4cout << "Outgoing primary energy: " << finalKinEnergy/keV << " keV" << G4endl;
334 G4cout << "Delta ray " << energyDelta/keV << " keV" << G4endl;
335 G4cout << "Fluorescence: " << (bindingEnergy-theEnergyDeposit)/keV << " keV" << G4endl;
336 G4cout << "Local energy deposit " << theEnergyDeposit/keV << " keV" << G4endl;
337 G4cout << "Total final state: " << (finalKinEnergy+energyDelta+bindingEnergy)
338 << " keV" << G4endl;
339 G4cout << "-----------------------------------------------------------" << G4endl;
340 }
341 return;
342}
Hep3Vector unit() const
G4double BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
const G4Material * GetMaterial() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
G4int SelectRandomShell(G4int Z, G4double e) const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
virtual G4double MaxEnergyOfSecondaries(G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=nullptr) const =0
virtual G4double SampleEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=nullptr) const =0
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SetVerboseLevel()

void G4LivermoreIonisationModel::SetVerboseLevel ( G4int  vl)
inline

Definition at line 85 of file G4LivermoreIonisationModel.hh.

85{verboseLevel = vl;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForLoss* G4LivermoreIonisationModel::fParticleChange
protected

Definition at line 92 of file G4LivermoreIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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