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

#include <G4BraggIonModel.hh>

+ Inheritance diagram for G4BraggIonModel:

Public Member Functions

 G4BraggIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="BraggIon")
 
 ~G4BraggIonModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
 
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss) override
 
G4BraggIonModeloperator= (const G4BraggIonModel &right)=delete
 
 G4BraggIonModel (const G4BraggIonModel &)=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
 

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
- 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 66 of file G4BraggIonModel.hh.

Constructor & Destructor Documentation

◆ G4BraggIonModel() [1/2]

G4BraggIonModel::G4BraggIonModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "BraggIon" 
)
explicit

Definition at line 79 of file G4BraggIonModel.cc.

81 : G4VEmModel(nam),
82 theElectron(G4Electron::Electron()),
83 HeMass(3.727417*CLHEP::GeV),
84 theZieglerFactor(CLHEP::eV*CLHEP::cm2*1.0e-15),
85 lowestKinEnergy(0.25*CLHEP::keV)
86{
87 SetHighEnergyLimit(2.0*CLHEP::MeV);
88
89 rateMassHe2p = HeMass/CLHEP::proton_mass_c2;
90 massFactor = 1000.*CLHEP::amu_c2/HeMass;
91
92 if(nullptr != p) { SetParticle(p); }
93 else { SetParticle(theElectron); }
94}
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746

◆ ~G4BraggIonModel()

G4BraggIonModel::~G4BraggIonModel ( )
override

Definition at line 98 of file G4BraggIonModel.cc.

99{
100 if(IsMaster()) { delete fASTAR; fASTAR = nullptr; }
101}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4BraggIonModel() [2/2]

G4BraggIonModel::G4BraggIonModel ( const G4BraggIonModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4BraggIonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 197 of file G4BraggIonModel.cc.

203{
204 G4double sigma =
205 Z*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
206 if(isAlpha) {
207 sigma *= (HeEffChargeSquare(Z, kinEnergy/CLHEP::MeV)/chargeSquare);
208 }
209 return sigma;
210}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

G4double G4BraggIonModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)

Definition at line 166 of file G4BraggIonModel.cc.

171{
172 G4double cross = 0.0;
173 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
174 const G4double maxEnergy = std::min(tmax, maxKinEnergy);
175 const G4double cutEnergy = std::max(lowestKinEnergy*massRate, minKinEnergy);
176
177 if(cutEnergy < tmax) {
178
179 const G4double energy = kineticEnergy + mass;
180 const G4double energy2 = energy*energy;
181 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
182
183 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
184 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
185 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
186
187 cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
188 cross = std::max(cross, 0.0);
189 }
190 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
191 // << " tmax= " << tmax << " cross= " << cross << G4endl;
192 return cross;
193}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double energy(const ThreeVector &p, const G4double m)

Referenced by ComputeCrossSectionPerAtom(), G4LindhardSorensenIonModel::ComputeCrossSectionPerElectron(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

G4double G4BraggIonModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 233 of file G4BraggIonModel.cc.

237{
238 const G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
239 const G4double tmin = std::max(lowestKinEnergy*massRate, minKinEnergy);
240 G4double dedx = 0.0;
241
242 // T is alpha energy
243 G4double T = kineticEnergy;
244 const G4double zeff = material->GetTotNbOfElectPerVolume()/
245 material->GetTotNbOfAtomsPerVolume();
246 heChargeSquare = HeEffChargeSquare(zeff, T/CLHEP::MeV);
247 if(!isAlpha) { T *= rateMassHe2p; }
248
249 if(T < lowestKinEnergy) {
250 dedx = DEDX(material, lowestKinEnergy)*std::sqrt(T/lowestKinEnergy);
251 } else {
252 dedx = DEDX(material, T);
253 }
254 if(!isAlpha) { dedx /= heChargeSquare; }
255 if (tmin < tmax) {
256 const G4double tau = kineticEnergy/mass;
257 const G4double x = tmin/tmax;
258
259 G4double del =
260 (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x) *
261 CLHEP::twopi_mc2_rcl2*material->GetElectronDensity();
262 if(isAlpha) { del *= heChargeSquare; }
263 dedx += del;
264 }
265 dedx = std::max(dedx, 0.0);
266 /*
267 G4cout << "BraggIon: tkin(MeV) = " << tkin/MeV << " dedx(MeV*cm^2/g) = "
268 << dedx*gram/(MeV*cm2*material->GetDensity())
269 << " q2 = " << chargeSquare << G4endl;
270 */
271 return dedx;
272}
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4double GetTotNbOfElectPerVolume() const
Definition: G4Material.hh:207
G4double GetElectronDensity() const
Definition: G4Material.hh:212

Referenced by G4BraggNoDeltaModel::ComputeDEDXPerVolume(), and G4LindhardSorensenIonModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

void G4BraggIonModel::CorrectionsAlongStep ( const G4MaterialCutsCouple couple,
const G4DynamicParticle dp,
const G4double length,
G4double eloss 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 276 of file G4BraggIonModel.cc.

280{
281 // no correction for alpha
282 if(isAlpha) { return; }
283
284 // no correction at a small step at the last step
285 const G4double preKinEnergy = dp->GetKineticEnergy();
286 if(eloss >= preKinEnergy || eloss < preKinEnergy*0.05) { return; }
287
288 // corrections only for ions
289 const G4ParticleDefinition* p = dp->GetDefinition();
290 if(p != particle) { SetParticle(p); }
291
292 // effective energy and charge at a step
293 const G4Material* mat = couple->GetMaterial();
294 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
295 const G4double q20 = corr->EffectiveChargeSquareRatio(p, mat, preKinEnergy);
296 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
297 const G4double qfactor = q2/q20;
298 /*
299 G4cout << "G4BraggIonModel::CorrectionsAlongStep: Epre(MeV)="
300 << preKinEnergy << " Eeff(MeV)=" << e
301 << " eloss=" << eloss << " elossnew=" << eloss*qfactor
302 << " qfactor=" << qfactor << " Qpre=" << q20
303 << p->GetParticleName() <<G4endl;
304 */
305 eloss *= qfactor;
306}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4Material * GetMaterial() const

◆ CrossSectionPerVolume()

G4double G4BraggIonModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 214 of file G4BraggIonModel.cc.

220{
221 G4double sigma = material->GetElectronDensity()*
222 ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
223 if(isAlpha) {
224 const G4double zeff = material->GetTotNbOfElectPerVolume()/
225 material->GetTotNbOfAtomsPerVolume();
226 sigma *= (HeEffChargeSquare(zeff, kinEnergy/CLHEP::MeV)/chargeSquare);
227 }
228 return sigma;
229}

◆ GetChargeSquareRatio()

G4double G4BraggIonModel::GetChargeSquareRatio ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 148 of file G4BraggIonModel.cc.

151{
152 return corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
153}

◆ GetParticleCharge()

G4double G4BraggIonModel::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 157 of file G4BraggIonModel.cc.

160{
161 return corr->GetParticleCharge(p,mat,kineticEnergy);
162}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ Initialise()

void G4BraggIonModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector  
)
overridevirtual

Implements G4VEmModel.

Definition at line 105 of file G4BraggIonModel.cc.

107{
108 if(p != particle) { SetParticle(p); }
109
110 // always false before the run
111 SetDeexcitationFlag(false);
112
113 // initialise once
114 if(nullptr == fParticleChange) {
115 const G4String& pname = particle->GetParticleName();
116 if(IsMaster()) {
117 if(pname == "proton" || pname == "GenericIon" || pname == "alpha") {
118 if(nullptr == fASTAR) { fASTAR = new G4ASTARStopping(); }
119 fASTAR->Initialise();
120
121 if(G4EmParameters::Instance()->UseICRU90Data()) {
123 fICRU90->Initialise();
124 }
125 }
126 }
127 if(pname == "alpha") { isAlpha = true; }
128
129 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
131 }
133
134 fParticleChange = GetParticleChangeForLoss();
135 }
136}
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
G4ICRU90StoppingData * GetICRU90StoppingData()
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:600
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:697
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ MaxSecondaryEnergy()

G4double G4BraggIonModel::MaxSecondaryEnergy ( const G4ParticleDefinition pd,
G4double  kinEnergy 
)
finalprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 388 of file G4BraggIonModel.cc.

390{
391 if(pd != particle) { SetParticle(pd); }
392 G4double tau = kinEnergy/mass;
393 G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
394 (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
395 return tmax;
396}

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ MinEnergyCut()

G4double G4BraggIonModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple couple 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 140 of file G4BraggIonModel.cc.

142{
144}
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 310 of file G4BraggIonModel.cc.

315{
316 const G4double tmax = MaxSecondaryKinEnergy(dp);
317 const G4double xmax = std::min(tmax, maxEnergy);
318 const G4double xmin = std::max(lowestKinEnergy*massRate, minEnergy);
319 if(xmin >= xmax) { return; }
320
321 G4double kineticEnergy = dp->GetKineticEnergy();
322 const G4double energy = kineticEnergy + mass;
323 const G4double energy2 = energy*energy;
324 const G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
325 const G4double grej = 1.0;
326 G4double deltaKinEnergy, f;
327
328 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
329 G4double rndm[2];
330
331 // sampling follows ...
332 do {
333 rndmEngineMod->flatArray(2, rndm);
334 deltaKinEnergy = xmin*xmax/(xmin*(1.0 - rndm[0]) + xmax*rndm[0]);
335
336 f = 1.0 - beta2*deltaKinEnergy/tmax;
337
338 if(f > grej) {
339 G4cout << "G4BraggIonModel::SampleSecondary Warning! "
340 << "Majorant " << grej << " < "
341 << f << " for e= " << deltaKinEnergy
342 << G4endl;
343 }
344
345 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
346 } while( grej*rndm[1] >= f );
347
348 G4ThreeVector deltaDirection;
349
351 const G4Material* mat = couple->GetMaterial();
353
354 deltaDirection =
355 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
356
357 } else {
358
359 G4double deltaMomentum =
360 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
361 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
362 (deltaMomentum * dp->GetTotalMomentum());
363 if(cost > 1.0) { cost = 1.0; }
364 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
365
366 G4double phi = twopi*rndmEngineMod->flat();
367
368 deltaDirection.set(sint*cos(phi),sint*sin(phi), cost) ;
369 deltaDirection.rotateUz(dp->GetMomentumDirection());
370 }
371
372 // create G4DynamicParticle object for delta ray
373 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
374
375 vdp->push_back(delta);
376
377 // Change kinematics of primary particle
378 kineticEnergy -= deltaKinEnergy;
379 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
380 finalP = finalP.unit();
381
382 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
383 fParticleChange->SetProposedMomentumDirection(finalP);
384}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const
Definition: G4VEmModel.cc:253
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:501

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