Geant4_10
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 66241 2012-12-13 18:34:42Z gunter $
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"
53 #include "G4Proton.hh"
54 #include "G4LossTableManager.hh"
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57 
58 using namespace std;
59 
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 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73 
75 {}
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80 {
81  return (p.GetPDGCharge() != 0.0 && !p.IsShortLived());
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
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 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
148 
150  const G4Material* mat)
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 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167 
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 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
tuple a
Definition: test.py:11
virtual G4bool IsApplicable(const G4ParticleDefinition &p)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void SetBuildTableFlag(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
const char * p
Definition: xmltok.h:285
void SetSplineFlag(G4bool val)
G4double MaxKinEnergy() const
G4double FactorForAngleLimit() const
const XML_Char * name
Definition: expat.h:151
void SetStartFromNullFlag(G4bool val)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
Float_t mat
Definition: plot.C:40
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
tuple degree
Definition: hepunit.py:69
G4double PolarAngleLimit() const
virtual void InitialiseProcess(const G4ParticleDefinition *)
bool G4bool
Definition: G4Types.hh:79
G4double GetInvA23() const
TString part[npart]
Definition: Style.C:32
const XML_Char XML_Content * model
Definition: expat.h:151
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
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
#define G4endl
Definition: G4ios.hh:61
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
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4CoulombScattering(const G4String &name="CoulombScat")