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

#include <G4CoulombScattering.hh>

Inheritance diagram for G4CoulombScattering:
Collaboration diagram for G4CoulombScattering:

Public Member Functions

 G4CoulombScattering (const G4String &name="CoulombScat")
 
virtual ~G4CoulombScattering ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) final
 
virtual void PrintInfo () override
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
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)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (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=nullptr)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4int GetNumberOfModels () const
 
G4int GetNumberOfRegionModels (size_t couple_index) const
 
G4VEmModelGetRegionModel (G4int idx, size_t couple_index) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4VEmModelGetCurrentModel () 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 SetIntegral (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 *) override
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *) final
 
- Protected Member Functions inherited from G4VEmProcess
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
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)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEmProcess
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
 

Detailed Description

Definition at line 55 of file G4CoulombScattering.hh.

Constructor & Destructor Documentation

G4CoulombScattering::G4CoulombScattering ( const G4String name = "CoulombScat")
explicit

Definition at line 61 of file G4CoulombScattering.cc.

62  : G4VEmProcess(name),q2Max(TeV*TeV),isInitialised(false)
63 {
64  // G4cout << "G4CoulombScattering constructor "<< G4endl;
65  SetBuildTableFlag(true);
66  SetStartFromNullFlag(false);
67  SetIntegral(true);
70 }
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:91
void SetBuildTableFlag(G4bool val)
void SetStartFromNullFlag(G4bool val)
static constexpr double TeV
Definition: G4SIunits.hh:218
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetIntegral(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)

Here is the call graph for this function:

G4CoulombScattering::~G4CoulombScattering ( )
virtual

Definition at line 74 of file G4CoulombScattering.cc.

75 {}

Member Function Documentation

void G4CoulombScattering::InitialiseProcess ( const G4ParticleDefinition p)
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 86 of file G4CoulombScattering.cc.

87 {
88  // second initialisation not allowed for the time being
89  // this means that polar angle limit change will not be appled
90  // after first initialisation
91  if(isInitialised) { return; }
92 
95  q2Max = 0.5*a*a;
96  G4double theta = param->MscThetaLimit();
97 
98  // restricted or non-restricted cross section table
99  G4bool yes = false;
100  if(theta == CLHEP::pi) { yes = true; }
102  /*
103  G4cout << "### G4CoulombScattering::InitialiseProcess: "
104  << p->GetParticleName()
105  << " Emin(MeV)= " << MinKinEnergy()/MeV
106  << " Emax(TeV)= " << MaxKinEnergy()/TeV
107  << " nbins= " << LambdaBinning()
108  << " theta= " << theta
109  << G4endl;
110  */
111 
112  isInitialised = true;
113  G4double mass = p->GetPDGMass();
115  //G4cout << name << " type: " << p->GetParticleType()
116  //<< " mass= " << mass << G4endl;
117  yes = true;
118  if (mass > GeV || p->GetParticleType() == "nucleus") {
119  SetBuildTableFlag(false);
120  yes = false;
121  if(name != "GenericIon") { SetVerboseLevel(0); }
122  } else {
123  if(name != "e-" && name != "e+" &&
124  name != "mu+" && name != "mu-" && name != "pi+" &&
125  name != "kaon+" && name != "proton" ) { SetVerboseLevel(0); }
126  }
127 
128  if(!EmModel(1)) {
129  if(yes) { SetEmModel(new G4eCoulombScatteringModel(), 1); }
130  else { SetEmModel(new G4IonCoulombScatteringModel(), 1); }
131  }
132  G4VEmModel* model = EmModel(1);
133  G4double emin = std::max(param->MinKinEnergy(),model->LowEnergyLimit());
134  G4double emax = std::min(param->MaxKinEnergy(),model->HighEnergyLimit());
135  model->SetPolarAngleLimit(theta);
136  model->SetLowEnergyLimit(emin);
137  model->SetHighEnergyLimit(emax);
138  AddEmModel(1, model);
139 }
const XML_Char * name
Definition: expat.h:151
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double MaxKinEnergy() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
void SetBuildTableFlag(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
static constexpr double hbarc
G4double MscThetaLimit() const
void SetStartFromNullFlag(G4bool val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void SetEmModel(G4VEmModel *, G4int index=1)
bool G4bool
Definition: G4Types.hh:79
const G4String & GetParticleType() const
G4double MinKinEnergy() const
static const G4double emax
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4EmParameters * Instance()
static constexpr double fermi
Definition: SystemOfUnits.h:83
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
const XML_Char XML_Content * model
Definition: expat.h:151
static constexpr double pi
Definition: SystemOfUnits.h:54
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:759
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4double FactorForAngleLimit() const

Here is the call graph for this function:

G4bool G4CoulombScattering::IsApplicable ( const G4ParticleDefinition p)
finalvirtual

Implements G4VEmProcess.

Definition at line 79 of file G4CoulombScattering.cc.

80 {
81  return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
82 }
G4double GetPDGCharge() const

Here is the call graph for this function:

G4double G4CoulombScattering::MinPrimaryEnergy ( const G4ParticleDefinition part,
const G4Material mat 
)
finalprotectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 143 of file G4CoulombScattering.cc.

145 {
146  // Pure Coulomb scattering
147  G4double emin = 0.0;
148 
149  // Coulomb scattering combined with multiple or hadronic scattering
151 
152  if(0.0 < theta) {
153  G4double p2 = q2Max*mat->GetIonisation()->GetInvA23()/(1.0 - cos(theta));
154  G4double mass = part->GetPDGMass();
155  emin = sqrt(p2 + mass*mass) - mass;
156  }
157 
158  return emin;
159 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double MscThetaLimit() const
G4double GetInvA23() const
G4double GetPDGMass() const
static G4EmParameters * Instance()
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4CoulombScattering::PrintInfo ( )
overridevirtual

Implements G4VEmProcess.

Definition at line 163 of file G4CoulombScattering.cc.

164 {
166  << " < Theta(degree) < 180";
167 
168  if(q2Max < DBL_MAX) { G4cout << "; pLimit(GeV^1)= " << sqrt(q2Max)/GeV; }
169  G4cout << G4endl;
170 }
G4double MscThetaLimit() const
G4GLOB_DLL std::ostream G4cout
static constexpr double degree
Definition: G4SIunits.hh:144
static G4EmParameters * Instance()
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:


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