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

#include <G4PAIPhotModel.hh>

+ Inheritance diagram for G4PAIPhotModel:

Public Member Functions

 G4PAIPhotModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PAI")
 
 ~G4PAIPhotModel () 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
 
G4PAIPhotDataGetPAIPhotData ()
 
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples ()
 
G4double ComputeMaxEnergy (G4double scaledEnergy)
 
void SetVerboseLevel (G4int verbose)
 
G4PAIPhotModeloperator= (const G4PAIPhotModel &right)=delete
 
 G4PAIPhotModel (const G4PAIPhotModel &)=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
 
- 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
 
- 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 63 of file G4PAIPhotModel.hh.

Constructor & Destructor Documentation

◆ G4PAIPhotModel() [1/2]

G4PAIPhotModel::G4PAIPhotModel ( const G4ParticleDefinition p = nullptr,
const G4String nam = "PAI" 
)
explicit

Definition at line 72 of file G4PAIPhotModel.cc.

74 fVerbose(0),
75 fModelData(nullptr),
76 fParticle(nullptr)
77{
78 fElectron = G4Electron::Electron();
79 fPositron = G4Positron::Positron();
80
81 fParticleChange = nullptr;
82
83 if(p) { SetParticle(p); }
84 else { SetParticle(fElectron); }
85
86 // default generator
88 fLowestTcut = 12.5*CLHEP::eV;
89}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4PAIPhotModel()

G4PAIPhotModel::~G4PAIPhotModel ( )
final

Definition at line 93 of file G4PAIPhotModel.cc.

94{
95 //G4cout << "G4PAIPhotModel::~G4PAIPhotModel() " << this << G4endl;
96 if(IsMaster()) { delete fModelData; fModelData = nullptr; }
97}
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4PAIPhotModel() [2/2]

G4PAIPhotModel::G4PAIPhotModel ( const G4PAIPhotModel )
delete

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 204 of file G4PAIPhotModel.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
214 G4double scaledTkin = kineticEnergy*fRatio;
215
216 return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
217}
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
Definition: G4VEmModel.hh:486

◆ ComputeMaxEnergy()

G4double G4PAIPhotModel::ComputeMaxEnergy ( G4double  scaledEnergy)
inline

Definition at line 160 of file G4PAIPhotModel.hh.

161{
162 return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
163}

Referenced by G4PAIPhotData::Initialise().

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 221 of file G4PAIPhotModel.cc.

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

◆ DefineForRegion()

void G4PAIPhotModel::DefineForRegion ( const G4Region r)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 444 of file G4PAIPhotModel.cc.

445{
446 fPAIRegionVector.push_back(r);
447}

◆ Dispersion()

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

Implements G4VEmFluctuationModel.

Definition at line 407 of file G4PAIPhotModel.cc.

412{
413 G4double particleMass = aParticle->GetMass();
414 G4double electronDensity = material->GetElectronDensity();
415 G4double kineticEnergy = aParticle->GetKineticEnergy();
416 G4double q = aParticle->GetCharge()/eplus;
417 G4double etot = kineticEnergy + particleMass;
418 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
419 G4double siga = (tmax/beta2 - 0.5*tcut) * twopi_mc2_rcl2 * step
420 * electronDensity * q * q;
421
422 return siga;
423}
G4double GetMass() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetElectronDensity() const
Definition: G4Material.hh:212

◆ GetPAIPhotData()

G4PAIPhotData * G4PAIPhotModel::GetPAIPhotData ( )
inline

Definition at line 149 of file G4PAIPhotModel.hh.

150{
151 return fModelData;
152}

Referenced by InitialiseLocal().

◆ GetVectorOfCouples()

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

Definition at line 155 of file G4PAIPhotModel.hh.

156{
157 return fMaterialCutsCoupleVector;
158}

Referenced by InitialiseLocal().

◆ Initialise()

void G4PAIPhotModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
finalvirtual

Implements G4VEmModel.

Definition at line 101 of file G4PAIPhotModel.cc.

103{
104 if(fVerbose > 0)
105 {
106 G4cout<<"G4PAIPhotModel::Initialise for "<<p->GetParticleName()<<G4endl;
107 }
108 SetParticle(p);
109 fParticleChange = GetParticleChangeForLoss();
110
111 if( IsMaster() )
112 {
113
115
116 delete fModelData;
117 fMaterialCutsCoupleVector.clear();
118
119 G4double tmin = LowEnergyLimit()*fRatio;
120 G4double tmax = HighEnergyLimit()*fRatio;
121 fModelData = new G4PAIPhotData(tmin, tmax, fVerbose);
122
123 // Prepare initialization
124 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
125 size_t numOfMat = G4Material::GetNumberOfMaterials();
126 size_t numRegions = fPAIRegionVector.size();
127
128 // protect for unit tests
129 if(0 == numRegions) {
130 G4Exception("G4PAIModel::Initialise()","em0106",JustWarning,
131 "no G4Regions are registered for the PAI model - World is used");
132 fPAIRegionVector.push_back(G4RegionStore::GetInstance()
133 ->GetRegion("DefaultRegionForTheWorld", false));
134 numRegions = 1;
135 }
136
137 for( size_t iReg = 0; iReg < numRegions; ++iReg )
138 {
139 const G4Region* curReg = fPAIRegionVector[iReg];
140 G4Region* reg = const_cast<G4Region*>(curReg);
141
142 for(size_t jMat = 0; jMat < numOfMat; ++jMat)
143 {
144 G4Material* mat = (*theMaterialTable)[jMat];
145 const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
146 //G4cout << "Couple <" << fCutCouple << G4endl;
147 if(cutCouple)
148 {
149 if(fVerbose>0)
150 {
151 G4cout << "Reg <" <<curReg->GetName() << "> mat <"
152 << mat->GetName() << "> fCouple= "
153 << cutCouple << ", idx= " << cutCouple->GetIndex()
154 <<" " << p->GetParticleName()
155 <<", cuts.size() = " << cuts.size() << G4endl;
156 }
157 // check if this couple is not already initialized
158
159 size_t n = fMaterialCutsCoupleVector.size();
160
161 G4bool isnew = true;
162 if( 0 < n )
163 {
164 for(size_t i=0; i<fMaterialCutsCoupleVector.size(); ++i)
165 {
166 if(cutCouple == fMaterialCutsCoupleVector[i]) {
167 isnew = false;
168 break;
169 }
170 }
171 }
172 // initialise data banks
173 if(isnew) {
174 fMaterialCutsCoupleVector.push_back(cutCouple);
175 G4double deltaCutInKinEnergy = cuts[cutCouple->GetIndex()];
176 fModelData->Initialise(cutCouple, deltaCutInKinEnergy, this);
177 }
178 }
179 }
180 }
181 }
182}
@ 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 *, G4double cut, G4PAIPhotModel *)
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 G4PAIPhotModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 186 of file G4PAIPhotModel.cc.

188{
189 fModelData = static_cast<G4PAIPhotModel*>(masterModel)->GetPAIPhotData();
190 fMaterialCutsCoupleVector = static_cast<G4PAIPhotModel*>(masterModel)->GetVectorOfCouples();
192}
G4PAIPhotData * GetPAIPhotData()
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 427 of file G4PAIPhotModel.cc.

429{
430 SetParticle(p);
431 G4double tmax = kinEnergy;
432 if(p == fElectron) { tmax *= 0.5; }
433 else if(p != fPositron) {
434 G4double ratio= electron_mass_c2/fMass;
435 G4double gamma= kinEnergy/fMass + 1.0;
436 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
437 (1. + 2.0*gamma*ratio + ratio*ratio);
438 }
439 return tmax;
440}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4PAIPhotModel.cc.

198{
199 return fLowestTcut;
200}

◆ operator=()

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

◆ SampleFluctuations()

G4double G4PAIPhotModel::SampleFluctuations ( const G4MaterialCutsCouple matCC,
const G4DynamicParticle aParticle,
const  G4double,
const  G4double,
const G4double  step,
const G4double  eloss 
)
finalvirtual

Implements G4VEmFluctuationModel.

Definition at line 366 of file G4PAIPhotModel.cc.

371{
372 // return 0.;
373 G4int coupleIndex = FindCoupleIndex(matCC);
374 if(0 > coupleIndex) { return eloss; }
375
376 SetParticle(aParticle->GetDefinition());
377
378
379 // G4cout << "G4PAIPhotModel::SampleFluctuations step(mm)= "<< step/mm
380 // << " Eloss(keV)= " << eloss/keV << " in "
381 // << matCC->GetMaterial()->GetName() << G4endl;
382
383
384 G4double Tkin = aParticle->GetKineticEnergy();
385 G4double scaledTkin = Tkin*fRatio;
386
387 G4double loss = fModelData->SampleAlongStepPhotonTransfer(coupleIndex, Tkin,
388 scaledTkin,
389 step*fChargeSquare);
390 loss += fModelData->SampleAlongStepPlasmonTransfer(coupleIndex, Tkin,
391 scaledTkin,
392 step*fChargeSquare);
393
394
395 // G4cout<<" PAIPhotModel::SampleFluctuations loss = "<<loss/keV<<" keV, on step = "
396 // <<step/mm<<" mm"<<G4endl;
397 return loss;
398
399}
G4ParticleDefinition * GetDefinition() const
G4double SampleAlongStepPlasmonTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double stepFactor) const
G4double SampleAlongStepPhotonTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double stepFactor) const

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 247 of file G4PAIPhotModel.cc.

252{
253 G4int coupleIndex = FindCoupleIndex(matCC);
254 if(0 > coupleIndex) { return; }
255
256 SetParticle(dp->GetDefinition());
257
258 G4double kineticEnergy = dp->GetKineticEnergy();
259
260 G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
261
262 if( maxEnergy < tmax) tmax = maxEnergy;
263 if( tmin >= tmax) return;
264
265 G4ThreeVector direction = dp->GetMomentumDirection();
266 G4double scaledTkin = kineticEnergy*fRatio;
267 G4double totalEnergy = kineticEnergy + fMass;
268 G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy + fMass));
269 G4double plRatio = fModelData->GetPlasmonRatio(coupleIndex, scaledTkin);
270
271 if( G4UniformRand() <= plRatio )
272 {
273 G4double deltaTkin = fModelData->SamplePostStepPlasmonTransfer(coupleIndex, scaledTkin);
274
275 // G4cout<<"G4PAIPhotModel::SampleSecondaries; dp "<<dp->GetParticleDefinition()->GetParticleName()
276 // <<"; Tkin = "<<kineticEnergy/keV<<" keV; transfer = "<<deltaTkin/keV<<" keV "<<G4endl;
277
278 if( deltaTkin <= 0. && fVerbose > 0)
279 {
280 G4cout<<"G4PAIPhotModel::SampleSecondary e- deltaTkin = "<<deltaTkin<<G4endl;
281 }
282 if( deltaTkin <= 0.) return;
283
284 if( deltaTkin > tmax) deltaTkin = tmax;
285
286 const G4Element* anElement = SelectTargetAtom(matCC,fParticle,kineticEnergy,
287 dp->GetLogKineticEnergy());
288 G4int Z = G4lrint(anElement->GetZ());
289
290 auto deltaRay = new G4DynamicParticle(fElectron,
291 GetAngularDistribution()->SampleDirection(dp, deltaTkin,
292 Z, matCC->GetMaterial()),
293 deltaTkin);
294
295 // primary change
296
297 kineticEnergy -= deltaTkin;
298
299 if( kineticEnergy <= 0. ) // kill primary as local? energy deposition
300 {
301 fParticleChange->SetProposedKineticEnergy(0.0);
302 fParticleChange->ProposeLocalEnergyDeposit(kineticEnergy+deltaTkin);
303 // fParticleChange->ProposeTrackStatus(fStopAndKill);
304 return;
305 }
306 else
307 {
308 G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
309 direction = dir.unit();
310 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
311 fParticleChange->SetProposedMomentumDirection(direction);
312 vdp->push_back(deltaRay);
313 }
314 }
315 else // secondary X-ray CR photon
316 {
317 G4double deltaTkin = fModelData->SamplePostStepPhotonTransfer(coupleIndex, scaledTkin);
318
319 // G4cout<<"PAIPhotonModel PhotonPostStepTransfer = "<<deltaTkin/keV<<" keV"<<G4endl;
320
321 if( deltaTkin <= 0. )
322 {
323 G4cout<<"G4PAIPhotonModel::SampleSecondary gamma deltaTkin = "<<deltaTkin<<G4endl;
324 }
325 if( deltaTkin <= 0.) return;
326
327 if( deltaTkin >= kineticEnergy ) // stop primary
328 {
329 deltaTkin = kineticEnergy;
330 kineticEnergy = 0.0;
331 }
332 G4double costheta = 0.; // G4UniformRand(); // VG: ??? for start only
333 G4double sintheta = sqrt((1.+costheta)*(1.-costheta));
334
335 // direction of the 'Cherenkov' photon
336 G4double phi = twopi*G4UniformRand();
337 G4double dirx = sintheta*cos(phi), diry = sintheta*sin(phi), dirz = costheta;
338
339 G4ThreeVector deltaDirection(dirx,diry,dirz);
340 deltaDirection.rotateUz(direction);
341
342 if( kineticEnergy > 0.) // primary change
343 {
344 kineticEnergy -= deltaTkin;
345 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
346 }
347 else // stop primary, but pass X-ray CR
348 {
349 // fParticleChange->ProposeLocalEnergyDeposit(deltaTkin);
350 fParticleChange->SetProposedKineticEnergy(0.0);
351 }
352 // create G4DynamicParticle object for photon ray
353
354 auto photonRay = new G4DynamicParticle;
355 photonRay->SetDefinition( G4Gamma::Gamma() );
356 photonRay->SetKineticEnergy( deltaTkin );
357 photonRay->SetMomentumDirection(deltaDirection);
358
359 vdp->push_back(photonRay);
360 }
361 return;
362}
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double GetLogKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4Material * GetMaterial() const
G4double SamplePostStepPlasmonTransfer(G4int coupleIndex, G4double scaledTkin) const
G4double GetPlasmonRatio(G4int coupleIndex, G4double scaledTkin) const
G4double SamplePostStepPhotonTransfer(G4int coupleIndex, G4double scaledTkin) 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
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
int G4lrint(double ad)
Definition: templates.hh:134

◆ SetVerboseLevel()

void G4PAIPhotModel::SetVerboseLevel ( G4int  verbose)
inline

Definition at line 165 of file G4PAIPhotModel.hh.

166{
167 fVerbose=verbose;
168}

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