Geant4  10.02.p03
G4IonCoulombScatteringModel Class Reference

#include <G4IonCoulombScatteringModel.hh>

Inheritance diagram for G4IonCoulombScatteringModel:
Collaboration diagram for G4IonCoulombScatteringModel:

Public Member Functions

 G4IonCoulombScatteringModel (const G4String &nam="IonCoulombScattering")
 
virtual ~G4IonCoulombScatteringModel ()
 
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)
 
void SetHeavyIonCorr (G4int b)
 
G4int GetHeavyIonCorr ()
 
- 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 GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 Member Functions

void DefineMaterial (const G4MaterialCutsCouple *)
 
void SetupParticle (const G4ParticleDefinition *)
 
G4IonCoulombScatteringModeloperator= (const G4IonCoulombScatteringModel &right)
 
 G4IonCoulombScatteringModel (const G4IonCoulombScatteringModel &)
 

Private Attributes

G4IonTabletheIonTable
 
G4ParticleChangeForGamma * fParticleChange
 
G4NistManagerfNistManager
 
G4IonCoulombCrossSectionioncross
 
const std::vector< G4double > * pCuts
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const G4ElementcurrentElement
 
G4int currentMaterialIndex
 
G4int heavycorr
 
G4double cosThetaMin
 
G4double recoilThreshold
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4double mass
 

Additional Inherited Members

- 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 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 71 of file G4IonCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4IonCoulombScatteringModel() [1/2]

G4IonCoulombScatteringModel::G4IonCoulombScatteringModel ( const G4String nam = "IonCoulombScattering")

Definition at line 75 of file G4IonCoulombScatteringModel.cc.

76  : G4VEmModel(nam),
77  cosThetaMin(1.0)
78 {
82 
83  pCuts=0;
84  currentMaterial = 0;
85  currentElement = 0;
86  currentCouple = 0;
87  fParticleChange = 0;
88 
89  recoilThreshold = 0.*eV;
90  heavycorr =0;
91  particle = 0;
92  mass=0;
94 
96 }
const std::vector< G4double > * pCuts
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
G4ParticleChangeForGamma * fParticleChange
G4IonTable * GetIonTable() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
const G4ParticleDefinition * theProton
static const double eV
Definition: G4SIunits.hh:212
static G4ParticleTable * GetParticleTable()
const G4ParticleDefinition * particle
const G4MaterialCutsCouple * currentCouple
Here is the call graph for this function:

◆ ~G4IonCoulombScatteringModel()

G4IonCoulombScatteringModel::~G4IonCoulombScatteringModel ( )
virtual

Definition at line 101 of file G4IonCoulombScatteringModel.cc.

102 {
103  delete ioncross;
104 }

◆ G4IonCoulombScatteringModel() [2/2]

G4IonCoulombScatteringModel::G4IonCoulombScatteringModel ( const G4IonCoulombScatteringModel )
private

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 125 of file G4IonCoulombScatteringModel.cc.

130 {
131  SetupParticle(p);
132 
133  G4double cross = 0.0;
134 
136 
137  G4int iz = G4int(Z);
138 
139  //from lab to pCM & mu_rel of effective particle
140  G4double tmass = proton_mass_c2;
141  if(1 < iz) {
142  tmass = fNistManager->GetAtomicMassAmu(iz)*amu_c2;
143  }
144  ioncross->SetupKinematic(kinEnergy, tmass);
145  ioncross->SetupTarget(Z, kinEnergy, heavycorr);
146  cross = ioncross->NuclearCrossSection();
147 
148  //cout<< "..........cross "<<G4BestUnit(cross,"Surface") <<endl;
149  return cross;
150 }
void SetupTarget(G4double Z, G4double kinEnergy, G4int heavycorr)
void SetupKinematic(G4double kinEnergy, G4double tmass)
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
float proton_mass_c2
Definition: hepunit.py:275
void DefineMaterial(const G4MaterialCutsCouple *)
void SetupParticle(const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
float amu_c2
Definition: hepunit.py:277
G4double GetAtomicMassAmu(const G4String &symb) const
Here is the call graph for this function:

◆ DefineMaterial()

void G4IonCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple cup)
inlineprivate

Definition at line 143 of file G4IonCoulombScatteringModel.hh.

144 {
145  if(cup != currentCouple) {
146  currentCouple = cup;
147  currentMaterial = cup->GetMaterial();
149  }
150 }
const G4Material * GetMaterial() const
const G4MaterialCutsCouple * currentCouple
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetHeavyIonCorr()

G4int G4IonCoulombScatteringModel::GetHeavyIonCorr ( )
inline

Definition at line 180 of file G4IonCoulombScatteringModel.hh.

181 {
182  return heavycorr;
183 }

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 108 of file G4IonCoulombScatteringModel.cc.

110 {
111  SetupParticle(p);
112  currentCouple = 0;
115 
116  pCuts = &cuts;
117  // G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
118  if(!fParticleChange) {
120  }
121 }
const std::vector< G4double > * pCuts
G4ParticleChangeForGamma * fParticleChange
void SetupParticle(const G4ParticleDefinition *)
const G4MaterialCutsCouple * currentCouple
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
Here is the call graph for this function:

◆ operator=()

G4IonCoulombScatteringModel& G4IonCoulombScatteringModel::operator= ( const G4IonCoulombScatteringModel right)
private

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 154 of file G4IonCoulombScatteringModel.cc.

159 {
160  G4double kinEnergy = dp->GetKineticEnergy();
161 
162  DefineMaterial(couple);
163 
165 
166  // Choose nucleus
167  currentElement = SelectRandomAtom(couple, particle, kinEnergy);
168 
169  G4double Z = currentElement->GetZ();
170  G4int iz = G4int(Z);
171  G4int ia = SelectIsotopeNumber(currentElement);
173 
174  ioncross->SetupKinematic(kinEnergy, mass2);
175 
176  ioncross->SetupTarget(Z, kinEnergy, heavycorr);
177 
178  //scattering angle, z1 == (1-cost)
180  if(z1 > 2.0) { z1 = 2.0; }
181  else if(z1 < 0.0) { z1 = 0.0; }
182 
183  G4double cost = 1.0 - z1;
184  G4double sint = sqrt(z1*(1.0 + cost));
185  G4double phi = twopi * G4UniformRand();
186 
187  // kinematics in the Lab system
188  G4double ptot = dp->GetTotalMomentum();
189  G4double e1 = dp->GetTotalEnergy();
190 
191  // Lab. system kinematics along projectile direction
192  G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1);
193  G4double bet = ptot/(e1 + mass2);
194  G4double gam = 1.0/sqrt((1.0 - bet)*(1.0 + bet));
195 
196  // CM projectile
197  G4double momCM = gam*(ptot - bet*e1);
198  G4double eCM = gam*(e1 - bet*ptot);
199 
200  // Momentum after scattering of incident particle
201  G4double pxCM = momCM*sint*cos(phi);
202  G4double pyCM = momCM*sint*sin(phi);
203  G4double pzCM = momCM*cost;
204 
205  // CM--->Lab
206  G4LorentzVector v1(pxCM , pyCM, gam*(pzCM + bet*eCM), gam*(eCM + bet*pzCM));
207 
208  // Rotate to global system
210  G4ThreeVector newDirection = v1.vect().unit();
211  newDirection.rotateUz(dir);
212 
213  fParticleChange->ProposeMomentumDirection(newDirection);
214 
215  // recoil v0 energy is kinetic
216  v0 -= v1;
217  G4double trec = v0.e();
218  G4double edep = 0.0;
219 
220  G4double tcut = recoilThreshold;
221  if(pCuts) {
222  tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]);
223  //G4cout<<" tcut eV "<<tcut/eV<<endl;
224  }
225 
226  // Recoil
227  if(trec > tcut) {
229  newDirection = v0.vect().unit();
230  newDirection.rotateUz(dir);
231  G4DynamicParticle* newdp = new G4DynamicParticle(ion, newDirection, trec);
232  fvect->push_back(newdp);
233  } else if(trec > 0.0) {
234  edep = trec;
235  fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
236  }
237 
238  // finelize primary energy and energy balance
239  G4double finalT = v1.e() - mass;
240  if(finalT < 0.0) {
241  edep += finalT;
242  finalT = 0.0;
243  }
244  edep = std::max(edep, 0.0);
245  fParticleChange->SetProposedKineticEnergy(finalT);
246  fParticleChange->ProposeLocalEnergyDeposit(edep);
247 }
void SetupTarget(G4double Z, G4double kinEnergy, G4int heavycorr)
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetupKinematic(G4double kinEnergy, G4double tmass)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4double GetTotalMomentum() const
const std::vector< G4double > * pCuts
TDirectory * dir
int G4int
Definition: G4Types.hh:78
Double_t edep
G4double GetTotalEnergy() const
G4ParticleChangeForGamma * fParticleChange
Hep3Vector vect() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
Float_t Z
G4double iz
Definition: TRTMaterials.hh:39
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static const double twopi
Definition: G4SIunits.hh:75
G4int SelectIsotopeNumber(const G4Element *)
Definition: G4VEmModel.hh:585
static const G4double e1
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * particle
void DefineMaterial(const G4MaterialCutsCouple *)
G4ParticleDefinition * GetDefinition() const
void SetupParticle(const G4ParticleDefinition *)
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:544
CLHEP::HepLorentzVector G4LorentzVector
Here is the call graph for this function:

◆ SetHeavyIonCorr()

void G4IonCoulombScatteringModel::SetHeavyIonCorr ( G4int  b)
inline

Definition at line 173 of file G4IonCoulombScatteringModel.hh.

◆ SetRecoilThreshold()

void G4IonCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 166 of file G4IonCoulombScatteringModel.hh.

167 {
168  recoilThreshold = eth;
169 }

◆ SetupParticle()

void G4IonCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition p)
inlineprivate

Definition at line 155 of file G4IonCoulombScatteringModel.hh.

156 {
157  if(p != particle) {
158  particle = p;
159  mass = particle->GetPDGMass();
161  }
162 }
void SetupParticle(const G4ParticleDefinition *)
const G4ParticleDefinition * particle
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ cosThetaMin

G4double G4IonCoulombScatteringModel::cosThetaMin
private

Definition at line 130 of file G4IonCoulombScatteringModel.hh.

◆ currentCouple

const G4MaterialCutsCouple* G4IonCoulombScatteringModel::currentCouple
private

Definition at line 123 of file G4IonCoulombScatteringModel.hh.

◆ currentElement

const G4Element* G4IonCoulombScatteringModel::currentElement
private

Definition at line 125 of file G4IonCoulombScatteringModel.hh.

◆ currentMaterial

const G4Material* G4IonCoulombScatteringModel::currentMaterial
private

Definition at line 124 of file G4IonCoulombScatteringModel.hh.

◆ currentMaterialIndex

G4int G4IonCoulombScatteringModel::currentMaterialIndex
private

Definition at line 126 of file G4IonCoulombScatteringModel.hh.

◆ fNistManager

G4NistManager* G4IonCoulombScatteringModel::fNistManager
private

Definition at line 119 of file G4IonCoulombScatteringModel.hh.

◆ fParticleChange

G4ParticleChangeForGamma* G4IonCoulombScatteringModel::fParticleChange
private

Definition at line 118 of file G4IonCoulombScatteringModel.hh.

◆ heavycorr

G4int G4IonCoulombScatteringModel::heavycorr
private

Definition at line 128 of file G4IonCoulombScatteringModel.hh.

◆ ioncross

G4IonCoulombCrossSection* G4IonCoulombScatteringModel::ioncross
private

Definition at line 120 of file G4IonCoulombScatteringModel.hh.

◆ mass

G4double G4IonCoulombScatteringModel::mass
private

Definition at line 136 of file G4IonCoulombScatteringModel.hh.

◆ particle

const G4ParticleDefinition* G4IonCoulombScatteringModel::particle
private

Definition at line 134 of file G4IonCoulombScatteringModel.hh.

◆ pCuts

const std::vector<G4double>* G4IonCoulombScatteringModel::pCuts
private

Definition at line 122 of file G4IonCoulombScatteringModel.hh.

◆ recoilThreshold

G4double G4IonCoulombScatteringModel::recoilThreshold
private

Definition at line 131 of file G4IonCoulombScatteringModel.hh.

◆ theIonTable

G4IonTable* G4IonCoulombScatteringModel::theIonTable
private

Definition at line 117 of file G4IonCoulombScatteringModel.hh.

◆ theProton

const G4ParticleDefinition* G4IonCoulombScatteringModel::theProton
private

Definition at line 135 of file G4IonCoulombScatteringModel.hh.


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