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

#include <G4EmMultiModel.hh>

Inheritance diagram for G4EmMultiModel:
Collaboration diagram for G4EmMultiModel:

Public Member Functions

 G4EmMultiModel (const G4String &nam="MultiModel")
 
virtual ~G4EmMultiModel ()
 
void AddModel (G4VEmModel *)
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) final
 
- 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 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 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 * > *)
 
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 G4EmMultiModel.hh.

Constructor & Destructor Documentation

G4EmMultiModel::G4EmMultiModel ( const G4String nam = "MultiModel")
explicit

Definition at line 53 of file G4EmMultiModel.cc.

54  : G4VEmModel(nam), nModels(0)
55 {
56  model.clear();
57  cross_section.clear();
58 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
const XML_Char XML_Content * model
Definition: expat.h:151
G4EmMultiModel::~G4EmMultiModel ( )
virtual

Definition at line 62 of file G4EmMultiModel.cc.

63 {}

Member Function Documentation

void G4EmMultiModel::AddModel ( G4VEmModel p)

Definition at line 67 of file G4EmMultiModel.cc.

68 {
69  cross_section.push_back(0.0);
70  model.push_back(p);
71  ++nModels;
72 }
const XML_Char XML_Content * model
Definition: expat.h:151
G4double G4EmMultiModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 112 of file G4EmMultiModel.cc.

118 {
119  G4double cross = 0.0;
120  if(nModels>0) {
121  for(G4int i=0; i<nModels; ++i) {
123  cross += (model[i])->ComputeCrossSectionPerAtom(p, kinEnergy, Z, A,
124  cutEnergy, maxEnergy);
125  }
126  }
127  return cross;
128 }
int G4int
Definition: G4Types.hh:78
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
double A(double temperature)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) final
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
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 G4EmMultiModel::ComputeDEDX ( const G4MaterialCutsCouple couple,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 93 of file G4EmMultiModel.cc.

97 {
98  SetCurrentCouple(couple);
99  G4double dedx = 0.0;
100 
101  if(nModels > 0) {
102  for(G4int i=0; i<nModels; i++) {
103  dedx += (model[i])->ComputeDEDX(couple, p, cutEnergy, kineticEnergy);
104  }
105  }
106 
107  return dedx;
108 }
int G4int
Definition: G4Types.hh:78
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
virtual G4double ComputeDEDX(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
double G4double
Definition: G4Types.hh:76
const XML_Char XML_Content * model
Definition: expat.h:151

Here is the call graph for this function:

void G4EmMultiModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
finalvirtual

Implements G4VEmModel.

Definition at line 76 of file G4EmMultiModel.cc.

78 {
79  if(nModels > 0) {
80  G4cout << "### Initialisation of EM MultiModel " << GetName()
81  << " including following list of models:" << G4endl;
82  for(G4int i=0; i<nModels; ++i) {
83  G4cout << " " << (model[i])->GetName();
85  (model[i])->Initialise(p, cuts);
86  }
87  G4cout << G4endl;
88  }
89 }
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:418
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:609
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4VParticleChange * pParticleChange
Definition: G4VEmModel.hh:422
#define G4endl
Definition: G4ios.hh:61
const G4String & GetName() const
Definition: G4VEmModel.hh:794
const XML_Char XML_Content * model
Definition: expat.h:151

Here is the call graph for this function:

void G4EmMultiModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  vdp,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  tmax 
)
finalvirtual

Implements G4VEmModel.

Definition at line 132 of file G4EmMultiModel.cc.

137 {
138  SetCurrentCouple(couple);
139  if(nModels > 0) {
140  G4int i;
141  G4double cross = 0.0;
142  for(i=0; i<nModels; ++i) {
143  cross += (model[i])->CrossSection(couple, dp->GetParticleDefinition(),
144  dp->GetKineticEnergy(), minEnergy, maxEnergy);
145  cross_section[i] = cross;
146  }
147 
148  cross *= G4UniformRand();
149 
150  for(i=0; i<nModels; ++i) {
151  if(cross <= cross_section[i]) {
152  (model[i])->SampleSecondaries(vdp, couple, dp, minEnergy, maxEnergy);
153  return;
154  }
155  }
156  }
157 }
G4double GetKineticEnergy() const
int G4int
Definition: G4Types.hh:78
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) final
#define G4UniformRand()
Definition: Randomize.hh:97
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:500
const G4ParticleDefinition * GetParticleDefinition() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:444
double G4double
Definition: G4Types.hh:76
const XML_Char XML_Content * model
Definition: expat.h:151

Here is the call graph for this function:


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