Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4eeToHadronsMultiModel Class Reference

#include <G4eeToHadronsMultiModel.hh>

Inheritance diagram for G4eeToHadronsMultiModel:
G4VEmModel

Public Member Functions

 G4eeToHadronsMultiModel (G4int ver=0, const G4String &nam="eeToHadrons")
 
virtual ~G4eeToHadronsMultiModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, 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)
 
virtual void PrintInfo ()
 
void SetCrossSecFactor (G4double fac)
 
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
- 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 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 *)
 
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)
 
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=0)
 
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
 

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
 

Detailed Description

Definition at line 62 of file G4eeToHadronsMultiModel.hh.

Constructor & Destructor Documentation

G4eeToHadronsMultiModel::G4eeToHadronsMultiModel ( G4int  ver = 0,
const G4String nam = "eeToHadrons" 
)

Definition at line 65 of file G4eeToHadronsMultiModel.cc.

References DBL_MAX, and python.hepunit::GeV.

66  : G4VEmModel(mname),
67  csFactor(1.0),
68  nModels(0),
69  verbose(ver),
70  isInitialised(false)
71 {
72  thKineticEnergy = DBL_MAX;
73  maxKineticEnergy = 1.2*GeV;
74  fParticleChange = 0;
75  cross = 0;
76 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
#define DBL_MAX
Definition: templates.hh:83
G4eeToHadronsMultiModel::~G4eeToHadronsMultiModel ( )
virtual

Definition at line 80 of file G4eeToHadronsMultiModel.cc.

References n.

81 {
82  G4int n = models.size();
83  if(n>0) {
84  for(G4int i=0; i<n; i++) {
85  delete models[i];
86  }
87  }
88  delete cross;
89 }
int G4int
Definition: G4Types.hh:78
const G4int n

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 168 of file G4eeToHadronsMultiModel.cc.

References ComputeCrossSectionPerElectron().

173 {
174  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
175 }
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4double G4eeToHadronsMultiModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
inline

Definition at line 132 of file G4eeToHadronsMultiModel.hh.

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

136 {
137  G4double res = 0.0;
138  if (kineticEnergy > thKineticEnergy) {
139  for(G4int i=0; i<nModels; i++) {
140  if(kineticEnergy >= ekinMin[i] && kineticEnergy <= ekinMax[i])
141  res += (models[i])->ComputeCrossSectionPerElectron(0,kineticEnergy);
142  cumSum[i] = res;
143  }
144  }
145  return res*csFactor;
146 }
int G4int
Definition: G4Types.hh:78
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
double G4double
Definition: G4Types.hh:76
G4double G4eeToHadronsMultiModel::CrossSectionPerVolume ( const G4Material mat,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 156 of file G4eeToHadronsMultiModel.cc.

References ComputeCrossSectionPerElectron(), and G4Material::GetElectronDensity().

161 {
162  return mat->GetElectronDensity()*
163  ComputeCrossSectionPerElectron(p, kineticEnergy);
164 }
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void G4eeToHadronsMultiModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 93 of file G4eeToHadronsMultiModel.cc.

References G4VEmModel::GetParticleChangeForGamma(), python.hepunit::GeV, G4Vee2hadrons::SetHighEnergy(), and G4Vee2hadrons::SetLowEnergy().

95 {
96  if(!isInitialised) {
97  isInitialised = true;
98 
99  cross = new G4eeCrossSections();
100 
101  G4eeToTwoPiModel* m2pi = new G4eeToTwoPiModel(cross);
102  m2pi->SetHighEnergy(maxKineticEnergy);
103  AddEEModel(m2pi);
104 
105  G4eeTo3PiModel* m3pi1 = new G4eeTo3PiModel(cross);
106  m3pi1->SetHighEnergy(0.95*GeV);
107  AddEEModel(m3pi1);
108 
109  G4eeTo3PiModel* m3pi2 = new G4eeTo3PiModel(cross);
110  m3pi2->SetLowEnergy(0.95*GeV);
111  m3pi2->SetHighEnergy(maxKineticEnergy);
112  AddEEModel(m3pi2);
113 
114  G4ee2KChargedModel* m2kc = new G4ee2KChargedModel(cross);
115  m2kc->SetHighEnergy(maxKineticEnergy);
116  AddEEModel(m2kc);
117 
118  G4ee2KNeutralModel* m2kn = new G4ee2KNeutralModel(cross);
119  m2kn->SetHighEnergy(maxKineticEnergy);
120  AddEEModel(m2kn);
121 
122  G4eeToPGammaModel* mpg1 = new G4eeToPGammaModel(cross,"pi0");
123  mpg1->SetLowEnergy(0.7*GeV);
124  mpg1->SetHighEnergy(maxKineticEnergy);
125  AddEEModel(mpg1);
126 
127  G4eeToPGammaModel* mpg2 = new G4eeToPGammaModel(cross,"eta");
128  mpg2->SetLowEnergy(0.7*GeV);
129  mpg2->SetHighEnergy(maxKineticEnergy);
130  AddEEModel(mpg2);
131 
132  nModels = models.size();
133 
134  fParticleChange = GetParticleChangeForGamma();
135  }
136 }
void SetLowEnergy(G4double val)
void SetHighEnergy(G4double val)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4eeToHadronsMultiModel::PrintInfo ( )
virtual

Definition at line 200 of file G4eeToHadronsMultiModel.cc.

References python.hepunit::electron_mass_c2, G4cout, G4endl, and python.hepunit::GeV.

Referenced by G4eeToHadrons::PrintInfo().

201 {
202  if(verbose > 0) {
203  G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
204  - 2.0*electron_mass_c2;
205  G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
206  - 2.0*electron_mass_c2;
207  G4cout << " e+ annihilation into hadrons active from "
208  << e1/GeV << " GeV to " << e2/GeV << " GeV"
209  << G4endl;
210  }
211 }
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4eeToHadronsMultiModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  newp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Implements G4VEmModel.

Definition at line 180 of file G4eeToHadronsMultiModel.cc.

References fStopAndKill, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), and G4VParticleChange::ProposeTrackStatus().

184 {
185  G4double kinEnergy = dp->GetKineticEnergy();
186  if (kinEnergy > thKineticEnergy) {
187  G4double q = cumSum[nModels-1]*G4UniformRand();
188  for(G4int i=0; i<nModels; i++) {
189  if(q <= cumSum[i]) {
190  (models[i])->SampleSecondaries(newp, couple,dp);
191  if(newp->size() > 0) fParticleChange->ProposeTrackStatus(fStopAndKill);
192  break;
193  }
194  }
195  }
196 }
G4double GetKineticEnergy() const
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX)
void G4eeToHadronsMultiModel::SetCrossSecFactor ( G4double  fac)

Definition at line 215 of file G4eeToHadronsMultiModel.cc.

References G4cout, and G4endl.

Referenced by G4eeToHadrons::InitialiseProcess(), and G4eeToHadrons::SetCrossSecFactor().

216 {
217  if(fac > 1.0) {
218  csFactor = fac;
219  if(verbose > 0)
220  G4cout << "### G4eeToHadronsMultiModel: The cross section for G4eeToHadronsMultiModel "
221  << " is increased by the Factor= " << csFactor << G4endl;
222  }
223 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

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