Geant4  10.02.p03
G4PAIModel Class Reference

#include <G4PAIModel.hh>

Inheritance diagram for G4PAIModel:
Collaboration diagram for G4PAIModel:

Public Member Functions

 G4PAIModel (const G4ParticleDefinition *p=0, const G4String &nam="PAI")
 
virtual ~G4PAIModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle *> *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double SampleFluctuations (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double, G4double, G4double)
 
virtual G4double Dispersion (const G4Material *, const G4DynamicParticle *, G4double, G4double)
 
void DefineForRegion (const G4Region *r)
 
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 InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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=0)
 
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

G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Private Member Functions

G4int FindCoupleIndex (const G4MaterialCutsCouple *)
 
void SetParticle (const G4ParticleDefinition *p)
 
G4PAIModeloperator= (const G4PAIModel &right)
 
 G4PAIModel (const G4PAIModel &)
 

Private Attributes

G4int fVerbose
 
G4PAIModelDatafModelData
 
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
 
std::vector< const G4Region * > fPAIRegionVector
 
const G4ParticleDefinitionfParticle
 
const G4ParticleDefinitionfElectron
 
const G4ParticleDefinitionfPositron
 
G4ParticleChangeForLoss * fParticleChange
 
G4double fMass
 
G4double fRatio
 
G4double fChargeSquare
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
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 69 of file G4PAIModel.hh.

Constructor & Destructor Documentation

◆ G4PAIModel() [1/2]

G4PAIModel::G4PAIModel ( const G4ParticleDefinition p = 0,
const G4String nam = "PAI" 
)

Definition at line 77 of file G4PAIModel.cc.

79  fVerbose(0),
80  fModelData(0),
81  fParticle(0)
82 {
85 
86  fParticleChange = 0;
87 
88  if(p) { SetParticle(p); }
89  else { SetParticle(fElectron); }
90 
91  // default generator
93 }
const G4ParticleDefinition * fPositron
Definition: G4PAIModel.hh:145
G4int fVerbose
Definition: G4PAIModel.hh:136
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:143
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:144
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:146
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4VEmFluctuationModel(const G4String &nam)
Here is the call graph for this function:

◆ ~G4PAIModel()

G4PAIModel::~G4PAIModel ( )
virtual

Definition at line 97 of file G4PAIModel.cc.

98 {
99  if(IsMaster()) { delete fModelData; }
100 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
Here is the call graph for this function:

◆ G4PAIModel() [2/2]

G4PAIModel::G4PAIModel ( const G4PAIModel )
private

Member Function Documentation

◆ ComputeDEDXPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 204 of file G4PAIModel.cc.

208 {
209  //G4cout << "===1=== " << CurrentCouple()
210  // << " idx= " << CurrentCouple()->GetIndex()
211  // << " " << fMaterialCutsCoupleVector[0]
212  // << G4endl;
213  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
214  //G4cout << "===2=== " << coupleIndex << G4endl;
215  if(0 > coupleIndex) { return 0.0; }
216 
217  G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
218 
219  G4double scaledTkin = kineticEnergy*fRatio;
220 
221  return fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin,
222  cut);
223 }
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:174
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:369
G4double fRatio
Definition: G4PAIModel.hh:149
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
G4double fChargeSquare
Definition: G4PAIModel.hh:150
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
G4double DEDXPerVolume(G4int coupleIndex, G4double scaledTkin, G4double cut) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ ComputeMaxEnergy()

G4double G4PAIModel::ComputeMaxEnergy ( G4double  scaledEnergy)
inline

Definition at line 164 of file G4PAIModel.hh.

165 {
166  return MaxSecondaryEnergy(fParticle, scaledEnergy/fRatio);
167 }
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:369
G4double fRatio
Definition: G4PAIModel.hh:149
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:143
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 227 of file G4PAIModel.cc.

232 {
233  //G4cout << "===3=== " << CurrentCouple()
234  // << " idx= " << CurrentCouple()->GetIndex()
235  // << " " << fMaterialCutsCoupleVector[0]
236  // << G4endl;
237  G4int coupleIndex = FindCoupleIndex(CurrentCouple());
238  //G4cout << "===4=== " << coupleIndex << G4endl;
239  if(0 > coupleIndex) { return 0.0; }
240 
241  G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
242  if(tmax <= cutEnergy) { return 0.0; }
243 
244  G4double scaledTkin = kineticEnergy*fRatio;
245 
246  return fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
247  scaledTkin,
248  cutEnergy,
249  tmax);
250 }
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:174
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:369
G4double fRatio
Definition: G4PAIModel.hh:149
const G4MaterialCutsCouple * CurrentCouple() const
Definition: G4VEmModel.hh:452
int G4int
Definition: G4Types.hh:78
G4double fChargeSquare
Definition: G4PAIModel.hh:150
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
G4double CrossSectionPerVolume(G4int coupleIndex, G4double scaledTkin, G4double tcut, G4double tmax) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ DefineForRegion()

void G4PAIModel::DefineForRegion ( const G4Region r)
virtual

Reimplemented from G4VEmModel.

Definition at line 386 of file G4PAIModel.cc.

387 {
388  fPAIRegionVector.push_back(r);
389 }
std::vector< const G4Region * > fPAIRegionVector
Definition: G4PAIModel.hh:141

◆ Dispersion()

G4double G4PAIModel::Dispersion ( const G4Material material,
const G4DynamicParticle aParticle,
G4double  tmax,
G4double  step 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 350 of file G4PAIModel.cc.

354 {
355  G4double particleMass = aParticle->GetMass();
356  G4double electronDensity = material->GetElectronDensity();
357  G4double kineticEnergy = aParticle->GetKineticEnergy();
358  G4double q = aParticle->GetCharge()/eplus;
359  G4double etot = kineticEnergy + particleMass;
360  G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
361  G4double siga = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * step
362  * electronDensity * q * q;
363 
364  return siga;
365 }
G4double GetMass() const
int twopi_mc2_rcl2
Definition: hepunit.py:294
G4double GetKineticEnergy() const
G4double GetCharge() const
G4double GetElectronDensity() const
Definition: G4Material.hh:217
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
Here is the call graph for this function:

◆ FindCoupleIndex()

G4int G4PAIModel::FindCoupleIndex ( const G4MaterialCutsCouple couple)
inlineprivate

Definition at line 174 of file G4PAIModel.hh.

175 {
176  G4int idx = -1;
177  size_t jMatMax = fMaterialCutsCoupleVector.size();
178  for(size_t jMat = 0;jMat < jMatMax; ++jMat) {
179  if(couple == fMaterialCutsCoupleVector[jMat]) {
180  idx = jMat;
181  break;
182  }
183  }
184  return idx;
185 }
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:140
int G4int
Definition: G4Types.hh:78
Here is the caller graph for this function:

◆ GetPAIModelData()

G4PAIModelData * G4PAIModel::GetPAIModelData ( )
inline

Definition at line 153 of file G4PAIModel.hh.

154 {
155  return fModelData;
156 }
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
Here is the caller graph for this function:

◆ GetVectorOfCouples()

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

Definition at line 159 of file G4PAIModel.hh.

160 {
162 }
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:140
Here is the caller graph for this function:

◆ Initialise()

void G4PAIModel::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 104 of file G4PAIModel.cc.

106 {
107  if(fVerbose > 0) {
108  G4cout<<"G4PAIModel::Initialise for "<<p->GetParticleName()<<G4endl;
109  }
110  SetParticle(p);
112 
113  if(IsMaster()) {
114 
115  delete fModelData;
116  fMaterialCutsCoupleVector.clear();
117  if(fVerbose > 0) {
118  G4cout << "G4PAIModel instantiates data for " << p->GetParticleName()
119  << G4endl;
120  }
121  G4double tmin = LowEnergyLimit()*fRatio;
123  fModelData = new G4PAIModelData(tmin, tmax, fVerbose);
124 
125  // Prepare initialization
126  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
127  size_t numOfMat = G4Material::GetNumberOfMaterials();
128  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");
135  ->GetRegion("DefaultRegionForTheWorld", false));
136  numRegions = 1;
137  }
138 
139  if(fVerbose > 0) {
140  G4cout << "G4PAIModel is defined for " << numRegions << " regions "
141  << G4endl;
142  G4cout << " total number of materials " << numOfMat << G4endl;
143  }
144  for(size_t iReg = 0; iReg<numRegions; ++iReg) {
145  const G4Region* curReg = fPAIRegionVector[iReg];
146  G4Region* reg = const_cast<G4Region*>(curReg);
147 
148  for(size_t jMat = 0; jMat<numOfMat; ++jMat) {
149  G4Material* mat = (*theMaterialTable)[jMat];
150  const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
151  size_t n = fMaterialCutsCoupleVector.size();
152  /*
153  G4cout << "Region: " << reg->GetName() << " " << reg
154  << " Couple " << cutCouple
155  << " PAI defined for " << n << " couples"
156  << " jMat= " << jMat << " " << mat->GetName()
157  << G4endl;
158  */
159  if(cutCouple) {
160  if(fVerbose > 0) {
161  G4cout << "Region <" << curReg->GetName() << "> mat <"
162  << mat->GetName() << "> CoupleIndex= "
163  << cutCouple->GetIndex()
164  << " " << p->GetParticleName()
165  << " cutsize= " << cuts.size() << G4endl;
166  }
167  // check if this couple is not already initialized
168  G4bool isnew = true;
169  if(0 < n) {
170  for(size_t i=0; i<n; ++i) {
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 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:140
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:120
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
std::vector< const G4Region * > fPAIRegionVector
Definition: G4PAIModel.hh:141
G4double fRatio
Definition: G4PAIModel.hh:149
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:589
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
std::vector< G4Material * > G4MaterialTable
G4int fVerbose
Definition: G4PAIModel.hh:136
static const G4double reg
Float_t mat
static G4RegionStore * GetInstance()
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
Char_t n[5]
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
bool G4bool
Definition: G4Types.hh:79
G4MaterialCutsCouple * FindCouple(G4Material *mat)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:596
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:146
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void Initialise(const G4MaterialCutsCouple *, G4PAIModel *)
const G4String & GetName() const
Definition: G4Material.hh:178
Here is the call graph for this function:

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 192 of file G4PAIModel.cc.

194 {
195  SetParticle(p);
196  fModelData = static_cast<G4PAIModel*>(masterModel)->GetPAIModelData();
198  static_cast<G4PAIModel*>(masterModel)->GetVectorOfCouples();
200 }
const std::vector< const G4MaterialCutsCouple * > & GetVectorOfCouples()
Definition: G4PAIModel.hh:159
G4PAIModelData * GetPAIModelData()
Definition: G4PAIModel.hh:153
std::vector< const G4MaterialCutsCouple * > fMaterialCutsCoupleVector
Definition: G4PAIModel.hh:140
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ MaxSecondaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 369 of file G4PAIModel.cc.

371 {
372  SetParticle(p);
373  G4double tmax = kinEnergy;
374  if(p == fElectron) { tmax *= 0.5; }
375  else if(p != fPositron) {
377  G4double gamma= kinEnergy/fMass + 1.0;
378  tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
379  (1. + 2.0*gamma*ratio + ratio*ratio);
380  }
381  return tmax;
382 }
const G4ParticleDefinition * fPositron
Definition: G4PAIModel.hh:145
G4double fMass
Definition: G4PAIModel.hh:148
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
float electron_mass_c2
Definition: hepunit.py:274
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:144
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

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

◆ SampleFluctuations()

G4double G4PAIModel::SampleFluctuations ( const G4MaterialCutsCouple matCC,
const G4DynamicParticle aParticle,
G4double  tmax,
G4double  step,
G4double  eloss 
)
virtual

Implements G4VEmFluctuationModel.

Definition at line 315 of file G4PAIModel.cc.

319 {
320  G4int coupleIndex = FindCoupleIndex(matCC);
321  if(0 > coupleIndex) { return eloss; }
322 
323  SetParticle(aParticle->GetDefinition());
324 
325  /*
326  G4cout << "G4PAIModel::SampleFluctuations step(mm)= "<< step/mm
327  << " Eloss(keV)= " << eloss/keV << " in "
328  << matCC->Getmaterial()->GetName() << G4endl;
329  */
330 
331  G4double Tkin = aParticle->GetKineticEnergy();
332  G4double scaledTkin = Tkin*fRatio;
333 
334  G4double loss = fModelData->SampleAlongStepTransfer(coupleIndex, Tkin,
335  scaledTkin, tmax,
336  step*fChargeSquare);
337 
338  // G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
339  //<<step/mm<<" mm"<<G4endl;
340  return loss;
341 
342 }
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:174
G4double fRatio
Definition: G4PAIModel.hh:149
int G4int
Definition: G4Types.hh:78
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
G4double GetKineticEnergy() const
G4double fChargeSquare
Definition: G4PAIModel.hh:150
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
G4double SampleAlongStepTransfer(G4int coupleIndex, G4double kinEnergy, G4double scaledTkin, G4double tmax, G4double stepFactor) const
G4ParticleDefinition * GetDefinition() const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 257 of file G4PAIModel.cc.

262 {
263  G4int coupleIndex = FindCoupleIndex(matCC);
264 
265  //G4cout << "G4PAIModel::SampleSecondaries: coupleIndex= "<<coupleIndex<<G4endl;
266  if(0 > coupleIndex) { return; }
267 
268  SetParticle(dp->GetDefinition());
269  G4double kineticEnergy = dp->GetKineticEnergy();
270 
271  G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
272  if(maxEnergy < tmax) { tmax = maxEnergy; }
273  if(tmin >= tmax) { return; }
274 
275  G4ThreeVector direction= dp->GetMomentumDirection();
276  G4double scaledTkin = kineticEnergy*fRatio;
277  G4double totalEnergy = kineticEnergy + fMass;
278  G4double totalMomentum = sqrt(kineticEnergy*(totalEnergy+fMass));
279 
280  G4double deltaTkin =
281  fModelData->SamplePostStepTransfer(coupleIndex, scaledTkin, tmin, tmax);
282 
283  //G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
284  // <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
285 
286  if( !(deltaTkin <= 0.) && !(deltaTkin > 0)) {
287  G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV
288  <<" keV "<< " Escaled(MeV)= " << scaledTkin << G4endl;
289  return;
290  }
291  if( deltaTkin <= 0.) { return; }
292 
293  if( deltaTkin > tmax) { deltaTkin = tmax; }
294 
295  const G4Element* anElement = SelectRandomAtom(matCC,fParticle,kineticEnergy);
296  G4int Z = G4lrint(anElement->GetZ());
297 
299  GetAngularDistribution()->SampleDirection(dp, deltaTkin,
300  Z, matCC->GetMaterial()),
301  deltaTkin);
302 
303  // primary change
304  kineticEnergy -= deltaTkin;
305  G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
306  direction = dir.unit();
307  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
308  fParticleChange->SetProposedMomentumDirection(direction);
309 
310  vdp->push_back(deltaRay);
311 }
G4int FindCoupleIndex(const G4MaterialCutsCouple *)
Definition: G4PAIModel.hh:174
G4double fMass
Definition: G4PAIModel.hh:148
const G4Material * GetMaterial() const
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
Definition: G4PAIModel.cc:369
G4double fRatio
Definition: G4PAIModel.hh:149
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
TDirectory * dir
int G4int
Definition: G4Types.hh:78
void SetParticle(const G4ParticleDefinition *p)
Definition: G4PAIModel.hh:187
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:143
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4PAIModelData * fModelData
Definition: G4PAIModel.hh:138
Hep3Vector unit() const
const G4ParticleDefinition * fElectron
Definition: G4PAIModel.hh:144
G4ParticleChangeForLoss * fParticleChange
Definition: G4PAIModel.hh:146
int G4lrint(double ad)
Definition: templates.hh:163
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
G4double SamplePostStepTransfer(G4int coupleIndex, G4double scaledTkin, G4double tmin, G4double tmax) const
Here is the call graph for this function:

◆ SetParticle()

void G4PAIModel::SetParticle ( const G4ParticleDefinition p)
inlineprivate

Definition at line 187 of file G4PAIModel.hh.

188 {
189  if(fParticle != p) {
190  fParticle = p;
194  fChargeSquare = q*q;
195  }
196 }
G4double fMass
Definition: G4PAIModel.hh:148
G4double fRatio
Definition: G4PAIModel.hh:149
const G4ParticleDefinition * fParticle
Definition: G4PAIModel.hh:143
G4double fChargeSquare
Definition: G4PAIModel.hh:150
static const double proton_mass_c2
static const double eplus
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerboseLevel()

void G4PAIModel::SetVerboseLevel ( G4int  verbose)
inline

Definition at line 169 of file G4PAIModel.hh.

170 {
171  fVerbose=verbose;
172 }
G4int fVerbose
Definition: G4PAIModel.hh:136

Member Data Documentation

◆ fChargeSquare

G4double G4PAIModel::fChargeSquare
private

Definition at line 150 of file G4PAIModel.hh.

◆ fElectron

const G4ParticleDefinition* G4PAIModel::fElectron
private

Definition at line 144 of file G4PAIModel.hh.

◆ fMass

G4double G4PAIModel::fMass
private

Definition at line 148 of file G4PAIModel.hh.

◆ fMaterialCutsCoupleVector

std::vector<const G4MaterialCutsCouple*> G4PAIModel::fMaterialCutsCoupleVector
private

Definition at line 140 of file G4PAIModel.hh.

◆ fModelData

G4PAIModelData* G4PAIModel::fModelData
private

Definition at line 138 of file G4PAIModel.hh.

◆ fPAIRegionVector

std::vector<const G4Region*> G4PAIModel::fPAIRegionVector
private

Definition at line 141 of file G4PAIModel.hh.

◆ fParticle

const G4ParticleDefinition* G4PAIModel::fParticle
private

Definition at line 143 of file G4PAIModel.hh.

◆ fParticleChange

G4ParticleChangeForLoss* G4PAIModel::fParticleChange
private

Definition at line 146 of file G4PAIModel.hh.

◆ fPositron

const G4ParticleDefinition* G4PAIModel::fPositron
private

Definition at line 145 of file G4PAIModel.hh.

◆ fRatio

G4double G4PAIModel::fRatio
private

Definition at line 149 of file G4PAIModel.hh.

◆ fVerbose

G4int G4PAIModel::fVerbose
private

Definition at line 136 of file G4PAIModel.hh.


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