Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PAIModel Class Referencefinal

#include <G4PAIModel.hh>

+ Inheritance diagram for G4PAIModel:

Public Member Functions

 G4PAIModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PAI")
 
 ~G4PAIModel () final
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) final
 
G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) final
 
G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
 
G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double, const G4double, const G4double, const G4double) final
 
G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double, const G4double, const G4double) final
 
void DefineForRegion (const G4Region *r) final
 
G4PAIModelDataGetPAIModelData ()
 
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples ()
 
G4double ComputeMaxEnergy (G4double scaledEnergy)
 
void SetVerboseLevel (G4int verbose)
 
G4PAIModeloperator= (const G4PAIModel &right)=delete
 
 G4PAIModel (const G4PAIModel &)=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 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 void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
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
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual void InitialiseMe (const G4ParticleDefinition *)
 
virtual void SetParticleAndCharge (const G4ParticleDefinition *, G4double q2)
 
const G4StringGetName () const
 
G4VEmFluctuationModeloperator= (const G4VEmFluctuationModel &right)=delete
 
 G4VEmFluctuationModel (const G4VEmFluctuationModel &)=delete
 

Protected Member Functions

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy) final
 
- 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 68 of file G4PAIModel.hh.

Constructor & Destructor Documentation

◆ G4PAIModel() [1/2]

G4PAIModel::G4PAIModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "PAI" )
explicit

Definition at line 76 of file G4PAIModel.cc.

78 fVerbose(0),
79 fModelData(nullptr),
80 fParticle(nullptr)
81{
82 fElectron = G4Electron::Electron();
83 fPositron = G4Positron::Positron();
84
85 fParticleChange = nullptr;
86
87 if(p) { SetParticle(p); }
88 else { SetParticle(fElectron); }
89
90 // default generator
92 fLowestTcut = 12.5*CLHEP::eV;
93}
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Positron * Positron()
Definition G4Positron.cc:90
G4VEmFluctuationModel(const G4String &nam)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4PAIModel()

G4PAIModel::~G4PAIModel ( )
final

Definition at line 97 of file G4PAIModel.cc.

98{
99 if(IsMaster()) { delete fModelData; }
100}
G4bool IsMaster() const

◆ G4PAIModel() [2/2]

G4PAIModel::G4PAIModel ( const G4PAIModel & )
delete

Member Function Documentation

◆ ComputeDEDXPerVolume()

G4double G4PAIModel::ComputeDEDXPerVolume ( const G4Material * ,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 204 of file G4PAIModel.cc.

208{
209 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
210 if(0 > coupleIndex) { return 0.0; }
211
212 G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
213 G4double scaledTkin = kineticEnergy*fRatio;
214 G4double dedx = fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
215 return dedx;
216}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
const G4MaterialCutsCouple * CurrentCouple() const

◆ ComputeMaxEnergy()

G4double G4PAIModel::ComputeMaxEnergy ( G4double scaledEnergy)
inline

Definition at line 166 of file G4PAIModel.hh.

167{
168 return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
169}

Referenced by G4PAIModelData::Initialise().

◆ CrossSectionPerVolume()

G4double G4PAIModel::CrossSectionPerVolume ( const G4Material * ,
const G4ParticleDefinition * p,
G4double kineticEnergy,
G4double cutEnergy,
G4double maxEnergy )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 220 of file G4PAIModel.cc.

225{
226 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
227 if(0 > coupleIndex) { return 0.0; }
228
229 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
230 if(tmax <= cutEnergy) { return 0.0; }
231
232 G4double scaledTkin = kineticEnergy*fRatio;
233 G4double xs = fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
234 scaledTkin, cutEnergy, tmax);
235 return xs;
236}
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const

◆ DefineForRegion()

void G4PAIModel::DefineForRegion ( const G4Region * r)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 377 of file G4PAIModel.cc.

378{
379 fPAIRegionVector.push_back(r);
380}

◆ Dispersion()

G4double G4PAIModel::Dispersion ( const G4Material * material,
const G4DynamicParticle * aParticle,
const G4double tcut,
const G4double tmax,
const G4double step )
finalvirtual

Implements G4VEmFluctuationModel.

Definition at line 340 of file G4PAIModel.cc.

345{
346 G4double particleMass = aParticle->GetMass();
347 G4double electronDensity = material->GetElectronDensity();
348 G4double kineticEnergy = aParticle->GetKineticEnergy();
349 G4double q = aParticle->GetCharge()/eplus;
350 G4double etot = kineticEnergy + particleMass;
351 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
352 G4double siga = (tmax/beta2 - 0.5*tcut) * twopi_mc2_rcl2 * step
353 * electronDensity * q * q;
354
355 return siga;
356}
G4double GetMass() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetElectronDensity() const

◆ GetPAIModelData()

G4PAIModelData * G4PAIModel::GetPAIModelData ( )
inline

Definition at line 155 of file G4PAIModel.hh.

156{
157 return fModelData;
158}

Referenced by InitialiseLocal().

◆ GetVectorOfCouples()

const std::vector< const G4MaterialCutsCouple * > & G4PAIModel::GetVectorOfCouples ( )
inline

Definition at line 161 of file G4PAIModel.hh.

162{
163 return fMaterialCutsCoupleVector;
164}

Referenced by InitialiseLocal().

◆ Initialise()

void G4PAIModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector & cuts )
finalvirtual

Implements G4VEmModel.

Definition at line 104 of file G4PAIModel.cc.

106{
107 if(fVerbose > 1) {
108 G4cout<<"G4PAIModel::Initialise for "<<p->GetParticleName()<<G4endl;
109 }
110 SetParticle(p);
111 fParticleChange = GetParticleChangeForLoss();
112
113 if(IsMaster()) {
114
115 delete fModelData;
116 fMaterialCutsCoupleVector.clear();
117 if(fVerbose > 1) {
118 G4cout << "G4PAIModel instantiates data for " << p->GetParticleName()
119 << G4endl;
120 }
121 G4double tmin = LowEnergyLimit()*fRatio;
122 G4double tmax = HighEnergyLimit()*fRatio;
123 fModelData = new G4PAIModelData(tmin, tmax, fVerbose);
124
125 // Prepare initialization
126 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
127 std::size_t numOfMat = G4Material::GetNumberOfMaterials();
128 std::size_t numRegions = fPAIRegionVector.size();
129
130 // protect for unit tests
131 if(0 == numRegions) {
132 G4Exception("G4PAIModel::Initialise()","em0106",JustWarning,
133 "no G4Regions are registered for the PAI model - World is used");
134 fPAIRegionVector.push_back(G4RegionStore::GetInstance()
135 ->GetRegion("DefaultRegionForTheWorld", false));
136 numRegions = 1;
137 }
138
139 if(fVerbose > 1) {
140 G4cout << "G4PAIModel is defined for " << numRegions << " regions "
141 << "; number of materials " << numOfMat << G4endl;
142 }
143 for(std::size_t iReg = 0; iReg<numRegions; ++iReg) {
144 const G4Region* curReg = fPAIRegionVector[iReg];
145 G4Region* reg = const_cast<G4Region*>(curReg);
146
147 for(std::size_t jMat = 0; jMat<numOfMat; ++jMat) {
148 G4Material* mat = (*theMaterialTable)[jMat];
149 const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
150 std::size_t n = fMaterialCutsCoupleVector.size();
151 if(nullptr != cutCouple) {
152 if(fVerbose > 1) {
153 G4cout << "Region <" << curReg->GetName() << "> mat <"
154 << mat->GetName() << "> CoupleIndex= "
155 << cutCouple->GetIndex()
156 << " " << p->GetParticleName()
157 << " cutsize= " << cuts.size() << G4endl;
158 }
159 // check if this couple is not already initialized
160 G4bool isnew = true;
161 if(0 < n) {
162 for(std::size_t i=0; i<n; ++i) {
163 G4cout << i << G4endl;
164 if(cutCouple == fMaterialCutsCoupleVector[i]) {
165 isnew = false;
166 break;
167 }
168 }
169 }
170 // initialise data banks
171 // G4cout << " isNew: " << isnew << " " << cutCouple << G4endl;
172 if(isnew) {
173 fMaterialCutsCoupleVector.push_back(cutCouple);
174 fModelData->Initialise(cutCouple, this);
175 }
176 }
177 }
178 }
180 }
181}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4Material * > G4MaterialTable
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const
void Initialise(const G4MaterialCutsCouple *, G4PAIModel *)
const G4String & GetParticleName() const
static G4RegionStore * GetInstance()
G4MaterialCutsCouple * FindCouple(G4Material *mat)
const G4String & GetName() const
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForLoss * GetParticleChangeForLoss()

◆ InitialiseLocal()

void G4PAIModel::InitialiseLocal ( const G4ParticleDefinition * ,
G4VEmModel * masterModel )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 185 of file G4PAIModel.cc.

187{
188 fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
189 fMaterialCutsCoupleVector =
190 static_cast<G4PAIModel*>(masterModel)->GetVectorOfCouples();
192}
G4PAIModelData * GetPAIModelData()
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ MaxSecondaryEnergy()

G4double G4PAIModel::MaxSecondaryEnergy ( const G4ParticleDefinition * p,
G4double kinEnergy )
finalprotectedvirtual

Reimplemented from G4VEmModel.

Definition at line 360 of file G4PAIModel.cc.

362{
363 SetParticle(p);
364 G4double tmax = kinEnergy;
365 if(p == fElectron) { tmax *= 0.5; }
366 else if(p != fPositron) {
367 G4double ratio= electron_mass_c2/fMass;
368 G4double gamma= kinEnergy/fMass + 1.0;
369 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
370 (1. + 2.0*gamma*ratio + ratio*ratio);
371 }
372 return tmax;
373}

Referenced by ComputeDEDXPerVolume(), ComputeMaxEnergy(), CrossSectionPerVolume(), and SampleSecondaries().

◆ MinEnergyCut()

G4double G4PAIModel::MinEnergyCut ( const G4ParticleDefinition * ,
const G4MaterialCutsCouple * couple )
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 196 of file G4PAIModel.cc.

198{
199 return fLowestTcut;
200}

◆ operator=()

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

◆ SampleFluctuations()

G4double G4PAIModel::SampleFluctuations ( const G4MaterialCutsCouple * matCC,
const G4DynamicParticle * aParticle,
const G4double tcut,
const G4double ,
const G4double step,
const G4double eloss )
finalvirtual

Implements G4VEmFluctuationModel.

Definition at line 303 of file G4PAIModel.cc.

309{
310 G4int coupleIndex = FindCoupleIndex(matCC);
311 if(0 > coupleIndex) { return eloss; }
312
313 SetParticle(aParticle->GetDefinition());
314
315 /*
316 G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
317 << " Eloss(keV)= " << eloss/keV << " in "
318 << matCC->Getmaterial()->GetName() << G4endl;
319 */
320
321 G4double Tkin = aParticle->GetKineticEnergy();
322 G4double scaledTkin = Tkin*fRatio;
323
324 G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
325 scaledTkin, tcut,
326 step*fChargeSquare);
327
328 // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
329 //<<step/mm<<" mm"<<G4endl;
330 return loss;
331
332}
G4ParticleDefinition * GetDefinition() const
G4double SampleAlongStepTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double tmax, G4double stepFactor) const

◆ SampleSecondaries()

void G4PAIModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * vdp,
const G4MaterialCutsCouple * matCC,
const G4DynamicParticle * dp,
G4double tmin,
G4double maxEnergy )
finalvirtual

Implements G4VEmModel.

Definition at line 243 of file G4PAIModel.cc.

248{
249 G4int coupleIndex = FindCoupleIndex(matCC);
250
251 //G4cout << "G4PAIModel::SampleSecondaries: coupleIndex= "<<coupleIndex<<G4endl;
252 if(0 > coupleIndex) { return; }
253
254 SetParticle(dp->GetDefinition());
255 G4double kineticEnergy = dp->GetKineticEnergy();
256
257 G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
258 if(maxEnergy < tmax) { tmax = maxEnergy; }
259 if(tmin >= tmax) { return; }
260
261 G4ThreeVector direction= dp->GetMomentumDirection();
262 G4double scaledTkin = kineticEnergy*fRatio;
263 G4double totalEnergy = kineticEnergy + fMass;
264 G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
265
266 G4double deltaTkin =
267 fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin, tmin, tmax);
268
269 //G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
270 // <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
271
272 if( !(deltaTkin <= 0.) && !(deltaTkin > 0)) {
273 G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
274 <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
275 return;
276 }
277 if( deltaTkin <= 0.) { return; }
278
279 if( deltaTkin > tmax) { deltaTkin = tmax; }
280
281 const G4Element* anElement = SelectTargetAtom(matCC, fParticle, kineticEnergy,
282 dp->GetLogKineticEnergy());
283
284 G4int Z = anElement->GetZasInt();
285
286 auto deltaRay = new G4DynamicParticle(fElectron,
287 GetAngularDistribution()->SampleDirection(dp, deltaTkin,
288 Z, matCC->GetMaterial()),
289 deltaTkin);
290
291 // primary change
292 kineticEnergy -= deltaTkin;
293 G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
294 direction = dir.unit();
295 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
296 fParticleChange->SetProposedMomentumDirection(direction);
297
298 vdp->push_back(deltaRay);
299}
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
G4double SamplePostStepTransfer(G4int coupleIndex, G4double scaledTkin, G4double tmin, G4double tmax) const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedMomentumDirection(const G4ThreeVector &dir)
G4VEmAngularDistribution * GetAngularDistribution()
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)

◆ SetVerboseLevel()

void G4PAIModel::SetVerboseLevel ( G4int verbose)
inline

Definition at line 171 of file G4PAIModel.hh.

172{
173 fVerbose=verbose;
174}

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