Geant4  10.00.p02
G4CoulombScattering.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4CoulombScattering.cc 81865 2014-06-06 11:32:58Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4CoulombScattering
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 22.08.2004
38 //
39 // Modifications:
40 // 01.08.06 V.Ivanchenko add choice between G4eCoulombScatteringModel and
41 // G4CoulombScatteringModel
42 //
43 
44 //
45 // -------------------------------------------------------------------
46 //
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
50 #include "G4CoulombScattering.hh"
51 #include "G4SystemOfUnits.hh"
54 #include "G4Proton.hh"
55 #include "G4LossTableManager.hh"
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58 
59 using namespace std;
60 
62  : G4VEmProcess(name),q2Max(TeV*TeV),isInitialised(false)
63 {
64  // G4cout << "G4CoulombScattering constructor "<< G4endl;
65  SetBuildTableFlag(true);
66  SetStartFromNullFlag(false);
67  SetIntegral(true);
70  SetSplineFlag(true);
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
76 {}
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81 {
82  return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88 {
89  // second initialisation not allowed for the time being
90  // this means that polar angle limit change will not be appled
91  // after first initialisation
92  if(isInitialised) { return; }
93 
95  *CLHEP::hbarc/CLHEP::fermi;
96  q2Max = 0.5*a*a;
97  G4double theta = PolarAngleLimit();
98 
99  // restricted or non-restricted cross section table
100  G4bool yes = false;
101  if(theta == CLHEP::pi) { yes = true; }
103  /*
104  G4cout << "### G4CoulombScattering::InitialiseProcess: "
105  << p->GetParticleName()
106  << " Emin(MeV)= " << MinKinEnergy()/MeV
107  << " Emax(TeV)= " << MaxKinEnergy()/TeV
108  << " nbins= " << LambdaBinning()
109  << " theta= " << theta
110  << G4endl;
111  */
112  /*
113  // second initialisation
114  if(isInitialised) {
115  G4VEmModel* mod = EmModel(1);
116  mod->SetPolarAngleLimit(theta);
117  mod = GetModelByIndex(1);
118  if(mod) { mod->SetPolarAngleLimit(theta); }
119 
120  // first initialisation
121  } else {
122  */
123 
124  isInitialised = true;
125  G4double mass = p->GetPDGMass();
127  //G4cout << name << " type: " << p->GetParticleType()
128  //<< " mass= " << mass << G4endl;
129  yes = true;
130  if (mass > GeV || p->GetParticleType() == "nucleus") {
131  SetBuildTableFlag(false);
132  yes = false;
133  if(name != "GenericIon") { SetVerboseLevel(0); }
134  } else {
135  if(name != "e-" && name != "e+" &&
136  name != "mu+" && name != "mu-" && name != "pi+" &&
137  name != "kaon+" && name != "proton" ) { SetVerboseLevel(0); }
138  }
139 
140  if(!EmModel(1)) {
141  if(yes) { SetEmModel(new G4eCoulombScatteringModel(), 1); }
142  else { SetEmModel(new G4IonCoulombScatteringModel(), 1); }
143  }
144  G4VEmModel* model = EmModel(1);
145  G4double emin = std::max(MinKinEnergy(),model->LowEnergyLimit());
146  G4double emax = std::min(MaxKinEnergy(),model->HighEnergyLimit());
147  model->SetPolarAngleLimit(theta);
148  model->SetLowEnergyLimit(emin);
149  model->SetHighEnergyLimit(emax);
150  AddEmModel(1, model);
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154 
156  const G4Material* mat)
157 {
158  // Pure Coulomb scattering
159  G4double emin = 0.0;
160 
161  // Coulomb scattering combined with multiple or hadronic scattering
162  G4double theta = PolarAngleLimit();
163  if(0.0 < theta) {
164  G4double p2 = q2Max*mat->GetIonisation()->GetInvA23()/(1.0 - cos(theta));
165  G4double mass = part->GetPDGMass();
166  emin = sqrt(p2 + mass*mass) - mass;
167  }
168 
169  return emin;
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173 
175 {
176  G4cout << " " << PolarAngleLimit()/degree
177  << " < Theta(degree) < 180";
178 
179  if(q2Max < DBL_MAX) { G4cout << "; pLimit(GeV^1)= " << sqrt(q2Max)/GeV; }
180  G4cout << G4endl;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void SetBuildTableFlag(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
G4String name
Definition: TRTMaterials.hh:40
const G4double pi
void SetSplineFlag(G4bool val)
G4double MaxKinEnergy() const
G4double FactorForAngleLimit() const
G4double a
Definition: TRTMaterials.hh:39
void SetStartFromNullFlag(G4bool val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
G4double PolarAngleLimit() const
virtual void InitialiseProcess(const G4ParticleDefinition *)
bool G4bool
Definition: G4Types.hh:79
G4double GetInvA23() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static const double GeV
Definition: G4SIunits.hh:196
const G4String & GetParticleType() const
void SetIntegral(G4bool val)
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
void SetSecondaryParticle(const G4ParticleDefinition *p)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const double degree
Definition: G4SIunits.hh:125
#define G4endl
Definition: G4ios.hh:61
static const double TeV
Definition: G4SIunits.hh:197
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4double GetPDGCharge() const
G4double MinKinEnergy() const
#define DBL_MAX
Definition: templates.hh:83
static const double fermi
Definition: G4SIunits.hh:93
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4CoulombScattering(const G4String &name="CoulombScat")