Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions
G4eSingleCoulombScatteringModel Class Reference

#include <G4eSingleCoulombScatteringModel.hh>

Inheritance diagram for G4eSingleCoulombScatteringModel:
G4VEmModel

Public Member Functions

 G4eSingleCoulombScatteringModel (const G4String &nam="eSingleCoulombScattering")
 
virtual ~G4eSingleCoulombScatteringModel ()
 
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 SetRecoilThreshold (G4double eth)
 
- 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 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 *)
 
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
 

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 *)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 67 of file G4eSingleCoulombScatteringModel.hh.

Constructor & Destructor Documentation

G4eSingleCoulombScatteringModel::G4eSingleCoulombScatteringModel ( const G4String nam = "eSingleCoulombScattering")

Definition at line 74 of file G4eSingleCoulombScatteringModel.cc.

References python.hepunit::eV, G4ParticleTable::GetIonTable(), G4ParticleTable::GetParticleTable(), and G4NistManager::Instance().

75  : G4VEmModel(nam),
76 
77  cosThetaMin(1.0),
78  isInitialised(false)
79 {
80  fNistManager = G4NistManager::Instance();
82  fParticleChange = 0;
83 
84  pCuts=0;
85  currentMaterial = 0;
86  currentElement = 0;
87  currentCouple = 0;
88 
89  lowEnergyLimit = 0*eV;
90  recoilThreshold = 0.*eV;
91  particle = 0;
92  mass=0;
93  currentMaterialIndex = -1;
94 
95  Mottcross = new G4ScreeningMottCrossSection();
96 
97 }
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
G4eSingleCoulombScatteringModel::~G4eSingleCoulombScatteringModel ( )
virtual

Definition at line 101 of file G4eSingleCoulombScatteringModel.cc.

102 { delete Mottcross;}

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 127 of file G4eSingleCoulombScatteringModel.cc.

References G4VEmModel::CurrentCouple(), DefineMaterial(), G4ScreeningMottCrossSection::NuclearCrossSection(), G4ScreeningMottCrossSection::SetupKinematic(), and SetupParticle().

134 {
135 
136  SetupParticle(p);
137 
138  G4double cross =0.0;
139  if(kinEnergy < lowEnergyLimit) return cross;
140 
142 
143  //Total Cross section
144  Mottcross->SetupKinematic(kinEnergy, Z);
145  cross = Mottcross->NuclearCrossSection();
146 
147  //cout<< "....cross "<<G4BestUnit(cross,"Surface") << " cm2 "<< cross/cm2 <<endl;
148  return cross;
149 }
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:426
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
void SetupKinematic(G4double kinEnergy, G4double Z)
void G4eSingleCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprotected

Definition at line 136 of file G4eSingleCoulombScatteringModel.hh.

References G4Material::GetIndex(), and G4MaterialCutsCouple::GetMaterial().

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

137 {
138  if(cup != currentCouple) {
139  currentCouple = cup;
140  currentMaterial = cup->GetMaterial();
141  currentMaterialIndex = currentCouple->GetIndex();
142 
143  }
144 }
const G4Material * GetMaterial() const
void G4eSingleCoulombScatteringModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 106 of file G4eSingleCoulombScatteringModel.cc.

References G4VEmModel::GetParticleChangeForGamma(), G4ScreeningMottCrossSection::Initialise(), G4VEmModel::PolarAngleLimit(), and SetupParticle().

108 {
109  SetupParticle(p);
110  currentCouple = 0;
111  currentMaterialIndex = -1;
112  cosThetaMin = cos(PolarAngleLimit());
113  Mottcross->Initialise(p,cosThetaMin);
114 
115  pCuts = &cuts;
116  //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
117 
118 
119  if(!isInitialised) {
120  isInitialised = true;
121  fParticleChange = GetParticleChangeForGamma();
122  }
123 }
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
void SetupParticle(const G4ParticleDefinition *)
G4double PolarAngleLimit() const
Definition: G4VEmModel.hh:620
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4eSingleCoulombScatteringModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 153 of file G4eSingleCoulombScatteringModel.cc.

References DefineMaterial(), G4DynamicParticle::GetDefinition(), G4IonTable::GetIon(), G4DynamicParticle::GetKineticEnergy(), G4ScreeningMottCrossSection::GetMom2Lab(), G4DynamicParticle::GetMomentumDirection(), G4ScreeningMottCrossSection::GetNewDirection(), G4ScreeningMottCrossSection::GetTotalCross(), G4ScreeningMottCrossSection::GetTrec(), iz, G4INCL::Math::min(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeNonIonizingEnergyDeposit(), CLHEP::Hep3Vector::rotateUz(), G4VEmModel::SelectIsotopeNumber(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and SetupParticle().

159 {
160  G4double kinEnergy = dp->GetKineticEnergy();
161  //cout<<"--- kinEnergy "<<kinEnergy<<endl;
162 
163  if(kinEnergy < lowEnergyLimit) return;
164 
165  DefineMaterial(couple);
167 
168  // Choose nucleus
169  //last two :cutEnergy= min e kinEnergy=max
170  currentElement = SelectRandomAtom(couple,particle,
171  kinEnergy,cutEnergy,kinEnergy);
172 
173  G4double Z = currentElement->GetZ();
174  G4int iz = G4int(Z);
175  G4int ia = SelectIsotopeNumber(currentElement);
176 
177  //cout<<"Element "<<currentElement->GetName()<<endl;;
178 
179  G4double cross= Mottcross->GetTotalCross();
180 
181  if(cross == 0.0) { return; }
182 
183  G4ThreeVector dir = dp->GetMomentumDirection(); //old direction
184  G4ThreeVector newDirection=Mottcross->GetNewDirection();//new direction
185  newDirection.rotateUz(dir);
186 
187  fParticleChange->ProposeMomentumDirection(newDirection);
188 
189  //Recoil energy
190  G4double trec= Mottcross->GetTrec();
191 
192  //Energy after scattering
193  if(trec > kinEnergy) { trec = kinEnergy; }
194  G4double finalT = kinEnergy - trec;
195  G4double edep = 0.0;
196 
197  G4double tcut = recoilThreshold;
198  if(pCuts) {
199  tcut= std::min(tcut,(*pCuts)[currentMaterialIndex]);
200  }
201 
202  if(trec > tcut) {
203 
204  //cout<<"Trec "<<trec/eV<<endl;
205  G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
206 
207  //incident before scattering
208  G4double ptot=sqrt(Mottcross->GetMom2Lab());
209  //incident after scattering
210  G4double plab = sqrt(finalT*(finalT + 2.0*mass));
211  G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit();
212  //secondary particle
213  G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec);
214  fvect->push_back(newdp);
215  } else if(trec > 0.0) {
216  edep = trec;
217  fParticleChange->ProposeNonIonizingEnergyDeposit(trec);
218  }
219 
220  // finelize primary energy and energy balance
221  if(finalT <= lowEnergyLimit) {
222  edep += finalT;
223  finalT = 0.0;
224  }
225  fParticleChange->SetProposedKineticEnergy(finalT);
226  fParticleChange->ProposeLocalEnergyDeposit(edep);
227 
228 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:548
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
void G4eSingleCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 159 of file G4eSingleCoulombScatteringModel.hh.

160 {
161  recoilThreshold = eth;
162 }
void G4eSingleCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 149 of file G4eSingleCoulombScatteringModel.hh.

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

150 {
151  if(p != particle) {
152  particle = p;
153  mass = particle->GetPDGMass();
154  Mottcross->SetupParticle(p);
155  }
156 }
const char * p
Definition: xmltok.h:285
void SetupParticle(const G4ParticleDefinition *)
G4double GetPDGMass() const

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