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

#include <G4MuElecInelasticModel.hh>

+ Inheritance diagram for G4MuElecInelasticModel:

Public Member Functions

 G4MuElecInelasticModel (const G4ParticleDefinition *p=0, const G4String &nam="MuElecInelasticModel")
 
virtual ~G4MuElecInelasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 64 of file G4MuElecInelasticModel.hh.

Constructor & Destructor Documentation

◆ G4MuElecInelasticModel()

G4MuElecInelasticModel::G4MuElecInelasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuElecInelasticModel" 
)

Definition at line 60 of file G4MuElecInelasticModel.cc.

62:G4VEmModel(nam),fAtomDeexcitation(0),isInitialised(false)
63{
65
66 verboseLevel= 0;
67 // Verbosity scale:
68 // 0 = nothing
69 // 1 = warning for energy non-conservation
70 // 2 = details of energy budget
71 // 3 = calculation of cross sections, file openings, sampling of atoms
72 // 4 = entering in methods
73
74 if( verboseLevel>0 )
75 {
76 G4cout << "MuElec inelastic model is constructed " << G4endl;
77 }
78
79 //Mark this model as "applicable" for atomic deexcitation
82}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4MuElecInelasticModel()

G4MuElecInelasticModel::~G4MuElecInelasticModel ( )
virtual

Definition at line 86 of file G4MuElecInelasticModel.cc.

87{
88 // Cross section
89
90 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
91 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
92 {
93 G4MuElecCrossSectionDataSet* table = pos->second;
94 delete table;
95 }
96
97 // Final state
98
99 eVecm.clear();
100 pVecm.clear();
101
102}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4MuElecInelasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 256 of file G4MuElecInelasticModel.cc.

261{
262 if (verboseLevel > 3)
263 G4cout << "Calling CrossSectionPerVolume() of G4MuElecInelasticModel" << G4endl;
264
265 G4double density = material->GetTotNbOfAtomsPerVolume();
266
267 /* if (
268 particleDefinition != G4Proton::ProtonDefinition()
269 &&
270 particleDefinition != G4Electron::ElectronDefinition()
271 &&
272 particleDefinition != G4GenericIon::GenericIonDefinition()
273 )
274
275 return 0;*/
276
277 // Calculate total cross section for model
278
279 G4double lowLim = 0;
280 G4double highLim = 0;
281 G4double sigma=0;
282
283 const G4String& particleName = particleDefinition->GetParticleName();
284 G4String nameLocal = particleName ;
285
286 G4double Zeff2 = 1.0;
287 G4double Mion_c2 = particleDefinition->GetPDGMass();
288
289 if (Mion_c2 > proton_mass_c2)
290 {
291 G4ionEffectiveCharge EffCharge ;
292 G4double Zeff = EffCharge.EffectiveCharge(particleDefinition, material,ekin);
293 Zeff2 = Zeff*Zeff;
294
295 if (verboseLevel > 3)
296 G4cout << "Before scaling : " << G4endl
297 << "Particle : " << nameLocal << ", mass : " << Mion_c2/proton_mass_c2 << "*mp, charge " << Zeff
298 << ", Ekin (eV) = " << ekin/eV << G4endl ;
299
300 ekin *= proton_mass_c2/Mion_c2 ;
301 nameLocal = "proton" ;
302
303 if (verboseLevel > 3)
304 G4cout << "After scaling : " << G4endl
305 << "Particle : " << nameLocal << ", Ekin (eV) = " << ekin/eV << G4endl ;
306 }
307
308 if (material == nistSi || material->GetBaseMaterial() == nistSi)
309 {
310
311 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
312 pos1 = lowEnergyLimit.find(nameLocal);
313 if (pos1 != lowEnergyLimit.end())
314 {
315 lowLim = pos1->second;
316 }
317
318 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
319 pos2 = highEnergyLimit.find(nameLocal);
320 if (pos2 != highEnergyLimit.end())
321 {
322 highLim = pos2->second;
323 }
324
325 if (ekin >= lowLim && ekin < highLim)
326 {
327 std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
328 pos = tableData.find(nameLocal);
329
330 if (pos != tableData.end())
331 {
332 G4MuElecCrossSectionDataSet* table = pos->second;
333 if (table != 0)
334 {
335 sigma = table->FindValue(ekin);
336 }
337 }
338 else
339 {
340 G4Exception("G4MuElecInelasticModel::CrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
341 }
342 }
343 else
344 {
345 if (nameLocal!="e-")
346 {
347 // G4cout << "Particle : " << nameLocal << ", Ekin (eV) = " << ekin/eV << G4endl;
348 // G4cout << "### Warning: particle energy out of bounds! ###" << G4endl;
349 }
350 }
351
352 if (verboseLevel > 3)
353 {
354 G4cout << "---> Kinetic energy (eV)=" << ekin/eV << G4endl;
355 G4cout << " - Cross section per Si atom (cm^2)=" << sigma*Zeff2/cm2 << G4endl;
356 G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density*Zeff2/(1./cm) << G4endl;
357 }
358
359 } // if (SiMaterial)
360 return sigma*density*Zeff2;
361
362
363}
@ FatalException
double G4double
Definition: G4Types.hh:64
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:232
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
virtual G4double FindValue(G4double e, G4int componentId=0) const
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ DifferentialCrossSection()

double G4MuElecInelasticModel::DifferentialCrossSection ( G4ParticleDefinition aParticleDefinition,
G4double  k,
G4double  energyTransfer,
G4int  shell 
)

Definition at line 608 of file G4MuElecInelasticModel.cc.

612{
613 G4double sigma = 0.;
614
615 if (energyTransfer >= SiStructure.Energy(LevelIndex))
616 {
617 G4double valueT1 = 0;
618 G4double valueT2 = 0;
619 G4double valueE21 = 0;
620 G4double valueE22 = 0;
621 G4double valueE12 = 0;
622 G4double valueE11 = 0;
623
624 G4double xs11 = 0;
625 G4double xs12 = 0;
626 G4double xs21 = 0;
627 G4double xs22 = 0;
628
629 if (particleDefinition == G4Electron::ElectronDefinition())
630 {
631 // k should be in eV and energy transfer eV also
632
633 std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
634 std::vector<double>::iterator t1 = t2-1;
635 // SI : the following condition avoids situations where energyTransfer >last vector element
636 if (energyTransfer <= eVecm[(*t1)].back() && energyTransfer <= eVecm[(*t2)].back() )
637 {
638 std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), energyTransfer);
639 std::vector<double>::iterator e11 = e12-1;
640
641 std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), energyTransfer);
642 std::vector<double>::iterator e21 = e22-1;
643
644 valueT1 =*t1;
645 valueT2 =*t2;
646 valueE21 =*e21;
647 valueE22 =*e22;
648 valueE12 =*e12;
649 valueE11 =*e11;
650
651 xs11 = eDiffCrossSectionData[LevelIndex][valueT1][valueE11];
652 xs12 = eDiffCrossSectionData[LevelIndex][valueT1][valueE12];
653 xs21 = eDiffCrossSectionData[LevelIndex][valueT2][valueE21];
654 xs22 = eDiffCrossSectionData[LevelIndex][valueT2][valueE22];
655 }
656
657 }
658
659 if (particleDefinition == G4Proton::ProtonDefinition())
660 {
661 // k should be in eV and energy transfer eV also
662 std::vector<double>::iterator t2 = std::upper_bound(pTdummyVec.begin(),pTdummyVec.end(), k);
663 std::vector<double>::iterator t1 = t2-1;
664 if (energyTransfer <= pVecm[(*t1)].back() && energyTransfer <= pVecm[(*t2)].back() )
665 {
666 std::vector<double>::iterator e12 = std::upper_bound(pVecm[(*t1)].begin(),pVecm[(*t1)].end(), energyTransfer);
667 std::vector<double>::iterator e11 = e12-1;
668
669 std::vector<double>::iterator e22 = std::upper_bound(pVecm[(*t2)].begin(),pVecm[(*t2)].end(), energyTransfer);
670 std::vector<double>::iterator e21 = e22-1;
671
672 valueT1 =*t1;
673 valueT2 =*t2;
674 valueE21 =*e21;
675 valueE22 =*e22;
676 valueE12 =*e12;
677 valueE11 =*e11;
678
679 xs11 = pDiffCrossSectionData[LevelIndex][valueT1][valueE11];
680 xs12 = pDiffCrossSectionData[LevelIndex][valueT1][valueE12];
681 xs21 = pDiffCrossSectionData[LevelIndex][valueT2][valueE21];
682 xs22 = pDiffCrossSectionData[LevelIndex][valueT2][valueE22];
683 }
684 }
685
686 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
687 if (xsProduct != 0.)
688 {
689 sigma = QuadInterpolator( valueE11, valueE12,
690 valueE21, valueE22,
691 xs11, xs12,
692 xs21, xs22,
693 valueT1, valueT2,
694 k, energyTransfer);
695 }
696
697 }
698
699 return sigma;
700}
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double Energy(G4int level)
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88

◆ Initialise()

void G4MuElecInelasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 106 of file G4MuElecInelasticModel.cc.

108{
109
110 if (verboseLevel > 3)
111 G4cout << "Calling G4MuElecInelasticModel::Initialise()" << G4endl;
112
113 // Energy limits
114
115 G4String fileElectron("muelec/sigma_inelastic_e_Si");
116 G4String fileProton("muelec/sigma_inelastic_p_Si");
117
120
121 G4String electron;
123
124 G4double scaleFactor = 1e-18 * cm *cm;
125
126 char *path = getenv("G4LEDATA");
127
128 // *** ELECTRON
129 electron = electronDef->GetParticleName();
130
131 tableFile[electron] = fileElectron;
132
133 lowEnergyLimit[electron] = 16.7 * eV;
134 highEnergyLimit[electron] = 100.0 * MeV;
135
136 // Cross section
137
139 tableE->LoadData(fileElectron);
140
141 tableData[electron] = tableE;
142
143 // Final state
144
145 std::ostringstream eFullFileName;
146 eFullFileName << path << "/muelec/sigmadiff_inelastic_e_Si.dat";
147 std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
148
149 if (!eDiffCrossSection)
150 {
151 G4Exception("G4MuElecInelasticModel::Initialise","em0003",FatalException,"Missing data file:/muelec/sigmadiff_inelastic_e_Si.dat");
152 }
153
154 eTdummyVec.push_back(0.);
155 while(!eDiffCrossSection.eof())
156 {
157 double tDummy;
158 double eDummy;
159 eDiffCrossSection>>tDummy>>eDummy;
160 if (tDummy != eTdummyVec.back()) eTdummyVec.push_back(tDummy);
161 for (int j=0; j<6; j++)
162 {
163 eDiffCrossSection>>eDiffCrossSectionData[j][tDummy][eDummy];
164
165 // SI - only if eof is not reached !
166 if (!eDiffCrossSection.eof()) eDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
167
168 eVecm[tDummy].push_back(eDummy);
169
170 }
171 }
172 //
173
174 // *** PROTON
175
176 proton = protonDef->GetParticleName();
177
178 tableFile[proton] = fileProton;
179
180 lowEnergyLimit[proton] = 50. * keV;
181 highEnergyLimit[proton] = 1. * GeV;
182
183 // Cross section
184
186 tableP->LoadData(fileProton);
187
188 tableData[proton] = tableP;
189
190 // Final state
191
192 std::ostringstream pFullFileName;
193 pFullFileName << path << "/muelec/sigmadiff_inelastic_p_Si.dat";
194 std::ifstream pDiffCrossSection(pFullFileName.str().c_str());
195
196 if (!pDiffCrossSection)
197 {
198 G4Exception("G4MuElecInelasticModel::Initialise","em0003",FatalException,"Missing data file:/muelec/sigmadiff_inelastic_p_Si.dat");
199 }
200
201 pTdummyVec.push_back(0.);
202 while(!pDiffCrossSection.eof())
203 {
204 double tDummy;
205 double eDummy;
206 pDiffCrossSection>>tDummy>>eDummy;
207 if (tDummy != pTdummyVec.back()) pTdummyVec.push_back(tDummy);
208 for (int j=0; j<6; j++)
209 {
210 pDiffCrossSection>>pDiffCrossSectionData[j][tDummy][eDummy];
211
212 // SI - only if eof is not reached !
213 if (!pDiffCrossSection.eof()) pDiffCrossSectionData[j][tDummy][eDummy]*=scaleFactor;
214
215 pVecm[tDummy].push_back(eDummy);
216 }
217 }
218
219
220 if (particle==electronDef)
221 {
222 SetLowEnergyLimit(lowEnergyLimit[electron]);
223 SetHighEnergyLimit(highEnergyLimit[electron]);
224 }
225
226 if (particle==protonDef)
227 {
228 SetLowEnergyLimit(lowEnergyLimit[proton]);
229 SetHighEnergyLimit(highEnergyLimit[proton]);
230 }
231
232 if( verboseLevel>0 )
233 {
234 G4cout << "MuElec Inelastic model is initialized " << G4endl
235 << "Energy range: "
236 << LowEnergyLimit() / eV << " eV - "
237 << HighEnergyLimit() / keV << " keV for "
238 << particle->GetParticleName()
239 << " with mass (amu) " << particle->GetPDGMass()/proton_mass_c2
240 << " and charge " << particle->GetPDGCharge()
241 << G4endl << G4endl ;
242 }
243
244 //
245
246 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
247
248 if (isInitialised) { return; }
250 isInitialised = true;
251
252}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
virtual G4bool LoadData(const G4String &argFileName)
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ SampleSecondaries()

void G4MuElecInelasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 367 of file G4MuElecInelasticModel.cc.

372{
373
374 if (verboseLevel > 3)
375 G4cout << "Calling SampleSecondaries() of G4MuElecInelasticModel" << G4endl;
376
377 G4double lowLim = 0;
378 G4double highLim = 0;
379
380 G4double ekin = particle->GetKineticEnergy();
381 G4double k = ekin ;
382
383 G4ParticleDefinition* PartDef = particle->GetDefinition();
384 const G4String& particleName = PartDef->GetParticleName();
385 G4String nameLocal2 = particleName ;
386 G4double particleMass = particle->GetDefinition()->GetPDGMass();
387
388 if (particleMass > proton_mass_c2)
389 {
390 k *= proton_mass_c2/particleMass ;
391 PartDef = G4Proton::ProtonDefinition();
392 nameLocal2 = "proton" ;
393 }
394
395 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
396 pos1 = lowEnergyLimit.find(nameLocal2);
397
398 if (pos1 != lowEnergyLimit.end())
399 {
400 lowLim = pos1->second;
401 }
402
403 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
404 pos2 = highEnergyLimit.find(nameLocal2);
405
406 if (pos2 != highEnergyLimit.end())
407 {
408 highLim = pos2->second;
409 }
410
411 if (k >= lowLim && k < highLim)
412 {
413 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
414 G4double totalEnergy = ekin + particleMass;
415 G4double pSquare = ekin * (totalEnergy + particleMass);
416 G4double totalMomentum = std::sqrt(pSquare);
417
418 G4int Shell = RandomSelect(k,nameLocal2);
419 G4double bindingEnergy = SiStructure.Energy(Shell);
420 if (verboseLevel > 3)
421 {
422 G4cout << "---> Kinetic energy (eV)=" << k/eV << G4endl ;
423 G4cout << "Shell: " << Shell << ", energy: " << bindingEnergy/eV << G4endl;
424 }
425
426 // sample deexcitation
427
428 G4int secNumberInit = 0; // need to know at a certain point the energy of secondaries
429 G4int secNumberFinal = 0; // So I'll make the difference and then sum the energies
430
431 if(fAtomDeexcitation && Shell > 3) {
432 G4int Z = 14;
434
435 if (Shell == 5)
436 {
438 }
439 else if (Shell == 4)
440 {
442 }
443
444 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
445 secNumberInit = fvect->size();
446 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
447 secNumberFinal = fvect->size();
448 }
449
450 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(PartDef,k,Shell);
451
452 if (verboseLevel > 3)
453 {
454 G4cout << "Ionisation process" << G4endl;
455 G4cout << "Shell: " << Shell << " Kin. energy (eV)=" << k/eV
456 << " Sec. energy (eV)=" << secondaryKinetic/eV << G4endl;
457 }
458
459 G4double cosTheta = 0.;
460 G4double phi = 0.;
461 RandomizeEjectedElectronDirection(PartDef, k, secondaryKinetic, cosTheta, phi);
462
463 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
464 G4double dirX = sinTheta*std::cos(phi);
465 G4double dirY = sinTheta*std::sin(phi);
466 G4double dirZ = cosTheta;
467 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
468 deltaDirection.rotateUz(primaryDirection);
469
470 //if (particle->GetDefinition() == G4Electron::ElectronDefinition())
471 //{
472 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
473
474 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
475 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
476 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
477 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
478 finalPx /= finalMomentum;
479 finalPy /= finalMomentum;
480 finalPz /= finalMomentum;
481
482 G4ThreeVector direction;
483 direction.set(finalPx,finalPy,finalPz);
484
486 //}
487 //else fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection) ;
488
489 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
490 G4double deexSecEnergy = 0;
491 for (G4int j=secNumberInit; j < secNumberFinal; j++) {
492 deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();}
493
494 fParticleChangeForGamma->SetProposedKineticEnergy(ekin-bindingEnergy-secondaryKinetic);
495 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy-deexSecEnergy);
496
497 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
498 fvect->push_back(dp);
499
500 }
501
502}
G4AtomicShellEnumerator
int G4int
Definition: G4Types.hh:66
double z() const
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4MuElecInelasticModel::fParticleChangeForGamma
protected

Definition at line 92 of file G4MuElecInelasticModel.hh.

Referenced by G4MuElecInelasticModel(), Initialise(), and SampleSecondaries().


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