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

#include <G4IonCoulombScatteringModel.hh>

+ Inheritance diagram for G4IonCoulombScatteringModel:

Public Member Functions

 G4IonCoulombScatteringModel (const G4String &nam="IonCoulombScattering")
 
 ~G4IonCoulombScatteringModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) final
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
void SetRecoilThreshold (G4double eth)
 
void SetHeavyIonCorr (G4int b)
 
G4int GetHeavyIonCorr ()
 
G4IonCoulombScatteringModeloperator= (const G4IonCoulombScatteringModel &right)=delete
 
 G4IonCoulombScatteringModel (const G4IonCoulombScatteringModel &)=delete
 
- 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 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 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 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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 = 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 71 of file G4IonCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4IonCoulombScatteringModel() [1/2]

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

Definition at line 75 of file G4IonCoulombScatteringModel.cc.

76 : G4VEmModel(nam),
77 cosThetaMin(1.0)
78{
79 fNistManager = G4NistManager::Instance();
81 theProton = G4Proton::Proton();
82
83 pCuts = nullptr;
84 currentMaterial = nullptr;
85 currentElement = nullptr;
86 currentCouple = nullptr;
87 fParticleChange = nullptr;
88
89 recoilThreshold = 0.*eV;
90 heavycorr =0;
91 particle = nullptr;
92 mass=0;
93 currentMaterialIndex = -1;
94
95 ioncross = new G4IonCoulombCrossSection();
96}
static G4NistManager * Instance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ ~G4IonCoulombScatteringModel()

G4IonCoulombScatteringModel::~G4IonCoulombScatteringModel ( )
override

Definition at line 101 of file G4IonCoulombScatteringModel.cc.

102{
103 delete ioncross;
104}

◆ G4IonCoulombScatteringModel() [2/2]

G4IonCoulombScatteringModel::G4IonCoulombScatteringModel ( const G4IonCoulombScatteringModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 125 of file G4IonCoulombScatteringModel.cc.

130{
131 SetupParticle(p);
132
133 G4double cross = 0.0;
134
135 DefineMaterial(CurrentCouple());
136
137 G4int iz = G4lrint(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 return cross;
148}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
void SetupTarget(G4double Z, G4double kinEnergy, G4int heavycorr)
void SetupKinematic(G4double kinEnergy, G4double tmass)
G4double GetAtomicMassAmu(const G4String &symb) const
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486
int G4lrint(double ad)
Definition: templates.hh:134

◆ GetHeavyIonCorr()

G4int G4IonCoulombScatteringModel::GetHeavyIonCorr ( )
inline

Definition at line 178 of file G4IonCoulombScatteringModel.hh.

179{
180 return heavycorr;
181}

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 108 of file G4IonCoulombScatteringModel.cc.

110{
111 SetupParticle(p);
112 currentCouple = nullptr;
113 currentMaterialIndex = -1;
114 ioncross->Initialise(p,cosThetaMin);
115
116 pCuts = &cuts;
117 // G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
118 if(!fParticleChange) {
119 fParticleChange = GetParticleChangeForGamma();
120 }
121}
void Initialise(const G4ParticleDefinition *, G4double cosThetaLim)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 152 of file G4IonCoulombScatteringModel.cc.

157{
158 G4double kinEnergy = dp->GetKineticEnergy();
159 DefineMaterial(couple);
160 SetupParticle(dp->GetDefinition());
161
162 // Choose nucleus
163 currentElement = SelectTargetAtom(couple, particle, kinEnergy,
164 dp->GetLogKineticEnergy());
165
166 G4int iz = currentElement->GetZasInt();
167 G4int ia = SelectIsotopeNumber(currentElement);
169
170 ioncross->SetupKinematic(kinEnergy, mass2);
171 ioncross->SetupTarget(currentElement->GetZ(), kinEnergy, heavycorr);
172
173 //scattering angle, z1 == (1-cost)
174 G4double z1 = ioncross->SampleCosineTheta();
175 if(z1 > 2.0) { z1 = 2.0; }
176 else if(z1 < 0.0) { z1 = 0.0; }
177 /*
178 G4cout << "Sample: " << particle->GetParticleName()
179 << " mass(GeV)= " << mass/GeV
180 << " Ekin(MeV)= " << kinEnergy << " cost= " << 1. - z1 << G4endl;
181 G4cout << " Z= " << iz << " A= " << ia
182 << " mass(GeV)= " << mass2/GeV << G4endl;
183 */
184 G4double cost = 1.0 - z1;
185 G4double sint = sqrt(z1*(1.0 + cost));
186 G4double phi = twopi * G4UniformRand();
187
188 // kinematics in the Lab system
189 G4double ptot = sqrt(kinEnergy*(kinEnergy + 2.0*mass));
190 G4double e1 = mass + kinEnergy;
191
192 // Lab. system kinematics along projectile direction
193 G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1+mass2);
194 G4LorentzVector v1 = G4LorentzVector(0, 0, ptot, e1);
195 G4ThreeVector bst = v0.boostVector();
196 v1.boost(-bst);
197 // CM projectile
198 G4double momCM = v1.pz();
199
200 // Momentum after scattering of incident particle
201 v1.setX(momCM*sint*cos(phi));
202 v1.setY(momCM*sint*sin(phi));
203 v1.setZ(momCM*cost);
204
205 // CM--->Lab
206 v1.boost(bst);
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 = std::max(v0.e() - mass2, 0.0);
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) {
228 G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
229 newDirection = v0.vect().unit();
230 newDirection.rotateUz(dir);
231 auto 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 //G4cout << "Efinal(MeV)= " << finalT << " Edep(MeV)= " << edep
246 // << " Trec(MeV)= " << trec << G4endl;
247 fParticleChange->SetProposedKineticEnergy(finalT);
248 fParticleChange->ProposeLocalEnergyDeposit(edep);
249}
CLHEP::HepLorentzVector G4LorentzVector
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4int GetZasInt() const
Definition: G4Element.hh:132
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4double GetNuclearMass(const G4double A, const G4double Z)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int SelectIsotopeNumber(const G4Element *) const
Definition: G4VEmModel.cc:276
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetHeavyIonCorr()

void G4IonCoulombScatteringModel::SetHeavyIonCorr ( G4int  b)
inline

Definition at line 171 of file G4IonCoulombScatteringModel.hh.

172{
173 heavycorr = b;
174}

◆ SetRecoilThreshold()

void G4IonCoulombScatteringModel::SetRecoilThreshold ( G4double  eth)
inline

Definition at line 164 of file G4IonCoulombScatteringModel.hh.

165{
166 recoilThreshold = eth;
167}

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