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

#include <G4LindhardSorensenIonModel.hh>

Inheritance diagram for G4LindhardSorensenIonModel:

Public Member Functions

 G4LindhardSorensenIonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LindhardSorensen")
 ~G4LindhardSorensenIonModel () 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
G4double GetChargeSquareRatio (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
void CorrectionsAlongStep (const G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, const G4double &length, G4double &eloss) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4LindhardSorensenIonModeloperator= (const G4LindhardSorensenIonModel &right)=delete
 G4LindhardSorensenIonModel (const G4LindhardSorensenIonModel &)=delete
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 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 void StartTracking (G4Track *)
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
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)
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)
void SetMasterThread (G4bool)
G4VEmModeloperator= (const G4VEmModel &right)=delete
 G4VEmModel (const G4VEmModel &)=delete

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) override
Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
G4ParticleChangeForGammaGetParticleChangeForGamma ()
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 62 of file G4LindhardSorensenIonModel.hh.

Constructor & Destructor Documentation

◆ G4LindhardSorensenIonModel() [1/2]

G4LindhardSorensenIonModel::G4LindhardSorensenIonModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LindhardSorensen" )
explicit

Definition at line 72 of file G4LindhardSorensenIonModel.cc.

74 : G4VEmModel(nam),
75 twoln10(2.0*G4Log(10.0))
76{
77 theElectron = G4Electron::Electron();
80 fBraggModel = new G4BraggModel();
81 fBBModel = new G4BetheBlochModel();
82 fElimit = 2.0*CLHEP::MeV;
83}
G4double G4Log(G4double x)
Definition G4Log.hh:169
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

Referenced by G4LindhardSorensenIonModel(), and operator=().

◆ ~G4LindhardSorensenIonModel()

G4LindhardSorensenIonModel::~G4LindhardSorensenIonModel ( )
override

Definition at line 87 of file G4LindhardSorensenIonModel.cc.

87 {
88 if(isFirst) {
89 delete lsdata;
90 delete fIonData;
91 lsdata = nullptr;
92 fIonData = nullptr;
93 }
94}

◆ G4LindhardSorensenIonModel() [2/2]

G4LindhardSorensenIonModel::G4LindhardSorensenIonModel ( const G4LindhardSorensenIonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 204 of file G4LindhardSorensenIonModel.cc.

210{
211 return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
212}
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 184 of file G4LindhardSorensenIonModel.cc.

189{
190 // take into account formfactor
191 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
192 G4double emax = std::min(tmax, maxKinEnergy);
193 G4double escaled = kinEnergy*pRatio;
194 G4double cross = (escaled <= fElimit)
195 ? fBraggModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax)
196 : fBBModel->ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,emax);
197 // G4cout << "LS: e= " << kinEnergy << " tmin= " << cutEnergy
198 // << " tmax= " << maxEnergy << " cross= " << cross << G4endl;
199 return cross;
200}
double G4double
Definition G4Types.hh:83
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

G4double G4LindhardSorensenIonModel::ComputeDEDXPerVolume ( const G4Material * mat,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 230 of file G4LindhardSorensenIonModel.cc.

234{
235 // formfactor is taken into account in CorrectionsAlongStep(..)
236 G4double tmax = MaxSecondaryEnergy(p, kinEnergy);
237 G4double cutEnergy = std::min(std::min(cut,tmax), tlimit);
238
239 G4double escaled = kinEnergy*pRatio;
240 G4double dedx = (escaled <= fElimit)
241 ? fBraggModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy)
242 : fBBModel->ComputeDEDXPerVolume(mat, p, kinEnergy, cutEnergy);
243
244 //G4cout << "E(MeV)=" << kinEnergy/MeV << " dedx=" << dedx
245 // << " " << mat->GetName() << " Ecut(MeV)=" << cutEnergy << G4endl;
246 return dedx;
247}

◆ CorrectionsAlongStep()

void G4LindhardSorensenIonModel::CorrectionsAlongStep ( const G4Material * mat,
const G4ParticleDefinition * p,
const G4double kinEnergy,
const G4double cutEnergy,
const G4double & length,
G4double & eloss )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 251 of file G4LindhardSorensenIonModel.cc.

258{
259 if(eloss >= preKinEnergy) { return; }
260
261 SetParticle(p);
262 const G4double eDensity = mat->GetElectronDensity();
263 const G4double e = std::max(preKinEnergy - eloss*0.5, preKinEnergy*0.5);
264 const G4double tmax = MaxSecondaryEnergy(p, e);
265 const G4double escaled = e*pRatio;
266 const G4double tau = e/mass;
267 const G4double q2 = corr->EffectiveChargeSquareRatio(p, mat, e);
268 const G4int Z = p->GetAtomicNumber();
269
270 G4double res = 0.0;
271 if (escaled <= fElimit) {
272 // data from ICRU73 or ICRU90
273 if(Z > 2 && Z <= 80) {
274 res = fIonData->GetDEDX(mat, Z, escaled, G4Log(escaled));
275 /*
276 G4cout << " GetDEDX for Z=" << Z << " in " << mat->GetName()
277 << " Escaled=" << escaled << " E="
278 << e << " dEdx=" << res << G4endl;
279 */
280 }
281 if (res > 0.0) {
282 if (cut < tmax) {
283 const G4double x = cut/tmax;
284 res += (G4Log(x)*(tau + 1.)*(tau + 1.)/(tau * (tau + 2.0)) + 1.0 - x)
285 *q2*CLHEP::twopi_mc2_rcl2*eDensity;
286 }
287 res *= length;
288 } else {
289 // simplified correction
290 res = eloss*q2/chargeSquareRatio;
291 }
292 } else {
293 // Lindhard-Sorensen model
294 const G4double gam = tau + 1.0;
295 const G4double beta2 = tau * (tau+2.0)/(gam*gam);
296 G4double deltaL0 = 2.0*corr->BarkasCorrection(p, mat, e)*(charge-1.)/charge;
297 G4double deltaL = lsdata->GetDeltaL(Zin, gam);
298
299 res = eloss +
300 CLHEP::twopi_mc2_rcl2*q2*eDensity*(deltaL+deltaL0)*length/beta2;
301 /*
302 G4cout << " E(GeV)=" << preKinEnergy/GeV << " eloss(MeV)=" << eloss
303 << " L= " << eloss*beta2/(twopi_mc2_rcl2*q2*eDensity*length)
304 << " dL0= " << deltaL0
305 << " dL= " << deltaL << " dE(MeV)=" << res - eloss << G4endl;
306 */
307 }
308 if(res > preKinEnergy || 2*res < eloss) { res = eloss; }
309 /*
310 G4cout << " G4LindhardSorensenIonModel::CorrectionsAlongStep: E(GeV)="
311 << preKinEnergy/GeV << " eloss(MeV)=" << eloss
312 << " res(MeV)=" << res << G4endl;
313 */
314 eloss = res;
315}
int G4int
Definition G4Types.hh:85
G4double GetElectronDensity() const
G4int GetAtomicNumber() const

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 216 of file G4LindhardSorensenIonModel.cc.

222{
223 return material->GetElectronDensity()
224 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
225}

◆ GetChargeSquareRatio()

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

Reimplemented from G4VEmModel.

Definition at line 131 of file G4LindhardSorensenIonModel.cc.

134{
135 chargeSquareRatio = corr->EffectiveChargeSquareRatio(p, mat, kinEnergy);
136 fBraggModel->SetChargeSquareRatio(chargeSquareRatio);
137 fBBModel->SetChargeSquareRatio(chargeSquareRatio);
138 return chargeSquareRatio;
139}

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 144 of file G4LindhardSorensenIonModel.cc.

147{
148 return corr->GetParticleCharge(p,mat,kinEnergy);
149}

◆ Initialise()

void G4LindhardSorensenIonModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector & ptr )
overridevirtual

Implements G4VEmModel.

Definition at line 98 of file G4LindhardSorensenIonModel.cc.

100{
101 fBraggModel->Initialise(p, ptr);
102 fBBModel->Initialise(p, ptr);
103 SetParticle(p);
104
105 // always false before the run
106 SetDeexcitationFlag(false);
107
108 if(nullptr == fParticleChange) {
109 fParticleChange = GetParticleChangeForLoss();
110 if(UseAngularGeneratorFlag() && nullptr == GetAngularDistribution()) {
111 SetAngularDistribution(new G4DeltaAngle());
112 }
113 }
114 if(nullptr == lsdata) {
115 G4AutoLock l(&ionXSMutex);
116 if(nullptr == lsdata) {
117 isFirst = true;
118 lsdata = new G4LindhardSorensenData();
119 fIonData = new G4IonICRU73Data();
120 }
121 l.unlock();
122 }
123 if(isFirst) {
124 fIonData->Initialise();
125 }
126}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4VEmAngularDistribution * GetAngularDistribution()
void SetDeexcitationFlag(G4bool val)
void SetAngularDistribution(G4VEmAngularDistribution *)
G4bool UseAngularGeneratorFlag() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ MaxSecondaryEnergy()

G4double G4LindhardSorensenIonModel::MaxSecondaryEnergy ( const G4ParticleDefinition * pd,
G4double kinEnergy )
overrideprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 437 of file G4LindhardSorensenIonModel.cc.

439{
440 // here particle type is checked for any method
441 SetParticle(pd);
442 G4double tau = kinEnergy/mass;
443 return 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.) /
444 (1. + 2.0*(tau + 1.)*eRatio + eRatio*eRatio);
445}

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), CorrectionsAlongStep(), and SampleSecondaries().

◆ MinEnergyCut()

G4double G4LindhardSorensenIonModel::MinEnergyCut ( const G4ParticleDefinition * ,
const G4MaterialCutsCouple * couple )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 175 of file G4LindhardSorensenIonModel.cc.

177{
179}
G4double GetMeanExcitationEnergy() const
const G4Material * GetMaterial() const
G4IonisParamMat * GetIonisation() const

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 319 of file G4LindhardSorensenIonModel.cc.

325{
326 G4double kineticEnergy = dp->GetKineticEnergy();
327 // take into account formfactor
328 const G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(), kineticEnergy);
329 const G4double minKinEnergy = std::min(cut, tmax);
330 const G4double maxKinEnergy = std::min(maxEnergy, tmax);
331 if(minKinEnergy >= maxKinEnergy) { return; }
332
333 //G4cout << "G4LindhardSorensenIonModel::SampleSecondaries Emin= "
334 // << minKinEnergy << " Emax= " << maxKinEnergy << G4endl;
335
336 G4double totEnergy = kineticEnergy + mass;
337 G4double etot2 = totEnergy*totEnergy;
338 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
339
340 G4double deltaKinEnergy, f;
341 G4double f1 = 0.0;
342 G4double fmax = 1.0;
343 if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; }
344
345 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
346 G4double rndm[2];
347
348 // sampling without nuclear size effect
349 do {
350 rndmEngineMod->flatArray(2, rndm);
351 deltaKinEnergy = minKinEnergy*maxKinEnergy
352 /(minKinEnergy*(1.0 - rndm[0]) + maxKinEnergy*rndm[0]);
353
354 f = 1.0 - beta2*deltaKinEnergy/tmax;
355 if( 0.0 < spin ) {
356 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
357 f += f1;
358 }
359
360 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
361 } while( fmax*rndm[1] > f);
362
363 // projectile formfactor - suppresion of high energy
364 // delta-electron production at high energy
365
366 G4double x = formfact*deltaKinEnergy;
367 if(x > 1.e-6) {
368
369 G4double x1 = 1.0 + x;
370 G4double grej = 1.0/(x1*x1);
371 if( 0.0 < spin ) {
372 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
373 grej *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
374 }
375 if(grej > 1.1) {
376 G4cout << "### G4LindhardSorensenIonModel WARNING: grej= " << grej
377 << " " << dp->GetDefinition()->GetParticleName()
378 << " Ekin(MeV)= " << kineticEnergy
379 << " delEkin(MeV)= " << deltaKinEnergy
380 << G4endl;
381 }
382 if(rndmEngineMod->flat() > grej) { return; }
383 }
384
385 G4ThreeVector deltaDirection;
386
388
389 const G4Material* mat = couple->GetMaterial();
391
392 deltaDirection =
393 GetAngularDistribution()->SampleDirection(dp, deltaKinEnergy, Z, mat);
394
395 } else {
396
397 G4double deltaMomentum =
398 std::sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
399 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
400 (deltaMomentum * dp->GetTotalMomentum());
401 cost = std::min(cost, 1.0);
402 G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
403
404 G4double phi = CLHEP::twopi*rndmEngineMod->flat();
405
406 deltaDirection.set(sint*std::cos(phi),sint*std::sin(phi), cost) ;
407 deltaDirection.rotateUz(dp->GetMomentumDirection());
408 }
409 /*
410 G4cout << "### G4LindhardSorensenIonModel "
411 << dp->GetDefinition()->GetParticleName()
412 << " Ekin(MeV)= " << kineticEnergy
413 << " delEkin(MeV)= " << deltaKinEnergy
414 << " tmin(MeV)= " << minKinEnergy
415 << " tmax(MeV)= " << maxKinEnergy
416 << " dir= " << dp->GetMomentumDirection()
417 << " dirDelta= " << deltaDirection
418 << G4endl;
419 */
420 // create G4DynamicParticle object for delta ray
421 auto delta = new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
422
423 vdp->push_back(delta);
424
425 // Change kinematics of primary particle
426 kineticEnergy -= deltaKinEnergy;
427 G4ThreeVector finalP = dp->GetMomentum() - delta->GetMomentum();
428 finalP = finalP.unit();
429
430 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
431 fParticleChange->SetProposedMomentumDirection(finalP);
432}
CLHEP::Hep3Vector G4ThreeVector
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetTotalMomentum() const
const G4String & GetParticleName() const
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
G4int SelectRandomAtomNumber(const G4Material *) const

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