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

#include <G4eSingleCoulombScatteringModel.hh>

Inheritance diagram for G4eSingleCoulombScatteringModel:
Collaboration diagram for G4eSingleCoulombScatteringModel:

Public Member Functions

 G4eSingleCoulombScatteringModel (const G4String &nam="eSingleCoulombScat")
 
virtual ~G4eSingleCoulombScatteringModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) final
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) final
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
void SetRecoilThreshold (G4double eth)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 * > *)
 
virtual 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=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 66 of file G4eSingleCoulombScatteringModel.hh.

Constructor & Destructor Documentation

G4eSingleCoulombScatteringModel::G4eSingleCoulombScatteringModel ( const G4String nam = "eSingleCoulombScat")
explicit

Definition at line 75 of file G4eSingleCoulombScatteringModel.cc.

76  : G4VEmModel(nam),
77  cosThetaMin(1.0)
78 {
79  fNistManager = G4NistManager::Instance();
81  fParticleChange = nullptr;
82 
83  pCuts=nullptr;
84  currentMaterial = nullptr;
85  currentElement = nullptr;
86  currentCouple = nullptr;
87 
88  lowEnergyLimit = 0*keV;
89  recoilThreshold = 0.*eV;
90  FormFactor = 0;
91  particle = nullptr;
92  mass=0.0;
93  currentMaterialIndex = -1;
94 
95  Mottcross = new G4ScreeningMottCrossSection();
96 }
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4IonTable * GetIonTable() const
static constexpr double eV
Definition: G4SIunits.hh:215
static G4ParticleTable * GetParticleTable()
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4eSingleCoulombScatteringModel::~G4eSingleCoulombScatteringModel ( )
virtual

Definition at line 100 of file G4eSingleCoulombScatteringModel.cc.

101 {
102  delete Mottcross;
103 }

Member Function Documentation

G4double G4eSingleCoulombScatteringModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 148 of file G4eSingleCoulombScatteringModel.cc.

155 {
156  SetupParticle(p);
157 
158  G4double cross =0.0;
159  if(kinEnergy < lowEnergyLimit) return cross;
160 
161  DefineMaterial(CurrentCouple());
162 
163  //Total Cross section
164  Mottcross->SetupKinematic(kinEnergy, Z);
165  cross = Mottcross->NuclearCrossSection(FormFactor);
166 
167  //cout<< "Compute Cross Section....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<" Z: "<<Z<<" kinEnergy: "<<kinEnergy<<endl;
168  return cross;
169 }
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
double G4double
Definition: G4Types.hh:76
void SetupKinematic(G4double kinEnergy, G4double Z)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 107 of file G4eSingleCoulombScatteringModel.cc.

109 {
111 
112  SetupParticle(p);
113  currentCouple = nullptr;
114  currentMaterialIndex = -1;
115  //cosThetaMin = cos(PolarAngleLimit());
116  Mottcross->Initialise(p,cosThetaMin);
117 
118  pCuts = &cuts;
119  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
120 
121  /*
122  G4cout << "!!! G4eSingleCoulombScatteringModel::Initialise for "
123  << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
124  << " cos(TetMax)= " << cosThetaMax <<G4endl;
125  G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
126  */
127 
128  if(!fParticleChange) {
129  fParticleChange = GetParticleChangeForGamma();
130  }
131 
132  if(IsMaster()) {
134  }
135 
136  FormFactor=param->NuclearFormfactorType();
137 
138  //G4cout<<"NUCLEAR FORM FACTOR: "<<FormFactor<<G4endl;
139 }
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
G4NuclearFormfactorType NuclearFormfactorType() const
static G4EmParameters * Instance()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

void G4eSingleCoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 141 of file G4eSingleCoulombScatteringModel.cc.

143 {
145 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 173 of file G4eSingleCoulombScatteringModel.cc.

179 {
180  G4double kinEnergy = dp->GetKineticEnergy();
181  //cout<<"--- kinEnergy "<<kinEnergy<<endl;
182 
183  if(kinEnergy < lowEnergyLimit) return;
184 
185  DefineMaterial(couple);
186  SetupParticle(dp->GetDefinition());
187 
188  // Choose nucleus
189  //last two :cutEnergy= min e kinEnergy=max
190  currentElement = SelectRandomAtom(couple,particle,
191  kinEnergy,cutEnergy,kinEnergy);
192 
193  G4double Z = currentElement->GetZ();
194  G4int iz = G4int(Z);
195  G4int ia = SelectIsotopeNumber(currentElement);
197 
198  //G4cout<<"..Z: "<<Z<<" ..iz: "<<iz<<" ..ia: "<<ia<<" ..mass2: "<<mass2<<G4endl;
199 
200  Mottcross->SetupKinematic(kinEnergy, Z);
201  G4double cross= Mottcross->NuclearCrossSection(FormFactor); //MODIFY TO LOAD TABLE
202  if(cross == 0.0) { return; }
203  //cout<< "Energy: "<<kinEnergy/MeV<<" Z: "<<Z<<"....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
204 
205  G4double z1 = Mottcross->GetScatteringAngle();
206  G4double sint = sin(z1);
207  G4double cost = sqrt(1.0 - sint*sint);
208  G4double phi = twopi* G4UniformRand();
209 
210  // kinematics in the Lab system
211  G4double ptot = dp->GetTotalMomentum();
212  G4double e1 = dp->GetTotalEnergy();
213  // Lab. system kinematics along projectile direction
214  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1);
215  G4double bet = ptot/(v0.e() + mass2);
216  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
217 
218  //CM Projectile
219  G4double momCM = gam*(ptot - bet*e1);
220  G4double eCM = gam*(e1 - bet*ptot);
221  //energy & momentum after scattering of incident particle
222  G4double pxCM = momCM*sint*cos(phi);
223  G4double pyCM = momCM*sint*sin(phi);
224  G4double pzCM = momCM*cost;
225 
226  //CM--->Lab
227  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
228 
229  // Rotate to global system
231  G4ThreeVector newDirection = v1.vect().unit();
232  newDirection.rotateUz(dir);
233 
234  fParticleChange->ProposeMomentumDirection(newDirection);
235 
236  // recoil
237  v0 -= v1;
238  G4double trec = v0.e();
239  G4double edep = 0.0;
240 
241  G4double tcut = recoilThreshold;
242 
243  //G4cout<<" Energy Transfered: "<<trec/eV<<G4endl;
244 
245  if(pCuts) {
246  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
247  }
248 
249  if(trec > tcut) {
250  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
251  newDirection = v0.vect().unit();
252  newDirection.rotateUz(dir);
253  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
254  fvect->push_back(newdp);
255  } else if(trec > 0.0) {
256  edep = trec;
257  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
258  }
259 
260  // finelize primary energy and energy balance
261  G4double finalT = v1.e() - mass;
262  //G4cout<<"Final Energy: "<<finalT/eV<<G4endl;
263  if(finalT <= lowEnergyLimit) {
264  edep += finalT;
265  finalT = 0.0;
266  }
267  edep = std::max(edep, 0.0);
268  fParticleChange->SetProposedKineticEnergy(finalT);
269  fParticleChange->ProposeLocalEnergyDeposit(edep);
270 
271 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
Definition: G4SIunits.hh:76
G4double GetTotalMomentum() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:584
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
void SetupKinematic(G4double kinEnergy, G4double Z)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:

void G4eSingleCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 155 of file G4eSingleCoulombScatteringModel.hh.

156 {
157  recoilThreshold = eth;
158 }

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