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

#include <G4LEPTSElasticModel.hh>

+ Inheritance diagram for G4LEPTSElasticModel:

Public Member Functions

 G4LEPTSElasticModel (const G4String &modelName="G4LEPTSElasticModel")
 
 ~G4LEPTSElasticModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
 
- Public Member Functions inherited from G4VLEPTSModel
 G4VLEPTSModel (const G4String &processName)
 
 ~G4VLEPTSModel () override
 
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 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VLEPTSModel
void Init ()
 
G4bool ReadParam (G4String fileName, const G4Material *aMaterial)
 
virtual std::map< G4int, std::vector< G4double > > ReadIXS (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
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
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 {false}
 
XSType theXSType
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 32 of file G4LEPTSElasticModel.hh.

Constructor & Destructor Documentation

◆ G4LEPTSElasticModel()

G4LEPTSElasticModel::G4LEPTSElasticModel ( const G4String & modelName = "G4LEPTSElasticModel")

Definition at line 29 of file G4LEPTSElasticModel.cc.

30 : G4VLEPTSModel( modelName )
31{
33 fParticleChangeForGamma = nullptr;
34} // constructor
@ XSElastic
G4VLEPTSModel(const G4String &processName)

◆ ~G4LEPTSElasticModel()

G4LEPTSElasticModel::~G4LEPTSElasticModel ( )
overridedefault

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4LEPTSElasticModel::CrossSectionPerVolume ( const G4Material * mate,
const G4ParticleDefinition * aParticle,
G4double kineticEnergy,
G4double cutEnergy = 0.0,
G4double maxEnergy = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 70 of file G4LEPTSElasticModel.cc.

75{
76 if( kineticEnergy < theLowestEnergyLimit ) return DBL_MAX;
77 return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
78
79}
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
G4double theLowestEnergyLimit
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4LEPTSElasticModel::Initialise ( const G4ParticleDefinition * aParticle,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 42 of file G4LEPTSElasticModel.cc.

44{
45 Init();
46 BuildPhysicsTable( *aParticle );
47
48 fParticleChangeForGamma = GetParticleChangeForGamma();
49
50 // static const G4double proton_mass_c2 = 938.272013 * MeV;
51 // static const G4double neutron_mass_c2 = 939.56536 * MeV;
52 // static const G4double h2o_mass_c2 = 8*neutron_mass_c2 + 10*(proton_mass_c2 + electron_mass_c2);
53 // G4cout << "mme " << h2o_mass_c2/MeV << " " << H2o_mass_c2/MeV << G4endl;
54
55 const G4MaterialTable * materialTable = G4Material::GetMaterialTable() ;
56 std::vector<G4Material*>::const_iterator matite;
57 for( matite = materialTable->begin(); matite != materialTable->end(); matite++ ) {
58 const G4Material * aMaterial = (*matite);
59 theMassTarget[aMaterial] = theMolecularMass[aMaterial] / (6.02214179e+23/CLHEP::mole) *CLHEP::c_light * CLHEP::c_light;
60 theMassProjectile[aMaterial] = CLHEP::electron_mass_c2;
61
62 if( verboseLevel >= 1) G4cout << "Material: " << aMaterial->GetName() << " MolecularMass: " << theMolecularMass[aMaterial]/(CLHEP::g/CLHEP::mole) << " g/mole "
63 << " MTarget: " << theMassTarget[aMaterial]/CLHEP::MeV << " MeV" << G4endl;
64 }
65
66
67}
std::vector< G4Material * > G4MaterialTable
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
std::map< const G4Material *, G4double > theMolecularMass
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 83 of file G4LEPTSElasticModel.cc.

88{
89 G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
90 G4ThreeVector P0Dir = aDynamicParticle->GetMomentumDirection();
91
92 if( P0KinEn < theLowestEnergyLimit ) {
93 fParticleChangeForGamma->ProposeMomentumDirection( P0Dir );
94 fParticleChangeForGamma->SetProposedKineticEnergy( 0.);
95 fParticleChangeForGamma->ProposeLocalEnergyDeposit( P0KinEn);
96 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
97 if( verboseLevel > 2 ) G4cout << " ENERGY LOW " << P0KinEn - theLowestEnergyLimit << G4endl;
98 return;
99 }
100
101 //- G4ParticleDefinition * particleDefDef = aTrack.GetDefinition();
102 //- G4String partName = particleDefDef->GetParticleName();
103
104 // G4ThreeVector pos, pos0, dpos;
105
106 //- G4StepPoG4int * PostPoG4int = aStep.GetPostStepPoG4int();
107 //- G4ThreeVector r = PostPoG4int->GetPosition();
108
109 //TypeOfInteraction=-10;
110
111 const G4Material* aMaterial = mateCuts->GetMaterial();
112 G4double ang = SampleAngle(aMaterial, P0KinEn/CLHEP::eV, 0.0);
113 G4ThreeVector P1Dir = SampleNewDirection(P0Dir, ang);
114#ifdef DEBUG_LEPTS
115 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( P1Dir " << P1Dir << " P0Dir " << P0Dir << " ang " << ang << G4endl;
116#endif
117
118 //G4ThreeVector P1Dir = SampleNewDirection(P0Dir, P0KinEn/eV, 0.0);
119 //G4double Energylost1= ElasticEnergyTransferWater2(P0KinEn, ang);
120 G4double Energylost = EnergyTransfer(P0KinEn, ang, theMassTarget[aMaterial], theMassProjectile[aMaterial]);
121 if( verboseLevel >= 3 ) G4cout << " ELASTIC Energylost "<< Energylost << " = " << P0KinEn << " " <<ang << " " << theMassTarget[aMaterial] << " " << theMassProjectile[aMaterial] << G4endl;
122
123 G4double P1KinEn = P0KinEn - Energylost;
124 if( verboseLevel >= 3 ) G4cout << " ELASTIC " << P1KinEn << " = " << P0KinEn << " - " << Energylost << G4endl;
125#ifdef DEBUG_LEPTS
126 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( SetProposedKineticEnergy " << P1KinEn << " " << P0KinEn << " - " << Energylost << G4endl;
127#endif
128 fParticleChangeForGamma->ProposeMomentumDirection( P1Dir );
129 fParticleChangeForGamma->SetProposedKineticEnergy( P1KinEn);
130 fParticleChangeForGamma->ProposeLocalEnergyDeposit( Energylost);
131 //G4cout << "elasticEnergyLost: " << Energylost << G4endl;
132
133#ifdef DEBUG_LEPTS
134 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( ProposeMomentumDirection " << fParticleChangeForGamma->GetProposedMomentumDirection() << G4endl;
135#endif
136}
@ fStopAndKill
double G4double
Definition G4Types.hh:83
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4ThreeVector & GetProposedMomentumDirection() const
G4ThreeVector SampleNewDirection(const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
G4double SampleAngle(const G4Material *aMaterial, G4double e, G4double el)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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