Geant4  10.02.p03
G4LEPTSExcitationModel Class Reference

#include <G4LEPTSExcitationModel.hh>

Inheritance diagram for G4LEPTSExcitationModel:
Collaboration diagram for G4LEPTSExcitationModel:

Public Member Functions

 G4LEPTSExcitationModel (const G4String &modelName="G4LEPTSExcitationModel")
 
 ~G4LEPTSExcitationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual std::map< G4int, std::vector< G4double > > ReadIXS (G4String, const G4Material *aMaterial)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
- Public Member Functions inherited from G4VLEPTSModel
 G4VLEPTSModel (const G4String &processName)
 
 ~G4VLEPTSModel ()
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
 
G4ThreeVector SampleNewDirection (const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
 
G4double SampleAngle (const G4Material *aMaterial, G4double e, G4double el)
 
G4ThreeVector SampleNewDirection (G4ThreeVector Dir, G4double ang)
 
G4VLEPTSModeloperator= (const G4VLEPTSModel &right)
 
 G4VLEPTSModel (const G4VLEPTSModel &)
 
- 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Private Attributes

G4ParticleChangeForGamma * fParticleChangeForGamma
 
G4double LowestExcitationEnergy
 
G4double LowestNeutralDisociationEnergy
 

Additional Inherited Members

- Protected Member Functions inherited from G4VLEPTSModel
void Init ()
 
G4bool ReadParam (G4String fileName, const G4Material *aMaterial)
 
G4double SampleEnergyLoss (const G4Material *aMaterial, G4double eMin, G4double eMax)
 
void BuildMeanFreePathTable (const G4Material *aMaterial, std::map< G4int, std::vector< G4double > > &integralXS)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VLEPTSModel
G4PhysicsTabletheMeanFreePathTable
 
G4double theLowestEnergyLimit
 
G4double theHighestEnergyLimit
 
G4int theNumbBinTable
 
std::map< const G4Material *, G4doubletheIonisPot
 
std::map< const G4Material *, G4doubletheIonisPotInt
 
std::map< const G4Material *, G4doubletheMolecularMass
 
std::map< const G4Material *, G4LEPTSDiffXS * > theDiffXS
 
std::map< const G4Material *, G4LEPTSDistribution * > theRMTDistr
 
std::map< const G4Material *, G4LEPTSElossDistr * > theElostDistr
 
std::map< const G4Material *, G4LEPTSDistribution * > theElostDistr2
 
std::map< const G4Material *, G4inttheNXSdat
 
std::map< const G4Material *, G4inttheNXSsub
 
G4bool isInitialised
 
XSType theXSType
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
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 31 of file G4LEPTSExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4LEPTSExcitationModel()

G4LEPTSExcitationModel::G4LEPTSExcitationModel ( const G4String modelName = "G4LEPTSExcitationModel")

Definition at line 30 of file G4LEPTSExcitationModel.cc.

31  : G4VLEPTSModel( modelName )
32 {
34 } // constructor
G4VLEPTSModel(const G4String &processName)

◆ ~G4LEPTSExcitationModel()

G4LEPTSExcitationModel::~G4LEPTSExcitationModel ( )

Definition at line 38 of file G4LEPTSExcitationModel.cc.

38  {
39 }

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4LEPTSExcitationModel::CrossSectionPerVolume ( const G4Material mate,
const G4ParticleDefinition aParticle,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 79 of file G4LEPTSExcitationModel.cc.

84 {
85  return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
86 
87 }
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
Here is the call graph for this function:

◆ Initialise()

void G4LEPTSExcitationModel::Initialise ( const G4ParticleDefinition aParticle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 43 of file G4LEPTSExcitationModel.cc.

45 {
46  Init();
47  BuildPhysicsTable( *aParticle );
48 
50 
53 
54 
55 }
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
G4ParticleChangeForGamma * fParticleChangeForGamma
Here is the call graph for this function:

◆ ReadIXS()

std::map< G4int, std::vector< G4double > > G4LEPTSExcitationModel::ReadIXS ( G4String  fileTXS,
const G4Material aMaterial 
)
virtual

Reimplemented from G4VLEPTSModel.

Definition at line 58 of file G4LEPTSExcitationModel.cc.

59 {
60  std::map<G4int,std::vector<G4double> > integralXS = G4VLEPTSModel::ReadIXS( fileTXS, aMaterial);
61 
62  if( integralXS.size() == 0 ) return integralXS;
63 
64  for (G4int jj=theNXSdat[aMaterial]; jj>=0; jj--) {
65  if( integralXS[XSDissociation][jj] > 0.001) LowestExcitationEnergy = integralXS[XSTotal][jj-1];
66  if( integralXS[XSVibration][jj] > 0.001) LowestNeutralDisociationEnergy = integralXS[XSTotal][jj-1]*CLHEP::eV;
67  // if( txs[5][j] > 0.001) LowestExcitationEnergy = txs[0][j-1];
68  // if( txs[6][j] > 0.001) LowestNeutralDisociationEnergy = txs[0][j-1]*CLHEP::eV;
69  }
70 
71  if( verboseLevel >= 1) G4cout << " LowestExcitationEnergy: " << LowestExcitationEnergy << G4endl
72  << "LowestNeutralDisociationEnergy: " << LowestNeutralDisociationEnergy/CLHEP::eV
73  << G4endl;
74 
75  return integralXS;
76 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
virtual std::map< G4int, std::vector< G4double > > ReadIXS(G4String fileName, const G4Material *aMaterial)
#define G4endl
Definition: G4ios.hh:61
static const double eV
std::map< const G4Material *, G4int > theNXSdat
Here is the call graph for this function:

◆ SampleSecondaries()

void G4LEPTSExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle *> *  ,
const G4MaterialCutsCouple mateCuts,
const G4DynamicParticle aDynamicParticle,
G4double  tmin = 0.0,
G4double  tmax = DBL_MAX 
)
virtual

Implements G4VEmModel.

Definition at line 91 of file G4LEPTSExcitationModel.cc.

96 {
97  G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
98 
99  G4double Edep=0;
100  G4double Energylost=0;
101  G4ThreeVector P0Dir = aDynamicParticle->GetMomentumDirection();
102 
103  G4double eMin = 0.0;
104  const G4Material* aMaterial = mateCuts->GetMaterial();
105  G4double eMax = std::min(theIonisPot[aMaterial], P0KinEn);
106  Energylost = SampleEnergyLoss(aMaterial, eMin, eMax);
107 
108  Edep = Energylost;
109 
110  G4ThreeVector P1Dir = SampleNewDirection(aMaterial, P0Dir, P0KinEn/CLHEP::eV, Energylost/CLHEP::eV);
111  G4double P1KinEn = P0KinEn - Edep;
112 
113  fParticleChangeForGamma->ProposeMomentumDirection( P1Dir);
114  fParticleChangeForGamma->SetProposedKineticEnergy( P1KinEn);
115  fParticleChangeForGamma->ProposeLocalEnergyDeposit( Edep);
116 
117 }
const G4Material * GetMaterial() const
std::map< const G4Material *, G4double > theIonisPot
G4double GetKineticEnergy() const
G4double SampleEnergyLoss(const G4Material *aMaterial, G4double eMin, G4double eMax)
const G4ThreeVector & GetMomentumDirection() const
G4ThreeVector SampleNewDirection(const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
double G4double
Definition: G4Types.hh:76
static const double eV
G4ParticleChangeForGamma * fParticleChangeForGamma
Here is the call graph for this function:

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4LEPTSExcitationModel::fParticleChangeForGamma
private

Definition at line 56 of file G4LEPTSExcitationModel.hh.

◆ LowestExcitationEnergy

G4double G4LEPTSExcitationModel::LowestExcitationEnergy
private

Definition at line 57 of file G4LEPTSExcitationModel.hh.

◆ LowestNeutralDisociationEnergy

G4double G4LEPTSExcitationModel::LowestNeutralDisociationEnergy
private

Definition at line 58 of file G4LEPTSExcitationModel.hh.


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