Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PAIPhotModel Class Reference

#include <G4PAIPhotModel.hh>

Inheritance diagram for G4PAIPhotModel:
Collaboration diagram for G4PAIPhotModel:

Public Member Functions

 G4PAIPhotModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="PAI")
 
virtual ~G4PAIPhotModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) final
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) final
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) final
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) final
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double) final
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, G4double, G4double) final
 
virtual void DefineForRegion (const G4Region *r) final
 
G4PAIPhotDataGetPAIPhotData ()
 
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 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
- 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
 

Protected Member Functions

virtual 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
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 64 of file G4PAIPhotModel.hh.

Constructor & Destructor Documentation

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

Definition at line 73 of file G4PAIPhotModel.cc.

75  fVerbose(0),
76  fModelData(nullptr),
77  fParticle(nullptr)
78 {
79  fElectron = G4Electron::Electron();
80  fPositron = G4Positron::Positron();
81 
82  fParticleChange = nullptr;
83 
84  if(p) { SetParticle(p); }
85  else { SetParticle(fElectron); }
86 
87  // default generator
89 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4VEmFluctuationModel(const G4String &nam)

Here is the call graph for this function:

G4PAIPhotModel::~G4PAIPhotModel ( )
virtual

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

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4PAIPhotModel.cc.

200 {
201  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
202  if(0 > coupleIndex) { return 0.0; }
203 
204  G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
205 
206  G4double scaledTkin = kineticEnergy*fRatio;
207 
208  return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
209 }
int G4int
Definition: G4Types.hh:78
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4PAIPhotModel::ComputeMaxEnergy ( G4double  scaledEnergy)
inline

Definition at line 161 of file G4PAIPhotModel.hh.

162 {
163  return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
164 }
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 213 of file G4PAIPhotModel.cc.

218 {
219  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
220 
221  if(0 > coupleIndex) return 0.0;
222 
223  G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
224 
225  if(tmax <= cutEnergy) return 0.0;
226 
227  G4double scaledTkin = kineticEnergy*fRatio;
228  G4double xsc = fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
229  scaledTkin,
230  cutEnergy, tmax);
231  return xsc;
232 }
int G4int
Definition: G4Types.hh:78
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) final
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:451
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4PAIPhotModel::DefineForRegion ( const G4Region r)
finalvirtual

Reimplemented from G4VEmModel.

Definition at line 445 of file G4PAIPhotModel.cc.

446 {
447  fPAIRegionVector.push_back(r);
448 }
G4double G4PAIPhotModel::Dispersion ( const G4Material material,
const G4DynamicParticle aParticle,
G4double  tmax,
G4double  step 
)
finalvirtual

Implements G4VEmFluctuationModel.

Definition at line 397 of file G4PAIPhotModel.cc.

401 {
402  G4double particleMass = aParticle->GetMass();
403  G4double electronDensity = material->GetElectronDensity();
404  G4double kineticEnergy = aParticle->GetKineticEnergy();
405  G4double q = aParticle->GetCharge()/eplus;
406  G4double etot = kineticEnergy + particleMass;
407  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
408  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * step
409  * electronDensity * q * q;
410 
411  return siga;
412  /*
413  G4double loss, sumLoss=0., sumLoss2=0., sigma2, meanLoss=0.;
414  for(G4int i = 0; i < fMeanNumber; i++)
415  {
416  loss = SampleFluctuations(material,aParticle,tmax,step,meanLoss);
417  sumLoss += loss;
418  sumLoss2 += loss*loss;
419  }
420  meanLoss = sumLoss/fMeanNumber;
421  sigma2 = meanLoss*meanLoss + (sumLoss2-2*sumLoss*meanLoss)/fMeanNumber;
422  return sigma2;
423  */
424 }
G4double GetKineticEnergy() const
static constexpr double twopi_mc2_rcl2
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double GetMass() const
G4double GetCharge() const
static constexpr double eplus
Definition: G4SIunits.hh:199
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4PAIPhotData * G4PAIPhotModel::GetPAIPhotData ( )
inline

Definition at line 150 of file G4PAIPhotModel.hh.

151 {
152  return fModelData;
153 }

Here is the caller graph for this function:

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

Definition at line 156 of file G4PAIPhotModel.hh.

157 {
158  return fMaterialCutsCoupleVector;
159 }

Here is the caller graph for this function:

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 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
const G4String & GetName() const
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:118
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
const G4String & GetName() const
Definition: G4Material.hh:178
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4String & GetParticleName() const
void Initialise(const G4MaterialCutsCouple *, G4double cut, G4PAIPhotModel *)
static const G4double reg
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4MaterialCutsCouple * FindCouple(G4Material *mat)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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();
191  SetElementSelectors( masterModel->GetElementSelectors() );
192 }
G4PAIPhotData * GetPAIPhotData()
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 428 of file G4PAIPhotModel.cc.

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

Here is the caller graph for this function:

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

Implements G4VEmFluctuationModel.

Definition at line 357 of file G4PAIPhotModel.cc.

361 {
362  // return 0.;
363  G4int coupleIndex = FindCoupleIndex(matCC);
364  if(0 > coupleIndex) { return eloss; }
365 
366  SetParticle(aParticle->GetDefinition());
367 
368 
369  // G4cout << "G4PAIPhotModel::SampleFluctuations step(mm)= "<< step/mm
370  // << " Eloss(keV)= " << eloss/keV << " in "
371  // << matCC->GetMaterial()->GetName() << G4endl;
372 
373 
374  G4double Tkin = aParticle->GetKineticEnergy();
375  G4double scaledTkin = Tkin*fRatio;
376 
377  G4double loss = fModelData->SampleAlongStepPhotonTransfer(coupleIndex, Tkin,
378  scaledTkin,
379  step*fChargeSquare);
380  loss += fModelData->SampleAlongStepPlasmonTransfer(coupleIndex, Tkin,
381  scaledTkin,
382  step*fChargeSquare);
383 
384 
385  // G4cout<<" PAIPhotModel::SampleFluctuations loss = "<<loss/keV<<" keV, on step = "
386  // <<step/mm<<" mm"<<G4endl;
387  return loss;
388 
389 }
G4double GetKineticEnergy() const
G4double SampleAlongStepPhotonTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double stepFactor) const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double SampleAlongStepPlasmonTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double stepFactor) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 239 of file G4PAIPhotModel.cc.

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

Here is the call graph for this function:

void G4PAIPhotModel::SetVerboseLevel ( G4int  verbose)
inline

Definition at line 166 of file G4PAIPhotModel.hh.

167 {
168  fVerbose=verbose;
169 }

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