Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4VEmProcess Class Referenceabstract

#include <G4VEmProcess.hh>

Inheritance diagram for G4VEmProcess:
G4VDiscreteProcess G4VProcess G4ComptonScattering G4CoulombScattering G4DNAAttachment G4DNAChargeDecrease G4DNAChargeIncrease G4DNAElastic G4DNAElectronSolvatation G4DNAExcitation G4DNAIonisation G4DNAVibExcitation G4eeToHadrons G4eplusAnnihilation G4eplusPolarizedAnnihilation G4GammaConversion G4MicroElecElastic G4MicroElecInelastic G4MuElecElastic G4MuElecInelastic G4NuclearStopping G4PhotoElectricEffect G4PolarizedCompton G4PolarizedGammaConversion G4PolarizedPhotoElectricEffect G4RayleighScattering

Public Member Functions

 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void StartTracking (G4Track *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double &kinEnergy, const G4MaterialCutsCouple *couple)
 
void SetLambdaBinning (G4int nbins)
 
G4int LambdaBinning () const
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMinKinEnergyPrim (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4ElementGetCurrentElement () const
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
G4double CrossSectionBiasingFactor () const
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void SetPolarAngleLimit (G4double a)
 
G4double PolarAngleLimit () const
 
void SetLambdaFactor (G4double val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetApplyCuts (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForGamma fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 92 of file G4VEmProcess.hh.

Constructor & Destructor Documentation

G4VEmProcess::G4VEmProcess ( const G4String name,
G4ProcessType  type = fElectromagnetic 
)

Definition at line 90 of file G4VEmProcess.cc.

References DBL_MAX, G4Electron::Electron(), fParticleChange, G4Gamma::Gamma(), G4VProcess::GetProcessName(), G4LossTableManager::Instance(), python.hepunit::keV, G4Positron::Positron(), G4VProcess::pParticleChange, G4LossTableManager::Register(), G4VParticleChange::SetSecondaryWeightByProcess(), G4VProcess::SetVerboseLevel(), and python.hepunit::TeV.

90  :
91  G4VDiscreteProcess(name, type),
92  secondaryParticle(0),
93  buildLambdaTable(true),
94  numberOfModels(0),
95  theLambdaTable(0),
96  theLambdaTablePrim(0),
97  theDensityFactor(0),
98  theDensityIdx(0),
99  integral(false),
100  applyCuts(false),
101  startFromNull(false),
102  splineFlag(true),
103  currentModel(0),
104  particle(0),
105  currentParticle(0),
106  currentCouple(0)
107 {
108  SetVerboseLevel(1);
109 
110  // Size of tables assuming spline
111  minKinEnergy = 0.1*keV;
112  maxKinEnergy = 10.0*TeV;
113  nLambdaBins = 77;
114  minKinEnergyPrim = DBL_MAX;
115 
116  // default lambda factor
117  lambdaFactor = 0.8;
118 
119  // default limit on polar angle
120  polarAngleLimit = 0.0;
121  biasFactor = 1.0;
122 
123  // particle types
124  theGamma = G4Gamma::Gamma();
125  theElectron = G4Electron::Electron();
126  thePositron = G4Positron::Positron();
127 
130  secParticles.reserve(5);
131 
132  preStepLambda = 0.0;
133  mfpKinEnergy = DBL_MAX;
134 
135  idxLambda = idxLambdaPrim = 0;
136 
137  modelManager = new G4EmModelManager();
138  biasManager = 0;
139  biasFlag = false;
140  weightFlag = false;
141  lManager = G4LossTableManager::Instance();
142  lManager->Register(this);
143  secID = fluoID = augerID = biasID = -1;
144  mainSecondaries = 100;
145  if("phot" == GetProcessName() || "compt" == GetProcessName()) {
146  mainSecondaries = 1;
147  }
148 }
static G4LossTableManager * Instance()
G4ParticleChangeForGamma fParticleChange
void SetSecondaryWeightByProcess(G4bool)
void Register(G4VEnergyLossProcess *p)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define DBL_MAX
Definition: templates.hh:83
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4VEmProcess::~G4VEmProcess ( )
virtual

Definition at line 152 of file G4VEmProcess.cc.

References G4LossTableManager::DeRegister(), G4cout, G4endl, G4VProcess::GetProcessName(), G4LossTableManager::IsMaster(), and G4VProcess::verboseLevel.

153 {
154  if(1 < verboseLevel) {
155  G4cout << "G4VEmProcess destruct " << GetProcessName()
156  << " " << this << " " << theLambdaTable <<G4endl;
157  }
158  if(lManager->IsMaster()) {
159  delete theLambdaTable;
160  delete theLambdaTablePrim;
161  }
162  delete modelManager;
163  delete biasManager;
164  lManager->DeRegister(this);
165  //G4cout << "G4VEmProcess removed " << G4endl;
166 }
G4int verboseLevel
Definition: G4VProcess.hh:368
void DeRegister(G4VEnergyLossProcess *p)
G4GLOB_DLL std::ostream G4cout
G4bool IsMaster() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

void G4VEmProcess::ActivateForcedInteraction ( G4double  length = 0.0,
const G4String r = "",
G4bool  flag = true 
)

Definition at line 1081 of file G4VEmProcess.cc.

References G4EmBiasingManager::ActivateForcedInteraction(), G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4VProcess::GetProcessName(), python.hepunit::mm, and G4VProcess::verboseLevel.

Referenced by G4EmProcessOptions::ActivateForcedInteraction().

1083 {
1084  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1085  if(1 < verboseLevel) {
1086  G4cout << "### ActivateForcedInteraction: for "
1087  << particle->GetParticleName()
1088  << " and process " << GetProcessName()
1089  << " length(mm)= " << length/mm
1090  << " in G4Region <" << r
1091  << "> weightFlag= " << flag
1092  << G4endl;
1093  }
1094  weightFlag = flag;
1095  biasManager->ActivateForcedInteraction(length, r);
1096 }
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
void G4VEmProcess::ActivateSecondaryBiasing ( const G4String region,
G4double  factor,
G4double  energyLimit 
)

Definition at line 1101 of file G4VEmProcess.cc.

References G4EmBiasingManager::ActivateSecondaryBiasing(), G4Electron::Electron(), G4cout, G4endl, G4VProcess::GetProcessName(), python.hepunit::MeV, and G4VProcess::verboseLevel.

Referenced by G4EmProcessOptions::ActivateSecondaryBiasingForGamma().

1104 {
1105  if (0.0 <= factor) {
1106 
1107  // Range cut can be applied only for e-
1108  if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1109  { return; }
1110 
1111  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
1112  biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1113  if(1 < verboseLevel) {
1114  G4cout << "### ActivateSecondaryBiasing: for "
1115  << " process " << GetProcessName()
1116  << " factor= " << factor
1117  << " in G4Region <" << region
1118  << "> energyLimit(MeV)= " << energyLimit/MeV
1119  << G4endl;
1120  }
1121  }
1122 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void G4VEmProcess::AddEmModel ( G4int  order,
G4VEmModel p,
const G4Region region = 0 
)

Definition at line 188 of file G4VEmProcess.cc.

References G4EmModelManager::AddEmModel(), fm, G4VProcess::pParticleChange, and G4VEmModel::SetParticleChange().

Referenced by DicomPhysicsList::ConstructEM(), G4EmLowEPPhysics::ConstructProcess(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmDNAPhysics::ConstructProcess(), PhysListEmStandardSSM::ConstructProcess(), PhysListEmLivermore::ConstructProcess(), PhysListEmPenelope::ConstructProcess(), PhysListEmStandardSS::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4DNAVibExcitation::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAElectronSolvatation::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4MicroElecElastic::InitialiseProcess(), G4MuElecElastic::InitialiseProcess(), G4RayleighScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4MuElecInelastic::InitialiseProcess(), G4MicroElecInelastic::InitialiseProcess(), G4eeToHadrons::InitialiseProcess(), G4PolarizedGammaConversion::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4PolarizedCompton::InitialiseProcess(), G4NuclearStopping::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), G4PolarizedPhotoElectricEffect::InitialiseProcess(), G4eplusPolarizedAnnihilation::InitialiseProcess(), and G4EmConfigurator::PrepareModels().

190 {
192  modelManager->AddEmModel(order, p, fm, region);
193  if(p) { p->SetParticleChange(pParticleChange); }
194 }
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
#define fm
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void G4VEmProcess::BuildPhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 325 of file G4VEmProcess.cc.

References G4cout, G4endl, G4LossTableBuilder::GetCoupleIndexes(), G4LossTableBuilder::GetDensityFactors(), G4VProcess::GetMasterProcess(), GetModelByIndex(), G4ParticleDefinition::GetParticleName(), G4VProcess::GetProcessName(), G4LossTableManager::GetTableBuilder(), G4LossTableBuilder::InitialiseBaseMaterials(), G4VEmModel::InitialiseLocal(), LambdaTable(), LambdaTablePrim(), G4EmModelManager::NumberOfModels(), and G4VProcess::verboseLevel.

Referenced by G4PolarizedCompton::BuildPhysicsTable(), and G4eplusPolarizedAnnihilation::BuildPhysicsTable().

326 {
327  const G4VEmProcess* masterProc = 0;
328  if(GetMasterProcess() != this) {
329  masterProc = static_cast<const G4VEmProcess*>(GetMasterProcess());
330  }
331 
332  G4String num = part.GetParticleName();
333  if(1 < verboseLevel) {
334  G4cout << "G4VEmProcess::BuildPhysicsTable() for "
335  << GetProcessName()
336  << " and particle " << num
337  << " buildLambdaTable= " << buildLambdaTable
338  << G4endl;
339  }
340 
341  if(particle == &part) {
342 
343  G4LossTableBuilder* bld = lManager->GetTableBuilder();
344 
345  // worker initialisation
346  if(masterProc) {
347  theLambdaTable = masterProc->LambdaTable();
348  theLambdaTablePrim = masterProc->LambdaTablePrim();
349 
350  if(theLambdaTable) {
351  bld->InitialiseBaseMaterials(theLambdaTable);
352  } else if(theLambdaTablePrim) {
353  bld->InitialiseBaseMaterials(theLambdaTablePrim);
354  }
355  theDensityFactor = bld->GetDensityFactors();
356  theDensityIdx = bld->GetCoupleIndexes();
357  if(theLambdaTable) { FindLambdaMax(); }
358 
359  // local initialisation of models
360  G4bool printing = true;
361  numberOfModels = modelManager->NumberOfModels();
362  for(G4int i=0; i<numberOfModels; ++i) {
363  G4VEmModel* mod = GetModelByIndex(i, printing);
364  G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing);
365  mod->InitialiseLocal(particle, mod0);
366  }
367  // master thread
368  } else {
369  theDensityFactor = bld->GetDensityFactors();
370  theDensityIdx = bld->GetCoupleIndexes();
371  if(buildLambdaTable || minKinEnergyPrim < maxKinEnergy) {
372  BuildLambdaTable();
373  }
374  }
375  }
376 
377  // explicitly defined printout by particle name
378  if(1 < verboseLevel ||
379  (0 < verboseLevel && (num == "gamma" || num == "e-" ||
380  num == "e+" || num == "mu+" ||
381  num == "mu-" || num == "proton"||
382  num == "pi+" || num == "pi-" ||
383  num == "kaon+" || num == "kaon-" ||
384  num == "alpha" || num == "anti_proton" ||
385  num == "GenericIon")))
386  {
387  PrintInfoProcess(part);
388  }
389 
390  if(1 < verboseLevel) {
391  G4cout << "G4VEmProcess::BuildPhysicsTable() done for "
392  << GetProcessName()
393  << " and particle " << num
394  << G4endl;
395  }
396 }
G4PhysicsTable * LambdaTable() const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
G4PhysicsTable * LambdaTablePrim() const
const std::vector< G4double > * GetDensityFactors()
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:206
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4LossTableBuilder * GetTableBuilder()
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4int NumberOfModels() const
void InitialiseBaseMaterials(G4PhysicsTable *table)
#define G4endl
Definition: G4ios.hh:61
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
const std::vector< G4int > * GetCoupleIndexes()
G4double G4VEmProcess::ComputeCrossSectionPerAtom ( G4double  kineticEnergy,
G4double  Z,
G4double  A = 0.,
G4double  cut = 0.0 
)

Definition at line 977 of file G4VEmProcess.cc.

References G4VEmModel::ComputeCrossSectionPerAtom(), SelectModel(), and test::x.

979 {
980  SelectModel(kineticEnergy, currentCoupleIndex);
981  G4double x = 0.0;
982  if(currentModel) {
983  x = currentModel->ComputeCrossSectionPerAtom(currentParticle,kineticEnergy,
984  Z,A,cut);
985  }
986  return x;
987 }
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:300
double G4double
Definition: G4Types.hh:76
G4double G4VEmProcess::CrossSectionBiasingFactor ( ) const
inline

Definition at line 595 of file G4VEmProcess.hh.

596 {
597  return biasFactor;
598 }
G4double G4VEmProcess::CrossSectionPerVolume ( G4double  kineticEnergy,
const G4MaterialCutsCouple couple 
)

Definition at line 934 of file G4VEmProcess.cc.

References G4VEmModel::CrossSectionPerVolume(), and SelectModel().

936 {
937  // Cross section per atom is calculated
938  DefineMaterial(couple);
939  G4double cross = 0.0;
940  if(buildLambdaTable && theLambdaTable) {
941  cross = GetCurrentLambda(kineticEnergy);
942  } else {
943  SelectModel(kineticEnergy, currentCoupleIndex);
944  cross = fFactor*currentModel->CrossSectionPerVolume(currentMaterial,
945  currentParticle,
946  kineticEnergy);
947  }
948 
949  if(cross < 0.0) { cross = 0.0; }
950  return cross;
951 }
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
double G4double
Definition: G4Types.hh:76
size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex ( ) const
inlineprotected
G4VEmModel * G4VEmProcess::EmModel ( G4int  index = 1) const
const G4Element * G4VEmProcess::GetCurrentElement ( ) const

Definition at line 1054 of file G4VEmProcess.cc.

References G4VEmModel::GetCurrentElement().

1055 {
1056  const G4Element* elm = 0;
1057  if(currentModel) {elm = currentModel->GetCurrentElement(); }
1058  return elm;
1059 }
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440
G4double G4VEmProcess::GetElectronEnergyCut ( )
inlineprotected

Definition at line 422 of file G4VEmProcess.hh.

423 {
424  return (*theCutsElectron)[currentCoupleIndex];
425 }
G4double G4VEmProcess::GetGammaEnergyCut ( )
inlineprotected

Definition at line 415 of file G4VEmProcess.hh.

416 {
417  return (*theCutsGamma)[currentCoupleIndex];
418 }
G4double G4VEmProcess::GetLambda ( G4double kinEnergy,
const G4MaterialCutsCouple couple 
)
inline

Definition at line 501 of file G4VEmProcess.hh.

References SelectModel().

Referenced by PostStepDoIt(), and G4AdjointComptonModel::RapidSampleSecondaries().

503 {
504  DefineMaterial(couple);
505  SelectModel(kinEnergy, currentCoupleIndex);
506  return GetCurrentLambda(kinEnergy);
507 }
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4double G4VEmProcess::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 955 of file G4VEmProcess.cc.

References MeanFreePath(), and NotForced.

Referenced by G4eplusPolarizedAnnihilation::GetMeanFreePath(), and G4PolarizedCompton::GetMeanFreePath().

958 {
959  *condition = NotForced;
960  return G4VEmProcess::MeanFreePath(track);
961 }
G4double condition(const G4ErrorSymMatrix &m)
G4double MeanFreePath(const G4Track &track)
G4VEmModel * G4VEmProcess::GetModelByIndex ( G4int  idx = 0,
G4bool  ver = false 
) const

Definition at line 224 of file G4VEmProcess.cc.

References G4EmModelManager::GetModel().

Referenced by BuildPhysicsTable().

225 {
226  return modelManager->GetModel(idx, ver);
227 }
G4VEmModel * GetModel(G4int, G4bool ver=false)
G4ParticleChangeForGamma * G4VEmProcess::GetParticleChange ( )
inlineprotected

Definition at line 665 of file G4VEmProcess.hh.

References fParticleChange.

666 {
667  return &fParticleChange;
668 }
G4ParticleChangeForGamma fParticleChange
virtual void G4VEmProcess::InitialiseProcess ( const G4ParticleDefinition )
protectedpure virtual
virtual G4bool G4VEmProcess::IsApplicable ( const G4ParticleDefinition p)
pure virtual
G4bool G4VEmProcess::IsIntegral ( ) const
inline

Definition at line 644 of file G4VEmProcess.hh.

645 {
646  return integral;
647 }
G4int G4VEmProcess::LambdaBinning ( ) const
inline

Definition at line 551 of file G4VEmProcess.hh.

Referenced by G4PolarizedCompton::BuildAsymmetryTable(), and G4eplusPolarizedAnnihilation::BuildAsymmetryTable().

552 {
553  return nLambdaBins;
554 }
G4PhysicsVector * G4VEmProcess::LambdaPhysicsVector ( const G4MaterialCutsCouple )
protected

Definition at line 1044 of file G4VEmProcess.cc.

References G4PhysicsVector::SetSpline(), G4LossTableManager::SplineFlag(), and test::v.

Referenced by G4PolarizedCompton::BuildAsymmetryTable(), and G4eplusPolarizedAnnihilation::BuildAsymmetryTable().

1045 {
1046  G4PhysicsVector* v =
1047  new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nLambdaBins);
1048  v->SetSpline(lManager->SplineFlag());
1049  return v;
1050 }
G4bool SplineFlag() const
void SetSpline(G4bool)
G4PhysicsTable * G4VEmProcess::LambdaTable ( ) const
inline

Definition at line 602 of file G4VEmProcess.hh.

Referenced by BuildPhysicsTable().

603 {
604  return theLambdaTable;
605 }
G4PhysicsTable * G4VEmProcess::LambdaTablePrim ( ) const
inline

Definition at line 609 of file G4VEmProcess.hh.

Referenced by BuildPhysicsTable().

610 {
611  return theLambdaTablePrim;
612 }
G4double G4VEmProcess::MaxKinEnergy ( ) const
inline
G4double G4VEmProcess::MeanFreePath ( const G4Track track)

Definition at line 965 of file G4VEmProcess.cc.

References DBL_MAX, G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), and test::x.

Referenced by GetMeanFreePath().

966 {
967  DefineMaterial(track.GetMaterialCutsCouple());
968  preStepLambda = GetCurrentLambda(track.GetKineticEnergy());
969  G4double x = DBL_MAX;
970  if(0.0 < preStepLambda) { x = 1.0/preStepLambda; }
971  return x;
972 }
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double G4VEmProcess::MinKinEnergy ( ) const
inline
G4double G4VEmProcess::MinPrimaryEnergy ( const G4ParticleDefinition ,
const G4Material  
)
protectedvirtual

Reimplemented in G4GammaConversion, and G4CoulombScattering.

Definition at line 180 of file G4VEmProcess.cc.

182 {
183  return 0.0;
184 }
const G4ParticleDefinition * G4VEmProcess::Particle ( ) const
inline

Definition at line 616 of file G4VEmProcess.hh.

617 {
618  return particle;
619 }
G4double G4VEmProcess::PolarAngleLimit ( ) const
inline

Definition at line 588 of file G4VEmProcess.hh.

Referenced by G4CoulombScattering::InitialiseProcess(), G4CoulombScattering::MinPrimaryEnergy(), and G4CoulombScattering::PrintInfo().

589 {
590  return polarAngleLimit;
591 }
G4VParticleChange * G4VEmProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 652 of file G4VEmProcess.cc.

References G4VParticleChange::AddSecondary(), G4EmBiasingManager::ApplySecondaryBiasing(), G4VProcess::ClearNumberOfInteractionLengthLeft(), DBL_MAX, python.hepunit::electron_mass_c2, fAlive, G4EmBiasingManager::ForcedInteractionRegion(), fParticleChange, fStopAndKill, fStopButAlive, G4cout, G4endl, G4UniformRand, G4ProcessManager::GetAtRestProcessVector(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetKineticEnergy(), GetLambda(), G4VParticleChange::GetLocalEnergyDeposit(), G4VParticleChange::GetParentWeight(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), G4Step::GetPostStepPoint(), G4ParticleDefinition::GetProcessManager(), G4VProcess::GetProcessName(), G4ParticleChangeForGamma::GetProposedKineticEnergy(), G4StepPoint::GetSafety(), G4Track::GetTouchableHandle(), G4VParticleChange::GetTrackStatus(), G4Track::GetTrackStatus(), G4ParticleChangeForGamma::InitializeForPostStep(), G4VEmModel::IsActive(), python.hepunit::MeV, G4VProcess::pParticleChange, G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::ProposeWeight(), G4VEmModel::SampleSecondaries(), G4EmBiasingManager::SecondaryBiasingRegion(), SelectModel(), G4Track::SetCreatorModelIndex(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), G4Track::SetWeight(), G4ProcessVector::size(), G4VProcess::theNumberOfInteractionLengthLeft, and G4VProcess::verboseLevel.

654 {
655  // In all cases clear number of interaction lengths
657  mfpKinEnergy = DBL_MAX;
658 
660 
661  // Do not make anything if particle is stopped, the annihilation then
662  // should be performed by the AtRestDoIt!
663  if (track.GetTrackStatus() == fStopButAlive) { return &fParticleChange; }
664 
665  G4double finalT = track.GetKineticEnergy();
666 
667  // forced process - should happen only once per track
668  if(biasFlag) {
669  if(biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
670  biasFlag = false;
671  }
672  }
673 
674  // Integral approach
675  if (integral) {
676  G4double lx = GetLambda(finalT, currentCouple);
677  if(preStepLambda<lx && 1 < verboseLevel) {
678  G4cout << "WARNING: for " << currentParticle->GetParticleName()
679  << " and " << GetProcessName()
680  << " E(MeV)= " << finalT/MeV
681  << " preLambda= " << preStepLambda << " < "
682  << lx << " (postLambda) "
683  << G4endl;
684  }
685 
686  if(preStepLambda*G4UniformRand() > lx) {
688  return &fParticleChange;
689  }
690  }
691 
692  SelectModel(finalT, currentCoupleIndex);
693  if(!currentModel->IsActive(finalT)) { return &fParticleChange; }
694 
695  // define new weight for primary and secondaries
697  if(weightFlag) {
698  weight /= biasFactor;
700  }
701 
702  /*
703  if(0 < verboseLevel) {
704  G4cout << "G4VEmProcess::PostStepDoIt: Sample secondary; E= "
705  << finalT/MeV
706  << " MeV; model= (" << currentModel->LowEnergyLimit()
707  << ", " << currentModel->HighEnergyLimit() << ")"
708  << G4endl;
709  }
710  */
711 
712  // sample secondaries
713  secParticles.clear();
714  currentModel->SampleSecondaries(&secParticles,
715  currentCouple,
716  track.GetDynamicParticle(),
717  (*theCuts)[currentCoupleIndex]);
718 
719  G4int num0 = secParticles.size();
720 
721  // splitting or Russian roulette
722  if(biasManager) {
723  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
724  G4double eloss = 0.0;
725  weight *= biasManager->ApplySecondaryBiasing(
726  secParticles, track, currentModel, &fParticleChange, eloss,
727  currentCoupleIndex, (*theCuts)[currentCoupleIndex],
728  step.GetPostStepPoint()->GetSafety());
729  if(eloss > 0.0) {
732  }
733  }
734  }
735 
736  // save secondaries
737  G4int num = secParticles.size();
738  if(num > 0) {
739 
742  G4double time = track.GetGlobalTime();
743 
744  for (G4int i=0; i<num; ++i) {
745  if (secParticles[i]) {
746  G4DynamicParticle* dp = secParticles[i];
748  G4double e = dp->GetKineticEnergy();
749  G4bool good = true;
750  if(applyCuts) {
751  if (p == theGamma) {
752  if (e < (*theCutsGamma)[currentCoupleIndex]) { good = false; }
753 
754  } else if (p == theElectron) {
755  if (e < (*theCutsElectron)[currentCoupleIndex]) { good = false; }
756 
757  } else if (p == thePositron) {
758  if (electron_mass_c2 < (*theCutsGamma)[currentCoupleIndex] &&
759  e < (*theCutsPositron)[currentCoupleIndex]) {
760  good = false;
761  e += 2.0*electron_mass_c2;
762  }
763  }
764  // added secondary if it is good
765  }
766  if (good) {
767  G4Track* t = new G4Track(dp, time, track.GetPosition());
769  t->SetWeight(weight);
771 
772  // define type of secondary
773  if(i < mainSecondaries) { t->SetCreatorModelIndex(secID); }
774  else if(i < num0) {
775  if(p == theGamma) {
776  t->SetCreatorModelIndex(fluoID);
777  } else {
778  t->SetCreatorModelIndex(augerID);
779  }
780  } else {
781  t->SetCreatorModelIndex(biasID);
782  }
783 
784  //G4cout << "Secondary(post step) has weight " << t->GetWeight()
785  // << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
786  } else {
787  delete dp;
788  edep += e;
789  }
790  }
791  }
793  }
794 
797  if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
800  }
801 
802  // ClearNumberOfInteractionLengthLeft();
803  return &fParticleChange;
804 }
G4double GetLambda(G4double &kinEnergy, const G4MaterialCutsCouple *couple)
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4int verboseLevel
Definition: G4VProcess.hh:368
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4bool SecondaryBiasingRegion(G4int coupleIdx)
const char * p
Definition: xmltok.h:285
const G4ThreeVector & GetPosition() const
G4ParticleChangeForGamma fParticleChange
G4TrackStatus GetTrackStatus() const
G4bool ForcedInteractionRegion(G4int coupleIdx)
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetParentWeight() const
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetCreatorModelIndex(G4int idx)
Definition: G4Track.hh:245
void ProposeWeight(G4double finalWeight)
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
float electron_mass_c2
Definition: hepunit.py:274
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetGlobalTime() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void AddSecondary(G4Track *aSecondary)
const G4TouchableHandle & GetTouchableHandle() const
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
G4int size() const
G4double GetProposedKineticEnergy() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4double GetLocalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4double GetSafety() const
#define G4endl
Definition: G4ios.hh:61
G4TrackStatus GetTrackStatus() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void InitializeForPostStep(const G4Track &)
#define DBL_MAX
Definition: templates.hh:83
G4double G4VEmProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 570 of file G4VEmProcess.cc.

References python.hepunit::cm, G4VProcess::currentInteractionLength, DBL_MAX, G4EmBiasingManager::ForcedInteractionRegion(), G4cout, G4endl, G4Log(), G4UniformRand, G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4Material::GetName(), G4Track::GetParentID(), G4ParticleDefinition::GetParticleName(), G4VProcess::GetProcessName(), G4EmBiasingManager::GetStepLimit(), G4VEmModel::IsActive(), python.hepunit::MeV, NotForced, SelectModel(), G4VProcess::theInitialNumberOfInteractionLength, G4VProcess::theNumberOfInteractionLengthLeft, G4VProcess::verboseLevel, and test::x.

Referenced by G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), and G4PolarizedCompton::PostStepGetPhysicalInteractionLength().

574 {
575  *condition = NotForced;
576  G4double x = DBL_MAX;
577 
578  preStepKinEnergy = track.GetKineticEnergy();
579  DefineMaterial(track.GetMaterialCutsCouple());
580  SelectModel(preStepKinEnergy, currentCoupleIndex);
581 
582  if(!currentModel->IsActive(preStepKinEnergy)) {
584  return x;
585  }
586 
587  // forced biasing only for primary particles
588  if(biasManager) {
589  if(0 == track.GetParentID()) {
590  if(biasFlag &&
591  biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
592  return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
593  }
594  }
595  }
596 
597  // compute mean free path
598  if(preStepKinEnergy < mfpKinEnergy) {
599  if (integral) { ComputeIntegralLambda(preStepKinEnergy); }
600  else { preStepLambda = GetCurrentLambda(preStepKinEnergy); }
601 
602  // zero cross section
603  if(preStepLambda <= 0.0) {
606  }
607  }
608 
609  // non-zero cross section
610  if(preStepLambda > 0.0) {
611 
613 
614  // beggining of tracking (or just after DoIt of this process)
615  // ResetNumberOfInteractionLengthLeft();
616 
619 
620  } else if(currentInteractionLength < DBL_MAX) {
621 
622  // subtract NumberOfInteractionLengthLeft using previous step
624  previousStepSize/currentInteractionLength;
625  //SubtractNumberOfInteractionLengthLeft(previousStepSize);
628  }
629  }
630 
631  // new mean free path and step limit for the next step
632  currentInteractionLength = 1.0/preStepLambda;
634 #ifdef G4VERBOSE
635  if (verboseLevel>2){
636  G4cout << "G4VEmProcess::PostStepGetPhysicalInteractionLength ";
637  G4cout << "[ " << GetProcessName() << "]" << G4endl;
638  G4cout << " for " << currentParticle->GetParticleName()
639  << " in Material " << currentMaterial->GetName()
640  << " Ekin(MeV)= " << preStepKinEnergy/MeV
641  <<G4endl;
642  G4cout << " MeanFreePath = " << currentInteractionLength/cm << "[cm]"
643  << " InteractionLength= " << x/cm <<"[cm] " <<G4endl;
644  }
645 #endif
646  }
647  return x;
648 }
G4double condition(const G4ErrorSymMatrix &m)
G4int GetParentID() const
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetName() const
Definition: G4Material.hh:176
G4bool ForcedInteractionRegion(G4int coupleIdx)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
const G4String & GetParticleName() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
G4double G4Log(G4double x)
Definition: G4Log.hh:227
#define G4endl
Definition: G4ios.hh:61
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:300
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void G4VEmProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 231 of file G4VEmProcess.cc.

References G4LossTableManager::AtomDeexcitation(), DBL_MAX, G4cout, G4endl, G4GenericIon::GenericIon(), G4ProductionCutsTable::GetEnergyCutsVector(), G4VProcess::GetMasterProcess(), G4EmModelManager::GetModel(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleSubType(), G4ParticleDefinition::GetParticleType(), G4VProcess::GetProcessName(), G4ProductionCutsTable::GetProductionCutsTable(), G4LossTableManager::GetTableBuilder(), G4ProductionCutsTable::GetTableSize(), G4VEmModel::HighEnergyLimit(), idxG4ElectronCut, idxG4GammaCut, idxG4PositronCut, G4EmBiasingManager::Initialise(), G4EmModelManager::Initialise(), G4LossTableBuilder::InitialiseBaseMaterials(), InitialiseProcess(), n, G4EmModelManager::NumberOfModels(), eplot::pname, G4PhysicsTableHelper::PreparePhysicsTable(), G4LossTableManager::PreparePhysicsTable(), G4PhysicsModelCatalog::Register(), G4EmModelManager::SetFluoFlag(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetMasterThread(), SetParticle(), G4VEmModel::SetPolarAngleLimit(), and G4VProcess::verboseLevel.

Referenced by G4PolarizedCompton::PreparePhysicsTable(), and G4eplusPolarizedAnnihilation::PreparePhysicsTable().

232 {
233  G4bool isMaster = false;
234  if(GetMasterProcess() == this) { isMaster = true; }
235  if(!particle) { SetParticle(&part); }
236 
237  if(part.GetParticleType() == "nucleus" &&
238  part.GetParticleSubType() == "generic") {
239 
240  G4String pname = part.GetParticleName();
241  if(pname != "deuteron" && pname != "triton" &&
242  pname != "alpha" && pname != "He3" &&
243  pname != "alpha+" && pname != "helium" &&
244  pname != "hydrogen") {
245 
246  particle = G4GenericIon::GenericIon();
247  }
248  }
249 
250  if(1 < verboseLevel) {
251  G4cout << "G4VEmProcess::PreparePhysicsTable() for "
252  << GetProcessName()
253  << " and particle " << part.GetParticleName()
254  << " local particle " << particle->GetParticleName()
255  << G4endl;
256  }
257 
258  if(particle != &part) { return; }
259 
260  G4LossTableBuilder* bld = lManager->GetTableBuilder();
261 
262  lManager->PreparePhysicsTable(&part, this, isMaster);
263 
264  Clear();
265  InitialiseProcess(particle);
266 
267  const G4ProductionCutsTable* theCoupleTable=
269  size_t n = theCoupleTable->GetTableSize();
270 
271  theEnergyOfCrossSectionMax.resize(n, 0.0);
272  theCrossSectionMax.resize(n, DBL_MAX);
273 
274  // initialisation of models
275  numberOfModels = modelManager->NumberOfModels();
276  for(G4int i=0; i<numberOfModels; ++i) {
277  G4VEmModel* mod = modelManager->GetModel(i);
278  if(0 == i) { currentModel = mod; }
279  mod->SetPolarAngleLimit(polarAngleLimit);
280  mod->SetMasterThread(isMaster);
281  if(mod->HighEnergyLimit() > maxKinEnergy) {
282  mod->SetHighEnergyLimit(maxKinEnergy);
283  }
284  }
285 
286  if(lManager->AtomDeexcitation()) { modelManager->SetFluoFlag(true); }
287  theCuts = modelManager->Initialise(particle,secondaryParticle,
288  2.,verboseLevel);
289  theCutsGamma = theCoupleTable->GetEnergyCutsVector(idxG4GammaCut);
290  theCutsElectron = theCoupleTable->GetEnergyCutsVector(idxG4ElectronCut);
291  theCutsPositron = theCoupleTable->GetEnergyCutsVector(idxG4PositronCut);
292 
293  // prepare tables
294  if(buildLambdaTable && isMaster){
295  theLambdaTable =
297  bld->InitialiseBaseMaterials(theLambdaTable);
298  }
299  // high energy table
300  if(isMaster && minKinEnergyPrim < maxKinEnergy){
301  theLambdaTablePrim =
302  G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTablePrim);
303  bld->InitialiseBaseMaterials(theLambdaTablePrim);
304  }
305  // forced biasing
306  if(biasManager) {
307  biasManager->Initialise(part,GetProcessName(),verboseLevel);
308  biasFlag = false;
309  }
310  // defined ID of secondary particles
311  if(isMaster) {
312  G4String nam1 = GetProcessName();
313  G4String nam2 = nam1 + "_fluo" ;
314  G4String nam3 = nam1 + "_auger";
315  G4String nam4 = nam1 + "_split";
316  secID = G4PhysicsModelCatalog::Register(nam1);
317  fluoID = G4PhysicsModelCatalog::Register(nam2);
318  augerID = G4PhysicsModelCatalog::Register(nam3);
319  biasID = G4PhysicsModelCatalog::Register(nam4);
320  }
321 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4VEmModel * GetModel(G4int, G4bool ver=false)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
const G4String & GetParticleSubType() const
int G4int
Definition: G4Types.hh:78
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
const G4String & GetParticleName() const
G4LossTableBuilder * GetTableBuilder()
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetParticle(const G4ParticleDefinition *p)
const G4String & GetParticleType() const
const G4int n
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
string pname
Definition: eplot.py:33
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:669
static G4ProductionCutsTable * GetProductionCutsTable()
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
G4int NumberOfModels() const
void InitialiseBaseMaterials(G4PhysicsTable *table)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
#define G4endl
Definition: G4ios.hh:61
void SetFluoFlag(G4bool val)
G4VAtomDeexcitation * AtomDeexcitation()
static G4int Register(G4String &)
#define DBL_MAX
Definition: templates.hh:83
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
virtual void G4VEmProcess::PrintInfo ( )
pure virtual
G4double G4VEmProcess::RecalculateLambda ( G4double  kinEnergy,
const G4MaterialCutsCouple couple 
)
inlineprotected

Definition at line 512 of file G4VEmProcess.hh.

References SelectModel().

513 {
514  DefineMaterial(couple);
515  SelectModel(e, currentCoupleIndex);
516  return fFactor*ComputeCurrentLambda(e);
517 }
G4VEmModel * SelectModel(G4double &kinEnergy, size_t index)
G4bool G4VEmProcess::RetrievePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii 
)
virtual

Reimplemented from G4VProcess.

Definition at line 856 of file G4VEmProcess.cc.

References G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4VProcess::GetPhysicsTableFileName(), G4VProcess::GetProcessName(), G4PhysicsTable::length(), n, G4PhysicsTableHelper::RetrievePhysicsTable(), G4LossTableManager::SplineFlag(), and G4VProcess::verboseLevel.

859 {
860  if(1 < verboseLevel) {
861  G4cout << "G4VEmProcess::RetrievePhysicsTable() for "
862  << part->GetParticleName() << " and process "
863  << GetProcessName() << G4endl;
864  }
865  G4bool yes = true;
866 
867  if((!buildLambdaTable && minKinEnergyPrim > maxKinEnergy)
868  || particle != part) { return yes; }
869 
870  const G4String particleName = part->GetParticleName();
871  G4String filename;
872 
873  if(buildLambdaTable) {
874  filename = GetPhysicsTableFileName(part,directory,"Lambda",ascii);
875  yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTable,
876  filename,ascii);
877  if ( yes ) {
878  if (0 < verboseLevel) {
879  G4cout << "Lambda table for " << particleName
880  << " is Retrieved from <"
881  << filename << ">"
882  << G4endl;
883  }
884  if(lManager->SplineFlag()) {
885  size_t n = theLambdaTable->length();
886  for(size_t i=0; i<n; ++i) {
887  if((* theLambdaTable)[i]) {
888  (* theLambdaTable)[i]->SetSpline(true);
889  }
890  }
891  }
892  } else {
893  if (1 < verboseLevel) {
894  G4cout << "Lambda table for " << particleName << " in file <"
895  << filename << "> is not exist"
896  << G4endl;
897  }
898  }
899  }
900  if(minKinEnergyPrim < maxKinEnergy) {
901  filename = GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
902  yes = G4PhysicsTableHelper::RetrievePhysicsTable(theLambdaTablePrim,
903  filename,ascii);
904  if ( yes ) {
905  if (0 < verboseLevel) {
906  G4cout << "Lambda table prim for " << particleName
907  << " is Retrieved from <"
908  << filename << ">"
909  << G4endl;
910  }
911  if(lManager->SplineFlag()) {
912  size_t n = theLambdaTablePrim->length();
913  for(size_t i=0; i<n; ++i) {
914  if((* theLambdaTablePrim)[i]) {
915  (* theLambdaTablePrim)[i]->SetSpline(true);
916  }
917  }
918  }
919  } else {
920  if (1 < verboseLevel) {
921  G4cout << "Lambda table prim for " << particleName << " in file <"
922  << filename << "> is not exist"
923  << G4endl;
924  }
925  }
926  }
927 
928  return yes;
929 }
G4bool SplineFlag() const
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4int n
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
size_t length() const
#define G4endl
Definition: G4ios.hh:61
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
const G4ParticleDefinition * G4VEmProcess::SecondaryParticle ( ) const
inline

Definition at line 623 of file G4VEmProcess.hh.

624 {
625  return secondaryParticle;
626 }
G4VEmModel * G4VEmProcess::SelectModel ( G4double kinEnergy,
size_t  index 
)
inlineprotected

Definition at line 447 of file G4VEmProcess.hh.

References G4EmModelManager::SelectModel(), and G4VEmModel::SetCurrentCouple().

Referenced by G4NuclearStopping::AlongStepDoIt(), ComputeCrossSectionPerAtom(), CrossSectionPerVolume(), GetLambda(), PostStepDoIt(), PostStepGetPhysicalInteractionLength(), and RecalculateLambda().

448 {
449  if(1 < numberOfModels) {
450  currentModel = modelManager->SelectModel(kinEnergy, index);
451  }
452  currentModel->SetCurrentCouple(currentCouple);
453  return currentModel;
454 }
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:419
G4VEmModel * SelectModel(G4double &energy, size_t &index)
G4VEmModel * G4VEmProcess::SelectModelForMaterial ( G4double  kinEnergy,
size_t &  idxRegion 
) const
inline

Definition at line 459 of file G4VEmProcess.hh.

References G4EmModelManager::SelectModel().

461 {
462  return modelManager->SelectModel(kinEnergy, idxRegion);
463 }
G4VEmModel * SelectModel(G4double &energy, size_t &index)
void G4VEmProcess::SetApplyCuts ( G4bool  val)
inline

Definition at line 651 of file G4VEmProcess.hh.

Referenced by G4EmProcessOptions::SetApplyCuts().

652 {
653  applyCuts = val;
654 }
void G4VEmProcess::SetBuildTableFlag ( G4bool  val)
inline
void G4VEmProcess::SetCrossSectionBiasingFactor ( G4double  f,
G4bool  flag = true 
)

Definition at line 1063 of file G4VEmProcess.cc.

References G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.

Referenced by G4EmProcessOptions::SetProcessBiasingFactor().

1064 {
1065  if(f > 0.0) {
1066  biasFactor = f;
1067  weightFlag = flag;
1068  if(1 < verboseLevel) {
1069  G4cout << "### SetCrossSectionBiasingFactor: for "
1070  << particle->GetParticleName()
1071  << " and process " << GetProcessName()
1072  << " biasFactor= " << f << " weightFlag= " << flag
1073  << G4endl;
1074  }
1075  }
1076 }
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
void G4VEmProcess::SetEmModel ( G4VEmModel p,
G4int  index = 1 
)

Definition at line 198 of file G4VEmProcess.cc.

References n.

Referenced by DMXPhysicsList::ConstructEM(), G4EmLivermorePolarizedPhysics::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmDNAPhysics::ConstructProcess(), PhysListEmStandardNR::ConstructProcess(), PhysListEmStandardWVI::ConstructProcess(), PhysListEmStandard_GS::ConstructProcess(), PhysListEmStandard_WVI::ConstructProcess(), PhysListEmStandard_option3::ConstructProcess(), PhysListEmStandard_SS::ConstructProcess(), PhysListEmStandard_option0::ConstructProcess(), GammaRayTelEMlowePhysics::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4DNAVibExcitation::InitialiseProcess(), G4DNAChargeIncrease::InitialiseProcess(), G4DNAAttachment::InitialiseProcess(), G4DNAElastic::InitialiseProcess(), G4DNAElectronSolvatation::InitialiseProcess(), G4DNAChargeDecrease::InitialiseProcess(), G4DNAIonisation::InitialiseProcess(), G4DNAExcitation::InitialiseProcess(), G4MuElecElastic::InitialiseProcess(), G4MicroElecElastic::InitialiseProcess(), G4RayleighScattering::InitialiseProcess(), G4CoulombScattering::InitialiseProcess(), G4MicroElecInelastic::InitialiseProcess(), G4MuElecInelastic::InitialiseProcess(), G4ComptonScattering::InitialiseProcess(), G4eplusAnnihilation::InitialiseProcess(), G4NuclearStopping::InitialiseProcess(), G4GammaConversion::InitialiseProcess(), G4PhotoElectricEffect::InitialiseProcess(), and G4PolarizedPhotoElectricEffect::InitialiseProcess().

199 {
200  G4int n = emModels.size();
201  if(index >= n) { for(G4int i=n; i<=index; ++i) {emModels.push_back(0);} }
202  emModels[index] = p;
203 }
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
const G4int n
void G4VEmProcess::SetIntegral ( G4bool  val)
inline

Definition at line 637 of file G4VEmProcess.hh.

Referenced by G4CoulombScattering::G4CoulombScattering(), G4eplusAnnihilation::G4eplusAnnihilation(), and G4eeToHadrons::InitialiseProcess().

638 {
639  if(particle && particle != theGamma) { integral = val; }
640 }
void G4VEmProcess::SetLambdaBinning ( G4int  nbins)
inline
void G4VEmProcess::SetLambdaFactor ( G4double  val)
inline

Definition at line 630 of file G4VEmProcess.hh.

631 {
632  if(val > 0.0 && val <= 1.0) { lambdaFactor = val; }
633 }
void G4VEmProcess::SetMaxKinEnergy ( G4double  e)
void G4VEmProcess::SetMinKinEnergy ( G4double  e)
void G4VEmProcess::SetMinKinEnergyPrim ( G4double  e)
inline
void G4VEmProcess::SetParticle ( const G4ParticleDefinition p)
inlineprotected

Definition at line 672 of file G4VEmProcess.hh.

Referenced by G4eeToHadrons::InitialiseProcess(), and PreparePhysicsTable().

673 {
674  particle = p;
675  currentParticle = p;
676 }
const char * p
Definition: xmltok.h:285
void G4VEmProcess::SetPolarAngleLimit ( G4double  a)
inline

Definition at line 579 of file G4VEmProcess.hh.

Referenced by G4EmProcessOptions::SetPolarAngleLimit().

580 {
581  if(val < 0.0) { polarAngleLimit = 0.0; }
582  else if(val > CLHEP::pi) { polarAngleLimit = CLHEP::pi; }
583  else { polarAngleLimit = val; }
584 }
void G4VEmProcess::SetSecondaryParticle ( const G4ParticleDefinition p)
inlineprotected
void G4VEmProcess::SetSplineFlag ( G4bool  val)
inlineprotected
void G4VEmProcess::SetStartFromNullFlag ( G4bool  val)
inlineprotected
void G4VEmProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 549 of file G4VEmProcess.cc.

References DBL_MAX, G4Track::GetParentID(), G4Track::GetParticleDefinition(), G4EmBiasingManager::ResetForcedInteraction(), and G4VProcess::theNumberOfInteractionLengthLeft.

550 {
551  // reset parameters for the new track
552  currentParticle = track->GetParticleDefinition();
554  //currentInteractionLength = -1.0;
555  // theInitialNumberOfInteractionLength=-1.0;
556  mfpKinEnergy = DBL_MAX;
557 
558  // forced biasing only for primary particles
559  if(biasManager) {
560  if(0 == track->GetParentID()) {
561  // primary particle
562  biasFlag = true;
563  biasManager->ResetForcedInteraction();
564  }
565  }
566 }
G4int GetParentID() const
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
const G4ParticleDefinition * GetParticleDefinition() const
#define DBL_MAX
Definition: templates.hh:83
G4bool G4VEmProcess::StorePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii = false 
)
virtual

Reimplemented from G4VProcess.

Definition at line 808 of file G4VEmProcess.cc.

References G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4VProcess::GetPhysicsTableFileName(), G4VProcess::GetProcessName(), and G4PhysicsTable::StorePhysicsTable().

811 {
812  G4bool yes = true;
813 
814  if ( theLambdaTable && part == particle) {
815  const G4String name =
816  GetPhysicsTableFileName(part,directory,"Lambda",ascii);
817  yes = theLambdaTable->StorePhysicsTable(name,ascii);
818 
819  if ( yes ) {
820  G4cout << "Physics table is stored for " << particle->GetParticleName()
821  << " and process " << GetProcessName()
822  << " in the directory <" << directory
823  << "> " << G4endl;
824  } else {
825  G4cout << "Fail to store Physics Table for "
826  << particle->GetParticleName()
827  << " and process " << GetProcessName()
828  << " in the directory <" << directory
829  << "> " << G4endl;
830  }
831  }
832  if ( theLambdaTablePrim && part == particle) {
833  const G4String name =
834  GetPhysicsTableFileName(part,directory,"LambdaPrim",ascii);
835  yes = theLambdaTablePrim->StorePhysicsTable(name,ascii);
836 
837  if ( yes ) {
838  G4cout << "Physics table prim is stored for "
839  << particle->GetParticleName()
840  << " and process " << GetProcessName()
841  << " in the directory <" << directory
842  << "> " << G4endl;
843  } else {
844  G4cout << "Fail to store Physics Table Prim for "
845  << particle->GetParticleName()
846  << " and process " << GetProcessName()
847  << " in the directory <" << directory
848  << "> " << G4endl;
849  }
850  }
851  return yes;
852 }
const XML_Char * name
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
void G4VEmProcess::UpdateEmModel ( const G4String nam,
G4double  emin,
G4double  emax 
)

Definition at line 216 of file G4VEmProcess.cc.

References G4EmModelManager::UpdateEmModel().

218 {
219  modelManager->UpdateEmModel(nam, emin, emax);
220 }
void UpdateEmModel(const G4String &, G4double, G4double)

Field Documentation

G4ParticleChangeForGamma G4VEmProcess::fParticleChange
protected

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