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

#include <G4eeToHadronsModel.hh>

+ Inheritance diagram for G4eeToHadronsModel:

Public Member Functions

 G4eeToHadronsModel (G4Vee2hadrons *, G4int ver=0, const G4String &nam="eeToHadrons")
 
virtual ~G4eeToHadronsModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
 
G4DynamicParticleGenerateCMPhoton (G4double)
 
G4double PeakEnergy () const
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 58 of file G4eeToHadronsModel.hh.

Constructor & Destructor Documentation

◆ G4eeToHadronsModel()

G4eeToHadronsModel::G4eeToHadronsModel ( G4Vee2hadrons mod,
G4int  ver = 0,
const G4String nam = "eeToHadrons" 
)
explicit

Definition at line 69 of file G4eeToHadronsModel.cc.

71 : G4VEmModel(nam),
72 model(mod),
73 crossPerElectron(0),
74 crossBornPerElectron(0),
75 isInitialised(false),
76 nbins(100),
77 verbose(ver)
78{
79 theGamma = G4Gamma::Gamma();
80 highKinEnergy = HighEnergyLimit();
81 lowKinEnergy = LowEnergyLimit();
82 emin = lowKinEnergy;
83 emax = highKinEnergy;
84 peakKinEnergy = highKinEnergy;
85 epeak = emax;
86 //verbose = 1;
87}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ ~G4eeToHadronsModel()

G4eeToHadronsModel::~G4eeToHadronsModel ( )
virtual

Definition at line 91 of file G4eeToHadronsModel.cc.

92{
93 delete model;
94}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4eeToHadronsModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 161 of file G4eeToHadronsModel.cc.

166{
167 return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
168}
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)

◆ ComputeCrossSectionPerElectron()

G4double G4eeToHadronsModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Definition at line 172 of file G4eeToHadronsModel.cc.

176{
177 return (crossPerElectron) ? crossPerElectron->Value(energy) : 0.0;
178}
G4double Value(G4double theEnergy, std::size_t &lastidx) const

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ CrossSectionPerVolume()

G4double G4eeToHadronsModel::CrossSectionPerVolume ( const G4Material mat,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 149 of file G4eeToHadronsModel.cc.

154{
155 return mat->GetElectronDensity()*
156 ComputeCrossSectionPerElectron(p, kineticEnergy);
157}
G4double GetElectronDensity() const
Definition: G4Material.hh:215

◆ GenerateCMPhoton()

G4DynamicParticle * G4eeToHadronsModel::GenerateCMPhoton ( G4double  e)

Definition at line 274 of file G4eeToHadronsModel.cc.

275{
276 G4double x;
277 G4DynamicParticle* gamma = nullptr;
278 G4double LL = 2.0*G4Log(e/electron_mass_c2);
279 G4double bt = 2.0*fine_structure_const*(LL - 1.)/pi;
280 G4double btm1= bt - 1.0;
281 G4double del = 1. + fine_structure_const*(1.5*LL + pi*pi/3. -2.)/pi;
282
283 G4double s0 = crossBornPerElectron->Value(e);
284 G4double de = (emax - emin)/(G4double)nbins;
285 G4double xmax = 0.5*(1.0 - (emin*emin)/(e*e));
286 G4double xmin = std::min(de/e, xmax);
287 G4double ds = s0*(del*G4Exp(G4Log(xmin)*bt) - bt*(xmin - 0.25*xmin*xmin));
288 G4double e1 = e*(1. - xmin);
289
290 //G4cout << "e1= " << e1 << G4endl;
291 if(e1 < emax && s0*G4UniformRand()<ds) {
292 x = xmin*G4Exp(G4Log(G4UniformRand())/bt);
293 } else {
294
295 x = xmin;
296 G4double s1 = crossBornPerElectron->Value(e1);
297 G4double w1 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
298 G4double grej = s1*w1;
299 G4double f;
300 /*
301 G4cout << "e(GeV)= " << e/GeV << " epeak(GeV)= " << epeak/GeV
302 << " s1= " << s1 << " w1= " << w1
303 << " grej= " << grej << G4endl;
304 */
305 // Above emax cross section is const
306 if(e1 > emax) {
307 x = 0.5*(1. - (emax*emax)/(e*e));
308 G4double s2 = crossBornPerElectron->Value(emax);
309 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
310 grej = s2*w2;
311 //G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2
312 // << " grej= " << grej << G4endl;
313 }
314
315 if(e1 > epeak) {
316 x = 0.5*(1.0 - (epeak*epeak)/(e*e));
317 G4double s2 = crossBornPerElectron->Value(epeak);
318 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
319 grej = std::max(grej,s2*w2);
320 //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2
321 // << " grej= " << grej << G4endl;
322 }
323 G4int ii = 0;
324 const G4int iimax = 1000;
325 do {
326 x = xmin + G4UniformRand()*(xmax - xmin);
327
328 G4double s2 = crossBornPerElectron->Value(sqrt(1.0 - 2*x)*e);
329 G4double w2 = bt*(del*G4Exp(G4Log(x)*btm1) - 1.0 + 0.5*x);
330 /*
331 G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
332 << " s2= " << s2 << " w2= " << w2 << G4endl;
333 */
334 f = s2*w2;
335 if(f > grej) {
336 G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
337 << f << " > " << grej << " majorant is`small!"
338 << G4endl;
339 }
340 if(++ii >= iimax) { break; }
341 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
342 } while (f < grej*G4UniformRand());
343 }
344
345 G4ThreeVector dir(0.0,0.0,1.0);
346 if(G4UniformRand() > 0.5) { dir.set(0.0,0.0,-1.0); }
347 //G4cout << "Egamma(MeV)= " << x*e << " " << dir << G4endl;
348 gamma = new G4DynamicParticle(theGamma,dir,x*e);
349 return gamma;
350}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4double pi

Referenced by SampleSecondaries().

◆ Initialise()

void G4eeToHadronsModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 98 of file G4eeToHadronsModel.cc.

100{
101 if(isInitialised) { return; }
102 isInitialised = true;
103
104 // CM system
105 emin = model->LowEnergy();
106 emax = model->HighEnergy();
107
108 // peak energy
109 epeak = std::min(model->PeakEnergy(), emax);
110
111 if(verbose>0) {
112 G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
113 G4cout << "CM System: emin(MeV)= " << emin/MeV
114 << " epeak(MeV)= " << epeak/MeV
115 << " emax(MeV)= " << emax/MeV
116 << G4endl;
117 }
118
119 crossBornPerElectron = model->PhysicsVector();
120 crossPerElectron = model->PhysicsVector();
121 nbins = crossPerElectron->GetVectorLength();
122 for(G4int i=0; i<nbins; ++i) {
123 G4double e = crossPerElectron->Energy(i);
124 G4double cs = model->ComputeCrossSection(e);
125 crossBornPerElectron->PutValue(i, cs);
126 }
127 ComputeCMCrossSectionPerElectron();
128
129 if(verbose>1) {
130 G4cout << "G4eeToHadronsModel: Cross sections per electron"
131 << " nbins= " << nbins
132 << " emin(MeV)= " << emin/MeV
133 << " emax(MeV)= " << emax/MeV
134 << G4endl;
135 for(G4int i=0; i<nbins; ++i) {
136 G4double e = crossPerElectron->Energy(i);
137 G4double s1 = crossPerElectron->Value(e);
138 G4double s2 = crossBornPerElectron->Value(e);
139 G4cout << "E(MeV)= " << e/MeV
140 << " cross(nb)= " << s1/nanobarn
141 << " crossBorn(nb)= " << s2/nanobarn
142 << G4endl;
143 }
144 }
145}
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
std::size_t GetVectorLength() const
G4PhysicsVector * PhysicsVector() const
G4double LowEnergy() const
virtual G4double ComputeCrossSection(G4double) const =0
virtual G4double PeakEnergy() const =0
G4double HighEnergy() const

◆ PeakEnergy()

G4double G4eeToHadronsModel::PeakEnergy ( ) const
inline

Definition at line 127 of file G4eeToHadronsModel.hh.

128{
129 return peakKinEnergy;
130}

◆ SampleSecondaries()

void G4eeToHadronsModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  newp,
const G4MaterialCutsCouple ,
const G4DynamicParticle dParticle,
G4double  tmin = 0.0,
G4double  maxEnergy = DBL_MAX 
)
overridevirtual

Implements G4VEmModel.

Definition at line 182 of file G4eeToHadronsModel.cc.

187{
188 if(crossPerElectron) {
189 G4double t = dParticle->GetKineticEnergy() + 2*electron_mass_c2;
190 G4LorentzVector inlv = dParticle->Get4Momentum() +
191 G4LorentzVector(0.0,0.0,0.0,electron_mass_c2);
192 G4double e = inlv.m();
193 G4ThreeVector inBoost = inlv.boostVector();
194 //G4cout << "G4eeToHadronsModel::SampleSecondaries e= " << e
195 // << " " << inlv << " " << inBoost <<G4endl;
196 if(e > emin) {
198 G4LorentzVector gLv = gamma->Get4Momentum();
199 G4LorentzVector lv(0.0,0.0,0.0,e);
200 lv -= gLv;
201 G4double mass = lv.m();
202 //G4cout << "mass= " << mass << " " << lv << G4endl;
203 G4ThreeVector boost = lv.boostVector();
204 //G4cout << "mass= " << mass << " " << boost << G4endl;
205 const G4ThreeVector dir = gamma->GetMomentumDirection();
206 model->SampleSecondaries(newp, mass, dir);
207 G4int np = newp->size();
208 for(G4int j=0; j<np; ++j) {
209 G4DynamicParticle* dp = (*newp)[j];
211 v.boost(boost);
212 //G4cout << j << ". " << v << G4endl;
213 v.boost(inBoost);
214 //G4cout << " " << v << G4endl;
215 dp->Set4Momentum(v);
216 t -= v.e();
217 }
218 //G4cout << "Gamma " << gLv << G4endl;
219 gLv.boost(inBoost);
220 //G4cout << " " << gLv << G4endl;
221 gamma->Set4Momentum(gLv);
222 t -= gLv.e();
223 newp->push_back(gamma);
224 if(std::abs(t) > MeV) {
225 G4cout << "G4eeToHadronsModel::SampleSecondaries: Ebalance(MeV)= "
226 << t/MeV << " primary 4-momentum: " << inlv << G4endl;
227 }
228 }
229 }
230}
CLHEP::HepLorentzVector G4LorentzVector
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
const G4ThreeVector & GetMomentumDirection() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, G4double, const G4ThreeVector &)=0
G4DynamicParticle * GenerateCMPhoton(G4double)

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