Geant4 11.2.2
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 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 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 67 of file G4PAIPhotModel.cc.

69 fVerbose(0),
70 fModelData(nullptr),
71 fParticle(nullptr)
72{
73 fElectron = G4Electron::Electron();
74 fPositron = G4Positron::Positron();
75
76 fParticleChange = nullptr;
77
78 if(p) { SetParticle(p); }
79 else { SetParticle(fElectron); }
80
81 // default generator
83 fLowestTcut = 12.5*CLHEP::eV;
84}
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 *)

◆ ~G4PAIPhotModel()

G4PAIPhotModel::~G4PAIPhotModel ( )
final

Definition at line 88 of file G4PAIPhotModel.cc.

89{
90 if(IsMaster()) { delete fModelData; fModelData = nullptr; }
91}
G4bool IsMaster() const

◆ 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 195 of file G4PAIPhotModel.cc.

199{
200 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
201 if(0 > coupleIndex) { return 0.0; }
202
203 G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
204 G4double scaledTkin = kineticEnergy*fRatio;
205 G4double dedx = fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
206 return dedx;
207}
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 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 211 of file G4PAIPhotModel.cc.

216{
217 G4int coupleIndex = FindCoupleIndex(CurrentCouple());
218 if(0 > coupleIndex) { return 0.0; }
219
220 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
221 if(tmax <= cutEnergy) { return 0.0; }
222
223 G4double scaledTkin = kineticEnergy*fRatio;
224 G4double xs = fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
225 scaledTkin, cutEnergy, tmax);
226 return xs;
227}
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 426 of file G4PAIPhotModel.cc.

427{
428 fPAIRegionVector.push_back(r);
429}

◆ 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 389 of file G4PAIPhotModel.cc.

394{
395 G4double particleMass = aParticle->GetMass();
396 G4double electronDensity = material->GetElectronDensity();
397 G4double kineticEnergy = aParticle->GetKineticEnergy();
398 G4double q = aParticle->GetCharge()/eplus;
399 G4double etot = kineticEnergy + particleMass;
400 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
401 G4double siga = (tmax/beta2 - 0.5*tcut) * twopi_mc2_rcl2 * step
402 * electronDensity * q * q;
403
404 return siga;
405}
G4double GetMass() const
G4double GetCharge() const
G4double GetKineticEnergy() const
G4double GetElectronDensity() const

◆ 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 95 of file G4PAIPhotModel.cc.

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

Reimplemented from G4VEmModel.

Definition at line 177 of file G4PAIPhotModel.cc.

179{
180 fModelData = static_cast<G4PAIPhotModel*>(masterModel)->GetPAIPhotData();
181 fMaterialCutsCoupleVector = static_cast<G4PAIPhotModel*>(masterModel)->GetVectorOfCouples();
183}
G4PAIPhotData * GetPAIPhotData()
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 409 of file G4PAIPhotModel.cc.

411{
412 SetParticle(p);
413 G4double tmax = kinEnergy;
414 if(p == fElectron) { tmax *= 0.5; }
415 else if(p != fPositron) {
416 G4double ratio= electron_mass_c2/fMass;
417 G4double gamma= kinEnergy/fMass + 1.0;
418 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
419 (1. + 2.0*gamma*ratio + ratio*ratio);
420 }
421 return tmax;
422}

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

◆ MinEnergyCut()

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

Reimplemented from G4VEmModel.

Definition at line 187 of file G4PAIPhotModel.cc.

189{
190 return fLowestTcut;
191}

◆ 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 352 of file G4PAIPhotModel.cc.

357{
358 // return 0.;
359 G4int coupleIndex = FindCoupleIndex(matCC);
360 if(0 > coupleIndex) { return eloss; }
361
362 SetParticle(aParticle->GetDefinition());
363
364 // G4cout << "G4PAIPhotModel::SampleFluctuations step(mm)= "<< step/mm
365 // << " Eloss(keV)= " << eloss/keV << " in "
366 // << matCC->GetMaterial()->GetName() << G4endl;
367
368 G4double Tkin = aParticle->GetKineticEnergy();
369 G4double scaledTkin = Tkin*fRatio;
370
371 G4double loss = fModelData->SampleAlongStepPhotonTransfer(coupleIndex, Tkin,
372 scaledTkin,
373 step*fChargeSquare);
374 loss += fModelData->SampleAlongStepPlasmonTransfer(coupleIndex, Tkin,
375 scaledTkin, step*fChargeSquare);
376
377 // G4cout<<" PAIPhotModel::SampleFluctuations loss = "<<loss/keV<<" keV, on step = "
378 // <<step/mm<<" mm"<<G4endl;
379 return loss;
380
381}
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 234 of file G4PAIPhotModel.cc.

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

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