Geant4 11.2.2
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 (G4bool combined=true)
 
 ~G4hCoulombScatteringModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) final
 
void SetLowEnergyThreshold (G4double val)
 
void SetRecoilThreshold (G4double eth)
 
void SetFixedCut (G4double)
 
G4double GetFixedCut () const
 
G4hCoulombScatteringModeloperator= (const G4hCoulombScatteringModel &right)=delete
 
 G4hCoulombScatteringModel (const G4hCoulombScatteringModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 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 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)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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

Detailed Description

Definition at line 66 of file G4hCoulombScatteringModel.hh.

Constructor & Destructor Documentation

◆ G4hCoulombScatteringModel() [1/2]

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( G4bool combined = true)
explicit

Definition at line 64 of file G4hCoulombScatteringModel.cc.

65 : G4VEmModel("hCoulombScattering"),
66 cosThetaMin(1.0),
67 cosThetaMax(-1.0),
68 isCombined(combined)
69{
70 fParticleChange = nullptr;
71 fNistManager = G4NistManager::Instance();
73 theProton = G4Proton::Proton();
74 currentMaterial = nullptr;
75 fixedCut = -1.0;
76
77 pCuts = nullptr;
78
79 recoilThreshold = 0.0; // by default does not work
80
81 particle = nullptr;
82 currentCouple = nullptr;
83 wokvi = new G4WentzelVIRelXSection();
84
85 currentMaterialIndex = 0;
86 mass = CLHEP::proton_mass_c2;
87 elecRatio = 0.0;
88}
static G4NistManager * Instance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4hCoulombScatteringModel()

G4hCoulombScatteringModel::~G4hCoulombScatteringModel ( )
override

Definition at line 92 of file G4hCoulombScatteringModel.cc.

93{
94 delete wokvi;
95}

◆ G4hCoulombScatteringModel() [2/2]

G4hCoulombScatteringModel::G4hCoulombScatteringModel ( const G4hCoulombScatteringModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 180 of file G4hCoulombScatteringModel.cc.

185{
186 //G4cout << "### G4hCoulombScatteringModel::ComputeCrossSectionPerAtom for "
187 //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl;
188 G4double cross = 0.0;
189 elecRatio = 0.0;
190 if(p != particle) { SetupParticle(p); }
191
192 // cross section is set to zero to avoid problems in sample secondary
193 if(kinEnergy <= 0.0) { return cross; }
195
196 G4int iz = G4lrint(Z);
197 G4double tmass = (1 == iz) ? proton_mass_c2 :
198 fNistManager->GetAtomicMassAmu(iz)*amu_c2;
199 wokvi->SetTargetMass(tmass);
200
201 G4double costmin =
202 wokvi->SetupKinematic(kinEnergy, currentMaterial);
203
204 if(cosThetaMax < costmin) {
205 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
206 costmin = wokvi->SetupTarget(iz, cut);
207 G4double costmax =
208 (1 == iz && particle == theProton && cosThetaMax < 0.0)
209 ? 0.0 : cosThetaMax;
210 if(costmin > costmax) {
211 cross = wokvi->ComputeNuclearCrossSection(costmin, costmax)
212 + wokvi->ComputeElectronCrossSection(costmin, costmax);
213 }
214 /*
215 if(p->GetParticleName() == "mu+")
216 G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn
217 << " 1-costmin= " << 1-costmin
218 << " 1-costmax= " << 1-costmax
219 << " 1-cosThetaMax= " << 1-cosThetaMax
220 << G4endl;
221 */
222 }
223 return cross;
224}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4MaterialCutsCouple * CurrentCouple() const
G4double SetupTarget(G4int Z, G4double cut)
G4double ComputeElectronCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double ComputeNuclearCrossSection(G4double CosThetaMin, G4double CosThetaMax)
G4double SetupKinematic(G4double kinEnergy, const G4Material *mat) override
void SetupParticle(const G4ParticleDefinition *)
void DefineMaterial(const G4MaterialCutsCouple *)
int G4lrint(double ad)
Definition templates.hh:134

◆ DefineMaterial()

void G4hCoulombScatteringModel::DefineMaterial ( const G4MaterialCutsCouple * cup)
inlineprotected

Definition at line 150 of file G4hCoulombScatteringModel.hh.

151{
152 if(cup != currentCouple) {
153 currentCouple = cup;
154 currentMaterial = cup->GetMaterial();
155 currentMaterialIndex = currentCouple->GetIndex();
156 }
157}
const G4Material * GetMaterial() const

Referenced by ComputeCrossSectionPerAtom(), and SampleSecondaries().

◆ GetFixedCut()

G4double G4hCoulombScatteringModel::GetFixedCut ( ) const
inline

Definition at line 188 of file G4hCoulombScatteringModel.hh.

189{
190 return fixedCut;
191}

◆ Initialise()

void G4hCoulombScatteringModel::Initialise ( const G4ParticleDefinition * part,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 99 of file G4hCoulombScatteringModel.cc.

101{
102 SetupParticle(part);
103 currentCouple = nullptr;
104
105 // defined theta limit between single and multiple scattering
106 isCombined = true;
108
109 if(tet <= 0.0) {
110 cosThetaMin = 1.0;
111 isCombined = false;
112 } else if(tet >= CLHEP::pi) {
113 cosThetaMin = -1.0;
114 } else {
115 cosThetaMin = std::cos(tet);
116 }
117
118 wokvi->Initialise(part, cosThetaMin);
119 /*
120 G4cout << "G4hCoulombScatteringModel: " << particle->GetParticleName()
121 << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin
122 << " cos(thetaMax)= " << cosThetaMax
123 << G4endl;
124 */
125 pCuts = &cuts;
126 //G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
127 /*
128 G4cout << "!!! G4hCoulombScatteringModel::Initialise for "
129 << part->GetParticleName() << " cos(TetMin)= " << cosThetaMin
130 << " cos(TetMax)= " << cosThetaMax <<G4endl;
131 G4cout << "cut= " << (*pCuts)[0] << " cut1= " << (*pCuts)[1] << G4endl;
132 */
133 if(!fParticleChange) {
134 fParticleChange = GetParticleChangeForGamma();
135 }
136 if(IsMaster() && mass < CLHEP::GeV && part->GetParticleName() != "GenericIon") {
137 InitialiseElementSelectors(part, cuts);
138 }
139}
G4double PolarAngleLimit() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsMaster() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
void Initialise(const G4ParticleDefinition *, G4double CosThetaLim)

◆ InitialiseLocal()

void G4hCoulombScatteringModel::InitialiseLocal ( const G4ParticleDefinition * ,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 143 of file G4hCoulombScatteringModel.cc.

145{
147}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ MinPrimaryEnergy()

G4double G4hCoulombScatteringModel::MinPrimaryEnergy ( const G4Material * material,
const G4ParticleDefinition * part,
G4double  )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4hCoulombScatteringModel.cc.

155{
156 SetupParticle(part);
157
158 // define cut using cuts for proton
159 G4double cut =
160 std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
161
162 // find out lightest element
163 const G4ElementVector* theElementVector = material->GetElementVector();
164 std::size_t nelm = material->GetNumberOfElements();
165
166 // select lightest element
167 G4int Z = 300;
168 for (std::size_t j=0; j<nelm; ++j) {
169 Z = std::min(Z,(*theElementVector)[j]->GetZasInt());
170 }
171 G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
173 G4double t = std::max(cut, 0.5*(cut + std::sqrt(2*cut*targetMass)));
174
175 return t;
176}
std::vector< const G4Element * > G4ElementVector
const G4double A[17]
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const
G4double GetAtomicMassAmu(const G4String &symb) const
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 228 of file G4hCoulombScatteringModel.cc.

234{
235 G4double kinEnergy = dp->GetKineticEnergy();
237 DefineMaterial(couple);
238
239 // Choose nucleus
240 G4double cut = (0.0 < fixedCut) ? fixedCut : cutEnergy;
241
242 const G4Element* elm = SelectRandomAtom(couple,particle,
243 kinEnergy,cut,kinEnergy);
244
245 G4int iz = elm->GetZasInt();
246 G4int ia = SelectIsotopeNumber(elm);
248
249 wokvi->SetTargetMass(mass2);
250 wokvi->SetupKinematic(kinEnergy, currentMaterial);
251 G4double costmin = wokvi->SetupTarget(iz, cut);
252 G4double costmax = (1 == iz && particle == theProton && cosThetaMax < 0.0)
253 ? 0.0 : cosThetaMax;
254 if(costmin <= costmax) { return; }
255
256 G4double cross = wokvi->ComputeNuclearCrossSection(costmin, costmax);
257 G4double ecross = wokvi->ComputeElectronCrossSection(costmin, costmax);
258 G4double ratio = ecross/(cross + ecross);
259
260 G4ThreeVector newDirection =
261 wokvi->SampleSingleScattering(costmin, costmax, ratio);
262
263 // kinematics in the Lab system
264 G4double ptot = std::sqrt(kinEnergy*(kinEnergy + 2.0*mass));
265 G4double e1 = mass + kinEnergy;
266
267 // Lab. system kinematics along projectile direction
268 G4LorentzVector v0 = G4LorentzVector(0, 0, ptot, e1+mass2);
269 G4LorentzVector v1 = G4LorentzVector(0, 0, ptot, e1);
270 G4ThreeVector bst = v0.boostVector();
271 v1.boost(-bst);
272 // CM projectile
273 G4double momCM = v1.pz();
274
275 // Momentum after scattering of incident particle
276 v1.setX(momCM*newDirection.x());
277 v1.setY(momCM*newDirection.y());
278 v1.setZ(momCM*newDirection.z());
279
280 // CM--->Lab
281 v1.boost(bst);
282
284 newDirection = v1.vect().unit();
285 newDirection.rotateUz(dir);
286
287 fParticleChange->ProposeMomentumDirection(newDirection);
288
289 // recoil
290 v0 -= v1;
291 G4double trec = std::max(v0.e() - mass2, 0.0);
292 G4double edep = 0.0;
293
294 G4double tcut = recoilThreshold;
295 if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
296
297 if(trec > tcut) {
298 G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
299 newDirection = v0.vect().unit();
300 newDirection.rotateUz(dir);
301 auto newdp = new G4DynamicParticle(ion, newDirection, trec);
302 fvect->push_back(newdp);
303 } else if(trec > 0.0) {
304 edep = trec;
305 fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
306 }
307
308 // finelize primary energy and energy balance
309 G4double finalT = v1.e() - mass;
310 if(finalT < 0.0) {
311 edep += finalT;
312 finalT = 0.0;
313 }
314 edep = std::max(edep, 0.0);
315 fParticleChange->SetProposedKineticEnergy(finalT);
316 fParticleChange->ProposeLocalEnergyDeposit(edep);
317}
CLHEP::HepLorentzVector G4LorentzVector
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int SelectIsotopeNumber(const G4Element *) const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ThreeVector & SampleSingleScattering(G4double CosThetaMin, G4double CosThetaMax, G4double elecRatio)

◆ SetFixedCut()

void G4hCoulombScatteringModel::SetFixedCut ( G4double val)
inline

Definition at line 181 of file G4hCoulombScatteringModel.hh.

182{
183 fixedCut = val;
184}

◆ SetLowEnergyThreshold()

void G4hCoulombScatteringModel::SetLowEnergyThreshold ( G4double val)
inline

◆ SetRecoilThreshold()

void G4hCoulombScatteringModel::SetRecoilThreshold ( G4double eth)
inline

Definition at line 174 of file G4hCoulombScatteringModel.hh.

175{
176 recoilThreshold = eth;
177}

◆ SetupParticle()

void G4hCoulombScatteringModel::SetupParticle ( const G4ParticleDefinition * p)
inlineprotected

Definition at line 162 of file G4hCoulombScatteringModel.hh.

163{
164 // Initialise mass and charge
165 if(p != particle) {
166 particle = p;
167 mass = particle->GetPDGMass();
168 wokvi->SetupParticle(p);
169 }
170}
void SetupParticle(const G4ParticleDefinition *)

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


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