Geant4 10.7.0
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")
 
virtual ~G4BraggIonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
 
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) override
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *p, const G4Material *mat, G4double kineticEnergy) override
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) override
 
- 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 *, 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 * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

virtual 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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 68 of file G4BraggIonModel.hh.

Constructor & Destructor Documentation

◆ G4BraggIonModel()

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

Definition at line 79 of file G4BraggIonModel.cc.

81 : G4VEmModel(nam),
82 corr(nullptr),
83 particle(nullptr),
84 fParticleChange(nullptr),
85 fICRU90(nullptr),
86 currentMaterial(nullptr),
87 baseMaterial(nullptr),
88 iMolecula(-1),
89 iASTAR(-1),
90 iICRU90(-1),
91 isIon(false)
92{
93 SetHighEnergyLimit(2.0*MeV);
94
95 HeMass = 3.727417*GeV;
96 rateMassHe2p = HeMass/proton_mass_c2;
97 lowestKinEnergy = 1.0*keV/rateMassHe2p;
98 massFactor = 1000.*amu_c2/HeMass;
99 theZieglerFactor = eV*cm2*1.0e-15;
100 theElectron = G4Electron::Electron();
101 corrFactor = 1.0;
102 if(p) { SetParticle(p); }
103 else { SetParticle(theElectron); }
104}
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757

◆ ~G4BraggIonModel()

G4BraggIonModel::~G4BraggIonModel ( )
virtual

Definition at line 108 of file G4BraggIonModel.cc.

109{
110 if(IsMaster()) { delete fASTAR; fASTAR = nullptr; }
111}
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

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 216 of file G4BraggIonModel.cc.

222{
223 return
224 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
225}
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)

◆ ComputeCrossSectionPerElectron()

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

Definition at line 187 of file G4BraggIonModel.cc.

192{
193 G4double cross = 0.0;
194 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
195 G4double maxEnergy = std::min(tmax,maxKinEnergy);
196 if(cutEnergy < tmax) {
197
198 G4double energy = kineticEnergy + mass;
199 G4double energy2 = energy*energy;
200 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
201 cross = (maxEnergy - cutEnergy)/(cutEnergy*maxEnergy)
202 - beta2*G4Log(maxEnergy/cutEnergy)/tmax;
203
204 if( 0.0 < spin ) { cross += 0.5*(maxEnergy - cutEnergy)/energy2; }
205
206 cross *= twopi_mc2_rcl2*chargeSquare/beta2;
207 }
208 // G4cout << "BR: e= " << kineticEnergy << " tmin= " << cutEnergy
209 // << " tmax= " << tmax << " cross= " << cross << G4endl;
210
211 return cross;
212}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
G4double energy(const ThreeVector &p, const G4double m)

Referenced by ComputeCrossSectionPerAtom(), and CrossSectionPerVolume().

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4BraggNoDeltaModel.

Definition at line 242 of file G4BraggIonModel.cc.

246{
247 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
248 G4double tmin = min(cutEnergy, tmax);
249 G4double tkin = kineticEnergy/massRate;
250 G4double dedx = 0.0;
251
252 if(tkin < lowestKinEnergy) {
253 dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy);
254 } else {
255 dedx = DEDX(material, tkin);
256 }
257
258 if (cutEnergy < tmax) {
259
260 G4double tau = kineticEnergy/mass;
261 G4double gam = tau + 1.0;
262 G4double bg2 = tau * (tau+2.0);
263 G4double beta2 = bg2/(gam*gam);
264 G4double x = tmin/tmax;
265
266 dedx += (G4Log(x) + (1.0 - x)*beta2) * twopi_mc2_rcl2
267 * (material->GetElectronDensity())/beta2;
268 }
269
270 // now compute the total ionization loss
271
272 dedx = std::max(dedx, 0.0);
273 dedx *= chargeSquare;
274 /*
275 G4cout << "Bragg: tkin(MeV) = " << tkin/MeV << " dedx(MeVxcm^2/g) = "
276 << dedx*gram/(MeV*cm2*material->GetDensity())
277 << " q2 = " << chargeSquare << G4endl;
278 */
279 return dedx;
280}
G4double GetElectronDensity() const
Definition: G4Material.hh:215
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

Referenced by G4BraggNoDeltaModel::ComputeDEDXPerVolume().

◆ CorrectionsAlongStep()

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

Reimplemented from G4VEmModel.

Definition at line 284 of file G4BraggIonModel.cc.

289{
290 // this method is called only for ions
291 const G4ParticleDefinition* p = dp->GetDefinition();
292 const G4Material* mat = couple->GetMaterial();
293 G4double preKinEnergy = dp->GetKineticEnergy();
294 G4double e = preKinEnergy - eloss*0.5;
295 if(e < 0.0) { e = preKinEnergy*0.5; }
296
297 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e);
299 G4double qfactor = q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor;
300 eloss *= qfactor;
301
302 //G4cout << "G4BraggIonModel::CorrectionsAlongStep e= " << e
303 // << " qfactor= " << qfactor << " " << p->GetParticleName() <<G4endl;
304}
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double EffectiveChargeCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4Material * GetMaterial() const
virtual void SetParticleAndCharge(const G4ParticleDefinition *, G4double q2)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:604

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Reimplemented in G4BraggNoDeltaModel.

Definition at line 229 of file G4BraggIonModel.cc.

235{
236 return material->GetElectronDensity()
237 *ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
238}

◆ GetChargeSquareRatio()

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

Reimplemented from G4VEmModel.

Definition at line 162 of file G4BraggIonModel.cc.

165{
166 //G4cout<<"G4BraggIonModel::GetChargeSquareRatio e= "<<kineticEnergy<<G4endl;
167 // this method is called only for ions
168 G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy);
169 corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy);
170 return corrFactor;
171}

◆ GetParticleCharge()

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

Reimplemented from G4VEmModel.

Definition at line 175 of file G4BraggIonModel.cc.

178{
179 //G4cout<<"G4BraggIonModel::GetParticleCharge e= "<<kineticEnergy <<
180 // " q= " << corr->GetParticleCharge(p,mat,kineticEnergy) <<G4endl;
181 // this method is called only for ions
182 return corr->GetParticleCharge(p,mat,kineticEnergy);
183}
G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 115 of file G4BraggIonModel.cc.

117{
118 if(p != particle) { SetParticle(p); }
119
120 corrFactor = chargeSquare;
121
122 // always false before the run
123 SetDeexcitationFlag(false);
124
125 if(IsMaster()) {
126 if(nullptr == fASTAR) { fASTAR = new G4ASTARStopping(); }
127 if(particle->GetPDGMass() < GeV) { fASTAR->Initialise(); }
128 if(G4EmParameters::Instance()->UseICRU90Data()) {
129 if(!fICRU90) {
131 } else if(particle->GetPDGMass() < GeV) { fICRU90->Initialise(); }
132 }
133 }
134
135 if(nullptr == fParticleChange) {
136
139 }
140 G4String pname = particle->GetParticleName();
141 if(particle->GetParticleType() == "nucleus" &&
142 pname != "deuteron" && pname != "triton" &&
143 pname != "alpha+" && pname != "helium" &&
144 pname != "hydrogen") { isIon = true; }
145
147
148 fParticleChange = GetParticleChangeForLoss();
149 }
150}
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
G4ICRU90StoppingData * GetICRU90StoppingData()
static G4NistManager * Instance()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
G4bool UseAngularGeneratorFlag() const
Definition: G4VEmModel.hh:708
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 386 of file G4BraggIonModel.cc.

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

Referenced by ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 154 of file G4BraggIonModel.cc.

156{
158}
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 308 of file G4BraggIonModel.cc.

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

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