Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eeToHadronsModel Class Reference

#include <G4eeToHadronsModel.hh>

Inheritance diagram for G4eeToHadronsModel:
Collaboration 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 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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

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 std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 59 of file G4eeToHadronsModel.hh.

Constructor & Destructor Documentation

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

Definition at line 70 of file G4eeToHadronsModel.cc.

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

Here is the call graph for this function:

G4eeToHadronsModel::~G4eeToHadronsModel ( )
virtual

Definition at line 92 of file G4eeToHadronsModel.cc.

93 {
94  delete model;
95 }
const XML_Char XML_Content * model
Definition: expat.h:151

Member Function Documentation

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 162 of file G4eeToHadronsModel.cc.

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

Here is the call graph for this function:

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

Definition at line 173 of file G4eeToHadronsModel.cc.

177 {
178  G4double cross = 0.0;
179  if(crossPerElectron) {
180  cross = crossPerElectron->Value(energy);
181  }
182  return cross;
183 }
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 150 of file G4eeToHadronsModel.cc.

155 {
156  return mat->GetElectronDensity()*
157  ComputeCrossSectionPerElectron(p, kineticEnergy);
158 }
G4double GetElectronDensity() const
Definition: G4Material.hh:217
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)

Here is the call graph for this function:

G4DynamicParticle * G4eeToHadronsModel::GenerateCMPhoton ( G4double  e)

Definition at line 279 of file G4eeToHadronsModel.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VEmModel.

Definition at line 99 of file G4eeToHadronsModel.cc.

101 {
102  if(isInitialised) { return; }
103  isInitialised = true;
104 
105  // CM system
106  emin = model->LowEnergy();
107  emax = model->HighEnergy();
108 
109  // peak energy
110  epeak = std::min(model->PeakEnergy(), emax);
111 
112  if(verbose>0) {
113  G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
114  G4cout << "CM System: emin(MeV)= " << emin/MeV
115  << " epeak(MeV)= " << epeak/MeV
116  << " emax(MeV)= " << emax/MeV
117  << G4endl;
118  }
119 
120  crossBornPerElectron = model->PhysicsVector();
121  crossPerElectron = model->PhysicsVector();
122  nbins = crossPerElectron->GetVectorLength();
123  for(G4int i=0; i<nbins; ++i) {
124  G4double e = crossPerElectron->Energy(i);
125  G4double cs = model->ComputeCrossSection(e);
126  crossBornPerElectron->PutValue(i, cs);
127  }
128  ComputeCMCrossSectionPerElectron();
129 
130  if(verbose>1) {
131  G4cout << "G4eeToHadronsModel: Cross sections per electron"
132  << " nbins= " << nbins
133  << " emin(MeV)= " << emin/MeV
134  << " emax(MeV)= " << emax/MeV
135  << G4endl;
136  for(G4int i=0; i<nbins; ++i) {
137  G4double e = crossPerElectron->Energy(i);
138  G4double s1 = crossPerElectron->Value(e);
139  G4double s2 = crossBornPerElectron->Value(e);
140  G4cout << "E(MeV)= " << e/MeV
141  << " cross(nb)= " << s1/nanobarn
142  << " crossBorn(nb)= " << s2/nanobarn
143  << G4endl;
144  }
145  }
146 }
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static constexpr double nanobarn
Definition: G4SIunits.hh:108
void PutValue(size_t index, G4double theValue)
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
const XML_Char XML_Content * model
Definition: expat.h:151

Here is the call graph for this function:

G4double G4eeToHadronsModel::PeakEnergy ( ) const
inline

Definition at line 128 of file G4eeToHadronsModel.hh.

129 {
130  return peakKinEnergy;
131 }
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 187 of file G4eeToHadronsModel.cc.

192 {
193  if(crossPerElectron) {
194  G4double t = dParticle->GetKineticEnergy() + 2*electron_mass_c2;
195  G4LorentzVector inlv = dParticle->Get4Momentum() +
196  G4LorentzVector(0.0,0.0,0.0,electron_mass_c2);
197  G4double e = inlv.m();
198  G4ThreeVector inBoost = inlv.boostVector();
199  //G4cout << "G4eeToHadronsModel::SampleSecondaries e= " << e
200  // << " " << inlv << " " << inBoost <<G4endl;
201  if(e > emin) {
203  G4LorentzVector gLv = gamma->Get4Momentum();
204  G4LorentzVector lv(0.0,0.0,0.0,e);
205  lv -= gLv;
206  G4double mass = lv.m();
207  //G4cout << "mass= " << mass << " " << lv << G4endl;
208  G4ThreeVector boost = lv.boostVector();
209  //G4cout << "mass= " << mass << " " << boost << G4endl;
210  const G4ThreeVector dir = gamma->GetMomentumDirection();
211  model->SampleSecondaries(newp, mass, dir);
212  G4int np = newp->size();
213  for(G4int j=0; j<np; ++j) {
214  G4DynamicParticle* dp = (*newp)[j];
215  G4LorentzVector v = dp->Get4Momentum();
216  v.boost(boost);
217  //G4cout << j << ". " << v << G4endl;
218  v.boost(inBoost);
219  //G4cout << " " << v << G4endl;
220  dp->Set4Momentum(v);
221  t -= v.e();
222  }
223  //G4cout << "Gamma " << gLv << G4endl;
224  gLv.boost(inBoost);
225  //G4cout << " " << gLv << G4endl;
226  gamma->Set4Momentum(gLv);
227  t -= gLv.e();
228  newp->push_back(gamma);
229  if(fabs(t) > MeV) {
230  G4cout << "G4eeToHadronsModel::SampleSecondaries: Ebalance(MeV)= "
231  << t/MeV << " primary 4-momentum: " << inlv << G4endl;
232  }
233  }
234  }
235 }
Hep3Vector boostVector() const
G4double GetKineticEnergy() const
int G4int
Definition: G4Types.hh:78
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
HepLorentzVector & boost(double, double, double)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4DynamicParticle * GenerateCMPhoton(G4double)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
const XML_Char XML_Content * model
Definition: expat.h:151
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:


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