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

#include <G4CoulombScattering.hh>

Inheritance diagram for G4CoulombScattering:
G4VEmProcess G4VDiscreteProcess G4VProcess

Public Member Functions

 G4CoulombScattering (const G4String &name="CoulombScat")
 
virtual ~G4CoulombScattering ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
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 *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Protected Member Functions inherited from G4VEmProcess
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 ()
 

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")

Definition at line 60 of file G4CoulombScattering.cc.

References fCoulombScattering, G4Proton::Proton(), G4VEmProcess::SetBuildTableFlag(), G4VEmProcess::SetIntegral(), G4VProcess::SetProcessSubType(), G4VEmProcess::SetSecondaryParticle(), G4VEmProcess::SetSplineFlag(), and G4VEmProcess::SetStartFromNullFlag().

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

Definition at line 74 of file G4CoulombScattering.cc.

75 {}

Member Function Documentation

void G4CoulombScattering::InitialiseProcess ( const G4ParticleDefinition p)
protectedvirtual

Implements G4VEmProcess.

Definition at line 86 of file G4CoulombScattering.cc.

References test::a, G4VEmProcess::AddEmModel(), G4VEmProcess::EmModel(), G4LossTableManager::FactorForAngleLimit(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGMass(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4LossTableManager::Instance(), G4VEmModel::LowEnergyLimit(), G4INCL::Math::max(), G4VEmProcess::MaxKinEnergy(), G4INCL::Math::min(), G4VEmProcess::MinKinEnergy(), G4VEmProcess::PolarAngleLimit(), G4VEmProcess::SetBuildTableFlag(), G4VEmProcess::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), G4VEmModel::SetPolarAngleLimit(), G4VEmProcess::SetStartFromNullFlag(), and G4VProcess::SetVerboseLevel().

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 
94  *CLHEP::hbarc/CLHEP::fermi;
95  q2Max = 0.5*a*a;
96  G4double theta = PolarAngleLimit();
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  // second initialisation
113  if(isInitialised) {
114  G4VEmModel* mod = EmModel(1);
115  mod->SetPolarAngleLimit(theta);
116  mod = GetModelByIndex(1);
117  if(mod) { mod->SetPolarAngleLimit(theta); }
118 
119  // first initialisation
120  } else {
121  */
122 
123  isInitialised = true;
124  G4double mass = p->GetPDGMass();
126  //G4cout << name << " type: " << p->GetParticleType()
127  //<< " mass= " << mass << G4endl;
128  if (mass > GeV || p->GetParticleType() == "nucleus") {
129  SetBuildTableFlag(false);
130  if(name != "GenericIon") { SetVerboseLevel(0); }
131  } else {
132  if(name != "e-" && name != "e+" &&
133  name != "mu+" && name != "mu-" && name != "pi+" &&
134  name != "kaon+" && name != "proton" ) { SetVerboseLevel(0); }
135  }
136 
137  if(!EmModel(1)) { SetEmModel(new G4eCoulombScatteringModel(), 1); }
138  G4VEmModel* model = EmModel(1);
139  G4double emin = std::max(MinKinEnergy(),model->LowEnergyLimit());
140  G4double emax = std::min(MaxKinEnergy(),model->HighEnergyLimit());
141  model->SetPolarAngleLimit(theta);
142  model->SetLowEnergyLimit(emin);
143  model->SetHighEnergyLimit(emax);
144  AddEmModel(1, model);
145 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void SetBuildTableFlag(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
G4double MaxKinEnergy() const
G4double FactorForAngleLimit() const
const XML_Char * name
void SetStartFromNullFlag(G4bool val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetEmModel(G4VEmModel *, G4int index=1)
G4double PolarAngleLimit() const
bool G4bool
Definition: G4Types.hh:79
const XML_Char XML_Content * model
const G4String & GetParticleType() const
G4double GetPDGMass() const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4double MinKinEnergy() const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4bool G4CoulombScattering::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 79 of file G4CoulombScattering.cc.

References G4ParticleDefinition::GetPDGCharge(), and G4ParticleDefinition::IsShortLived().

80 {
81  return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
82 }
G4double GetPDGCharge() const
G4double G4CoulombScattering::MinPrimaryEnergy ( const G4ParticleDefinition part,
const G4Material mat 
)
protectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 149 of file G4CoulombScattering.cc.

References G4IonisParamMat::GetInvA23(), G4Material::GetIonisation(), G4ParticleDefinition::GetPDGMass(), and G4VEmProcess::PolarAngleLimit().

151 {
152  // Pure Coulomb scattering
153  G4double emin = 0.0;
154 
155  // Coulomb scattering combined with multiple or hadronic scattering
156  G4double theta = PolarAngleLimit();
157  if(0.0 < theta) {
158  G4double p2 = q2Max*mat->GetIonisation()->GetInvA23()/(1.0 - cos(theta));
159  G4double mass = part->GetPDGMass();
160  emin = sqrt(p2 + mass*mass) - mass;
161  }
162 
163  return emin;
164 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double PolarAngleLimit() const
G4double GetInvA23() const
G4double GetPDGMass() const
double G4double
Definition: G4Types.hh:76
void G4CoulombScattering::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 168 of file G4CoulombScattering.cc.

References DBL_MAX, python.hepunit::degree, G4cout, G4endl, python.hepunit::GeV, and G4VEmProcess::PolarAngleLimit().

169 {
170  G4cout << " " << PolarAngleLimit()/degree
171  << " < Theta(degree) < 180";
172 
173  if(q2Max < DBL_MAX) { G4cout << "; pLimit(GeV^1)= " << sqrt(q2Max)/GeV; }
174  G4cout << G4endl;
175 }
G4GLOB_DLL std::ostream G4cout
tuple degree
Definition: hepunit.py:69
G4double PolarAngleLimit() const
#define G4endl
Definition: G4ios.hh:61
#define DBL_MAX
Definition: templates.hh:83

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