Geant4 11.1.1
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)
 
- 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
 
- Public Member Functions inherited from G4VEmFluctuationModel
 G4VEmFluctuationModel (const G4String &nam)
 
virtual ~G4VEmFluctuationModel ()
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length)=0
 
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
 
G4PAIModeloperator= (const G4PAIModel &right)=delete
 
 G4PAIModel (const G4PAIModel &)=delete
 
- 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 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:93
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4PAIModel()

G4PAIModel::~G4PAIModel ( )
final

Definition at line 97 of file G4PAIModel.cc.

98{
99 if(IsMaster()) { delete fModelData; }
100}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4PAIModel() [2/2]

G4PAIModel::G4PAIModel ( const G4PAIModel )
protecteddelete

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 212 of file G4PAIModel.cc.

216{
217 //G4cout << "===1=== " << CurrentCouple()
218 // << " idx= " << CurrentCouple()->GetIndex()
219 // << " " << fMaterialCutsCoupleVector[0]
220 // << G4endl;
221 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
222 //G4cout << "===2=== " << coupleIndex << G4endl;
223 if(0 > coupleIndex) { return 0.0; }
224
225 G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
226
227 G4double scaledTkin = kineticEnergy*fRatio;
228
229 return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin,
230 cut);
231}
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
Definition: G4PAIModel.cc:382
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:486

◆ 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 235 of file G4PAIModel.cc.

240{
241 //G4cout << "===3=== " << CurrentCouple()
242 // << " idx= " << CurrentCouple()->GetIndex()
243 // << " " << fMaterialCutsCoupleVector[0]
244 // << G4endl;
245 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
246 //G4cout << "===4=== " << coupleIndex << G4endl;
247 if(0 > coupleIndex) { return 0.0; }
248
249 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
250 if(tmax <= cutEnergy) { return 0.0; }
251
252 G4double scaledTkin = kineticEnergy*fRatio;
253
254 return fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
255 scaledTkin,
256 cutEnergy,
257 tmax);
258}
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 399 of file G4PAIModel.cc.

400{
401 fPAIRegionVector.push_back(r);
402}

◆ 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 362 of file G4PAIModel.cc.

367{
368 G4double particleMass = aParticle->GetMass();
369 G4double electronDensity = material->GetElectronDensity();
370 G4double kineticEnergy = aParticle->GetKineticEnergy();
371 G4double q = aParticle->GetCharge()/eplus;
372 G4double etot = kineticEnergy + particleMass;
373 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
374 G4double siga = (tmax/beta2 - 0.5*tcut) * twopi_mc2_rcl2 * step
375 * electronDensity * q * q;
376
377 return siga;
378}
G4double GetMass() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:212

◆ 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 /*
152 G4cout << "Region: " << reg->GetName() << " " << reg
153 << " Couple " << cutCouple
154 << " PAI defined for " << n << " couples"
155 << " jMat= " << jMat << " " << mat->GetName()
156 << G4endl;
157 */
158 if(nullptr != cutCouple) {
159 if(fVerbose > 1) {
160 G4cout << "Region <" << curReg->GetName() << "> mat <"
161 << mat->GetName() << "> CoupleIndex= "
162 << cutCouple->GetIndex()
163 << " " << p->GetParticleName()
164 << " cutsize= " << cuts.size() << G4endl;
165 }
166 // check if this couple is not already initialized
167 G4bool isnew = true;
168 if(0 < n) {
169 for(std::size_t i=0; i<n; ++i) {
170 G4cout << i << G4endl;
171 if(cutCouple == fMaterialCutsCoupleVector[i]) {
172 isnew = false;
173 break;
174 }
175 }
176 }
177 // initialise data banks
178 // G4cout << " isNew: " << isnew << " " << cutCouple << G4endl;
179 if(isnew) {
180 fMaterialCutsCoupleVector.push_back(cutCouple);
181 fModelData->Initialise(cutCouple, this);
182 }
183 }
184 }
185 }
187 }
188}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< G4Material * > G4MaterialTable
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
void Initialise(const G4MaterialCutsCouple *, G4PAIModel *)
const G4String & GetParticleName() const
static G4RegionStore * GetInstance()
G4MaterialCutsCouple * FindCouple(G4Material *mat)
const G4String & GetName() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:109

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 192 of file G4PAIModel.cc.

194{
195 SetParticle(p);
196 fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
197 fMaterialCutsCoupleVector =
198 static_cast<G4PAIModel*>(masterModel)->GetVectorOfCouples();
200}
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:155
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
Definition: G4PAIModel.hh:161
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 382 of file G4PAIModel.cc.

384{
385 SetParticle(p);
386 G4double tmax = kinEnergy;
387 if(p == fElectron) { tmax *= 0.5; }
388 else if(p != fPositron) {
389 G4double ratio= electron_mass_c2/fMass;
390 G4double gamma= kinEnergy/fMass + 1.0;
391 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
392 (1. + 2.0*gamma*ratio + ratio*ratio);
393 }
394 return tmax;
395}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 204 of file G4PAIModel.cc.

206{
207 return fLowestTcut;
208}

◆ operator=()

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

◆ 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 325 of file G4PAIModel.cc.

331{
332 G4int coupleIndex = FindCoupleIndex(matCC);
333 if(0 > coupleIndex) { return eloss; }
334
335 SetParticle(aParticle->GetDefinition());
336
337 /*
338 G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
339 << " Eloss(keV)= " << eloss/keV << " in "
340 << matCC->Getmaterial()->GetName() << G4endl;
341 */
342
343 G4double Tkin = aParticle->GetKineticEnergy();
344 G4double scaledTkin = Tkin*fRatio;
345
346 G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
347 scaledTkin, tcut,
348 step*fChargeSquare);
349
350 // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
351 //<<step/mm<<" mm"<<G4endl;
352 return loss;
353
354}
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 265 of file G4PAIModel.cc.

270{
271 G4int coupleIndex = FindCoupleIndex(matCC);
272
273 //G4cout << "G4PAIModel::SampleSecondaries: coupleIndex= "<<coupleIndex<<G4endl;
274 if(0 > coupleIndex) { return; }
275
276 SetParticle(dp->GetDefinition());
277 G4double kineticEnergy = dp->GetKineticEnergy();
278
279 G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
280 if(maxEnergy < tmax) { tmax = maxEnergy; }
281 if(tmin >= tmax) { return; }
282
283 G4ThreeVector direction= dp->GetMomentumDirection();
284 G4double scaledTkin = kineticEnergy*fRatio;
285 G4double totalEnergy = kineticEnergy + fMass;
286 G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
287
288 G4double deltaTkin =
289 fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin, tmin, tmax);
290
291 //G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
292 // <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
293
294 if( !(deltaTkin <= 0.) && !(deltaTkin > 0)) {
295 G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
296 <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
297 return;
298 }
299 if( deltaTkin <= 0.) { return; }
300
301 if( deltaTkin > tmax) { deltaTkin = tmax; }
302
303 const G4Element* anElement = SelectTargetAtom(matCC, fParticle, kineticEnergy,
304 dp->GetLogKineticEnergy());
305
306 G4int Z = G4lrint(anElement->GetZ());
307
308 auto deltaRay = new G4DynamicParticle(fElectron,
309 GetAngularDistribution()->SampleDirection(dp, deltaTkin,
310 Z, matCC->GetMaterial()),
311 deltaTkin);
312
313 // primary change
314 kineticEnergy -= deltaTkin;
315 G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
316 direction = dir.unit();
317 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
318 fParticleChange->SetProposedMomentumDirection(direction);
319
320 vdp->push_back(deltaRay);
321}
const G4int Z[17]
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetLogKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
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()
Definition: G4VEmModel.hh:600
const G4Element * SelectTargetAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:577
int G4lrint(double ad)
Definition: templates.hh:134

◆ 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: