Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eeToHadronsMultiModel.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: G4eeToHadronsMultiModel.cc 97391 2016-06-02 10:08:45Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eeToHadronsMultiModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 02.08.2004
38 //
39 // Modifications:
40 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
41 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
42 //
43 
44 //
45 // -------------------------------------------------------------------
46 //
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
51 #include "G4PhysicalConstants.hh"
52 #include "G4SystemOfUnits.hh"
53 #include "G4eeToTwoPiModel.hh"
54 #include "G4eeTo3PiModel.hh"
55 #include "G4eeToPGammaModel.hh"
56 #include "G4ee2KNeutralModel.hh"
57 #include "G4ee2KChargedModel.hh"
58 #include "G4eeCrossSections.hh"
59 #include "G4Vee2hadrons.hh"
60 #include "G4Positron.hh"
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
64 using namespace std;
65 
67  const G4String& mname) : G4VEmModel(mname),
68  csFactor(1.0),
69  nModels(0),
70  verbose(ver),
71  isInitialised(false)
72 {
73  thKineticEnergy = DBL_MAX;
74  maxKineticEnergy = 4.521*GeV; //crresponding to 10TeV in lab
75  fParticleChange = nullptr;
76  cross = nullptr;
77  delta = 1.0*MeV; //for bin width
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83 {
84  delete cross;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
90  const G4DataVector& cuts)
91 {
92  if(!isInitialised) {
93  isInitialised = true;
94 
95  //G4cout<<"###Initialise in HadronMultiModel###"<<G4endl;
96 
97  cross = new G4eeCrossSections();
98 
99  G4eeToTwoPiModel* m2pi =
100  new G4eeToTwoPiModel(cross,maxKineticEnergy,delta);
101  AddEEModel(m2pi,cuts);
102 
103  G4eeTo3PiModel* m3pi =
104  new G4eeTo3PiModel(cross,maxKineticEnergy,delta);
105  AddEEModel(m3pi,cuts);
106 
107  G4ee2KChargedModel* m2kc =
108  new G4ee2KChargedModel(cross,maxKineticEnergy,delta);
109  AddEEModel(m2kc,cuts);
110 
111  G4ee2KNeutralModel* m2kn =
112  new G4ee2KNeutralModel(cross,maxKineticEnergy,delta);
113  AddEEModel(m2kn,cuts);
114 
115  G4eeToPGammaModel* mpg1 =
116  new G4eeToPGammaModel(cross,"pi0",maxKineticEnergy,delta);
117  AddEEModel(mpg1,cuts);
118 
119  G4eeToPGammaModel* mpg2 =
120  new G4eeToPGammaModel(cross,"eta",maxKineticEnergy,delta);
121  AddEEModel(mpg2,cuts);
122 
123  nModels = models.size();
124 
125  fParticleChange = GetParticleChangeForGamma();
126  }
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 
132 
133 void G4eeToHadronsMultiModel::AddEEModel(G4Vee2hadrons* mod,
134  const G4DataVector& cuts)
135 {
136  G4eeToHadronsModel* model = new G4eeToHadronsModel(mod, verbose);
137  models.push_back(model);
138  G4double elow = mod->LowEnergy();
139  ekinMin.push_back(elow);
140  if(thKineticEnergy > elow) { thKineticEnergy = elow; }
141  ekinMax.push_back(mod->HighEnergy());
142  ekinPeak.push_back(mod->PeakEnergy());
143  cumSum.push_back(0.0);
144 
146  model->Initialise(positron,cuts);
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150 
152  const G4Material* mat,
153  const G4ParticleDefinition* p,
154  G4double kineticEnergy,
156 {
157  return mat->GetElectronDensity()*
158  ComputeCrossSectionPerElectron(p, kineticEnergy);
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
164  const G4ParticleDefinition* p,
165  G4double kineticEnergy,
168 {
169  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
170 }
171 
172 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
176  std::vector<G4DynamicParticle*>* newp,
177  const G4MaterialCutsCouple* couple,
178  const G4DynamicParticle* dp,
180 {
181  G4double kinEnergy = dp->GetKineticEnergy();
182  G4double energy = LabToCM(kinEnergy);
183  if (energy > thKineticEnergy) {
184  G4double q = cumSum[nModels-1]*G4UniformRand();
185  for(G4int i=0; i<nModels; i++) {
186  if(q <= cumSum[i]) {
187  (models[i])->SampleSecondaries(newp, couple,dp);
188  if(newp->size() > 0) {
189  fParticleChange->ProposeTrackStatus(fStopAndKill);
190  }
191  break;
192  }
193  }
194  }
195 }
196 
197 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
198 
200 {
201  if(verbose > 0) {
202  G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
203  - 2.0*electron_mass_c2;
204  G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
205  - 2.0*electron_mass_c2;
206  G4cout << " e+ annihilation into hadrons active from "
207  << e1/GeV << " GeV to " << e2/GeV << " GeV"
208  << G4endl;
209  }
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213 
215 {
216  if(fac > 1.0) {
217  csFactor = fac;
218  if(verbose > 0) {
219  G4cout << "### G4eeToHadronsMultiModel: The cross section for "
220  << "G4eeToHadronsMultiModel is increased by "
221  << csFactor << " times" << G4endl;
222  }
223  }
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double LowEnergy() const
G4double GetKineticEnergy() const
const char * p
Definition: xmltok.h:285
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
static constexpr double electron_mass_c2
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
virtual G4double PeakEnergy() const =0
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4eeToHadronsMultiModel(G4int ver=0, const G4String &nam="eeToHadrons")
G4double energy(const ThreeVector &p, const G4double m)
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static const G4double fac
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double HighEnergy() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
#define DBL_MAX
Definition: templates.hh:83
const XML_Char XML_Content * model
Definition: expat.h:151
void SetCrossSecFactor(G4double fac)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132