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

#include <G4hCoulombScatteringModel.hh>

+ Inheritance diagram for G4hCoulombScatteringModel:

Public Member Functions

 G4hCoulombScatteringModel (const G4String &nam="eCoulombScattering")
 
virtual ~G4hCoulombScatteringModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetLowEnergyThreshold (G4double val)
 
void SetRecoilThreshold (G4double eth)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

void DefineMaterial (const G4MaterialCutsCouple *)
 
void SetupParticle (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleTabletheParticleTable
 
G4ParticleChangeForGammafParticleChange
 
G4WentzelVIRelXSectionwokvi
 
G4NistManagerfNistManager
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
G4int currentMaterialIndex
 
G4double cosThetaMin
 
G4double cosThetaMax
 
G4double cosTetMinNuc
 
G4double cosTetMaxNuc
 
G4double recoilThreshold
 
G4double elecRatio
 
G4double mass
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4double lowEnergyThreshold
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Detailed Description

Definition at line 67 of file G4hCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4hCoulombScatteringModel()

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( const G4String nam = "eCoulombScattering")

Definition at line 69 of file G4hCoulombScatteringModel.cc.

70 : G4VEmModel(nam),
71 cosThetaMin(1.0),
72 cosThetaMax(-1.0),
73 isInitialised(false)
74{
79 currentMaterial = 0;
80
81 pCuts = 0;
82
83 lowEnergyThreshold = 1*keV; // particle will be killed for lower energy
84 recoilThreshold = 0.*keV; // by default does not work
85
86 particle = 0;
87 currentCouple = 0;
89
91
92 cosTetMinNuc = 1.0;
93 cosTetMaxNuc = -1.0;
94 elecRatio = 0.0;
95 mass = proton_mass_c2;
96}
static G4NistManager * Instance()
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * theProton
const std::vector< G4double > * pCuts
const G4ParticleDefinition * particle
G4ParticleChangeForGamma * fParticleChange

◆ ~G4hCoulombScatteringModel()

G4hCoulombScatteringModel::~G4hCoulombScatteringModel ( )
virtual

Definition at line 100 of file G4hCoulombScatteringModel.cc.

101{
102 delete wokvi;
103}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4hCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 136 of file G4hCoulombScatteringModel.cc.

141{
142 //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
143 // << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
144 G4double cross = 0.0;
145 if(p != particle) { SetupParticle(p); }
146
147 // cross section is set to zero to avoid problems in sample secondary
148 if(kinEnergy <= 0.0) { return cross; }
152 G4int iz = G4int(Z);
153 cosTetMinNuc = wokvi->SetupTarget(iz, cutEnergy);
155 if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) {
156 cosTetMaxNuc = 0.0;
157 }
160 cross += elecRatio;
161 if(cross > 0.0) { elecRatio /= cross; }
162 }
163 /*
164 if(p->GetParticleName() == "e-")
165 G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
166 << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
167 << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
168 << G4endl;
169 */
170 return cross;
171}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:377
G4double SetupTarget(G4int Z, G4double cut=DBL_MAX)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
void SetupParticle(const G4ParticleDefinition *)
void DefineMaterial(const G4MaterialCutsCouple *)

Referenced by SampleSecondaries().

◆ DefineMaterial()

void G4hCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 145 of file G4hCoulombScatteringModel.hh.

146{
147 if(cup != currentCouple) {
148 currentCouple = cup;
151 }
152}
const G4Material * GetMaterial() const

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ Initialise()

void G4hCoulombScatteringModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 107 of file G4hCoulombScatteringModel.cc.

109{
110 SetupParticle(p);
111 currentCouple = 0;
114 /*
115 G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
116 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
117 << " cos(thetaMax)= " << cosThetaMax
118 << G4endl;
119 */
121 //G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
122 // << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin
123 // << " cos(TetMax)= " << cosThetaMax <<G4endl;
124 // G4cout << "cut0= " << cuts[0] << " cut1= " << cuts[1] << G4endl;
125 if(!isInitialised) {
126 isInitialised = true;
128 }
129 if(mass < GeV) {
131 }
132}
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:550
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)

◆ SampleSecondaries()

void G4hCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 175 of file G4hCoulombScatteringModel.cc.

181{
182 G4double kinEnergy = dp->GetKineticEnergy();
183
184 // absorb particle below low-energy limit to avoid situation
185 // when a particle has no energy loss
186 if(kinEnergy < lowEnergyThreshold) {
190 return;
191 }
193 DefineMaterial(couple);
194
195 //G4cout << "G4hCoulombScatteringModel::SampleSecondaries e(MeV)= "
196 // << kinEnergy << " " << particle->GetParticleName()
197 // << " cut= " << cutEnergy<< G4endl;
198
199 // Choose nucleus
200 const G4Element* currentElement =
201 SelectRandomAtom(couple,particle,kinEnergy,cutEnergy,kinEnergy);
202
203 G4double Z = currentElement->GetZ();
204
205 if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
206 kinEnergy, cutEnergy, kinEnergy) == 0.0)
207 { return; }
208
209 G4int iz = G4int(Z);
210 G4int ia = SelectIsotopeNumber(currentElement);
211 G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
212 wokvi->SetTargetMass(targetMass);
213
214 G4ThreeVector newDirection =
216 G4double cost = newDirection.z();
217
218 G4ThreeVector direction = dp->GetMomentumDirection();
219 newDirection.rotateUz(direction);
220
222
223 // recoil sampling assuming a small recoil
224 // and first order correction to primary 4-momentum
226 G4double trec = mom2*(1.0 - cost)/(targetMass + (mass + kinEnergy)*(1.0 - cost));
227 G4double finalT = kinEnergy - trec;
228 //G4cout<<"G4hCoulombScatteringModel: finalT= "<<finalT<<" Trec= "<<trec<<G4endl;
229 if(finalT <= lowEnergyThreshold) {
230 trec = kinEnergy;
231 finalT = 0.0;
232 }
233
236 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
237
238 if(trec > tcut) {
239 G4ParticleDefinition* ion = theParticleTable->GetIon(iz, ia, 0.0);
240 G4ThreeVector dir = (direction*sqrt(mom2) -
241 newDirection*sqrt(finalT*(2*mass + finalT))).unit();
242 G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
243 fvect->push_back(newdp);
244 } else {
247 }
248
249 return;
250}
double z() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ParticleDefinition * GetIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:478
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetTargetMass(G4double value)
G4ThreeVector SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio=0.0)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)

◆ SetLowEnergyThreshold()

void G4hCoulombScatteringModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 169 of file G4hCoulombScatteringModel.hh.

170{
171 lowEnergyThreshold = val;
172}

◆ SetRecoilThreshold()

void G4hCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 176 of file G4hCoulombScatteringModel.hh.

177{
178 recoilThreshold = eth;
179}

◆ SetupParticle()

void G4hCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 157 of file G4hCoulombScatteringModel.hh.

158{
159 // Initialise mass and charge
160 if(p != particle) {
161 particle = p;
164 }
165}
void SetupParticle(const G4ParticleDefinition *)

Referenced by ComputeCrossSectionPerAtom(), Initialise(), and SampleSecondaries().

Member Data Documentation

◆ cosTetMaxNuc

G4double G4hCoulombScatteringModel::cosTetMaxNuc
protected

◆ cosTetMinNuc

G4double G4hCoulombScatteringModel::cosTetMinNuc
protected

◆ cosThetaMax

G4double G4hCoulombScatteringModel::cosThetaMax
protected

Definition at line 124 of file G4hCoulombScatteringModel.hh.

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ cosThetaMin

G4double G4hCoulombScatteringModel::cosThetaMin
protected

Definition at line 123 of file G4hCoulombScatteringModel.hh.

Referenced by Initialise().

◆ currentCouple

const G4MaterialCutsCouple* G4hCoulombScatteringModel::currentCouple
protected

◆ currentMaterial

const G4Material* G4hCoulombScatteringModel::currentMaterial
protected

◆ currentMaterialIndex

G4int G4hCoulombScatteringModel::currentMaterialIndex
protected

◆ elecRatio

G4double G4hCoulombScatteringModel::elecRatio
protected

◆ fNistManager

G4NistManager* G4hCoulombScatteringModel::fNistManager
protected

Definition at line 115 of file G4hCoulombScatteringModel.hh.

Referenced by G4hCoulombScatteringModel().

◆ fParticleChange

G4ParticleChangeForGamma* G4hCoulombScatteringModel::fParticleChange
protected

◆ lowEnergyThreshold

G4double G4hCoulombScatteringModel::lowEnergyThreshold
protected

◆ mass

G4double G4hCoulombScatteringModel::mass
protected

◆ particle

const G4ParticleDefinition* G4hCoulombScatteringModel::particle
protected

◆ pCuts

const std::vector<G4double>* G4hCoulombScatteringModel::pCuts
protected

◆ recoilThreshold

G4double G4hCoulombScatteringModel::recoilThreshold
protected

◆ theParticleTable

G4ParticleTable* G4hCoulombScatteringModel::theParticleTable
protected

Definition at line 112 of file G4hCoulombScatteringModel.hh.

Referenced by G4hCoulombScatteringModel(), and SampleSecondaries().

◆ theProton

const G4ParticleDefinition* G4hCoulombScatteringModel::theProton
protected

◆ wokvi


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