Geant4  10.00.p02
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 66241 2012-12-13 18:34:42Z gunter $
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 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
63 using namespace std;
64 
66  : G4VEmModel(mname),
67  csFactor(1.0),
68  nModels(0),
69  verbose(ver),
70  isInitialised(false)
71 {
73  maxKineticEnergy = 1.2*GeV;
74  fParticleChange = 0;
75  cross = 0;
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81 {
82  G4int n = models.size();
83  if(n>0) {
84  for(G4int i=0; i<n; i++) {
85  delete models[i];
86  }
87  }
88  delete cross;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94  const G4DataVector&)
95 {
96  if(!isInitialised) {
97  isInitialised = true;
98 
99  cross = new G4eeCrossSections();
100 
103  AddEEModel(m2pi);
104 
105  G4eeTo3PiModel* m3pi1 = new G4eeTo3PiModel(cross);
106  m3pi1->SetHighEnergy(0.95*GeV);
107  AddEEModel(m3pi1);
108 
109  G4eeTo3PiModel* m3pi2 = new G4eeTo3PiModel(cross);
110  m3pi2->SetLowEnergy(0.95*GeV);
112  AddEEModel(m3pi2);
113 
116  AddEEModel(m2kc);
117 
120  AddEEModel(m2kn);
121 
122  G4eeToPGammaModel* mpg1 = new G4eeToPGammaModel(cross,"pi0");
123  mpg1->SetLowEnergy(0.7*GeV);
125  AddEEModel(mpg1);
126 
127  G4eeToPGammaModel* mpg2 = new G4eeToPGammaModel(cross,"eta");
128  mpg2->SetLowEnergy(0.7*GeV);
130  AddEEModel(mpg2);
131 
132  nModels = models.size();
133 
135  }
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
141 {
142  G4eeToHadronsModel* model = new G4eeToHadronsModel(mod, verbose);
145  models.push_back(model);
146  G4double elow = mod->ThresholdEnergy();
147  ekinMin.push_back(elow);
148  if(thKineticEnergy > elow) thKineticEnergy = elow;
149  ekinMax.push_back(mod->HighEnergy());
150  ekinPeak.push_back(mod->PeakEnergy());
151  cumSum.push_back(0.0);
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
157  const G4Material* mat,
158  const G4ParticleDefinition* p,
159  G4double kineticEnergy,
161 {
162  return mat->GetElectronDensity()*
163  ComputeCrossSectionPerElectron(p, kineticEnergy);
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
167 
169  const G4ParticleDefinition* p,
170  G4double kineticEnergy,
171  G4double Z, G4double,
173 {
174  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
175 }
176 
177 
178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179 
180 void G4eeToHadronsMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
181  const G4MaterialCutsCouple* couple,
182  const G4DynamicParticle* dp,
184 {
185  G4double kinEnergy = dp->GetKineticEnergy();
186  if (kinEnergy > thKineticEnergy) {
188  for(G4int i=0; i<nModels; i++) {
189  if(q <= cumSum[i]) {
190  (models[i])->SampleSecondaries(newp, couple,dp);
191  if(newp->size() > 0) fParticleChange->ProposeTrackStatus(fStopAndKill);
192  break;
193  }
194  }
195  }
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199 
201 {
202  if(verbose > 0) {
203  G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
204  - 2.0*electron_mass_c2;
205  G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
206  - 2.0*electron_mass_c2;
207  G4cout << " e+ annihilation into hadrons active from "
208  << e1/GeV << " GeV to " << e2/GeV << " GeV"
209  << G4endl;
210  }
211 }
212 
213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214 
216 {
217  if(fac > 1.0) {
218  csFactor = fac;
219  if(verbose > 0)
220  G4cout << "### G4eeToHadronsMultiModel: The cross section for G4eeToHadronsMultiModel "
221  << " is increased by the Factor= " << csFactor << G4endl;
222  }
223 }
224 
225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static const G4double fac
G4double GetKineticEnergy() const
void SetLowEnergy(G4double val)
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
static const G4double e2
void SetHighEnergy(G4double val)
std::vector< G4double > ekinMax
int G4int
Definition: G4Types.hh:78
std::vector< G4double > ekinPeak
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
virtual G4double ThresholdEnergy() const =0
std::vector< G4eeToHadronsModel * > models
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:215
std::vector< G4double > ekinMin
static const double GeV
Definition: G4SIunits.hh:196
std::vector< G4double > cumSum
const G4int n
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4ParticleChangeForGamma * fParticleChange
static const G4double e1
virtual G4double PeakEnergy() const =0
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4eeToHadronsMultiModel(G4int ver=0, const G4String &nam="eeToHadrons")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
#define G4endl
Definition: G4ios.hh:61
G4double HighEnergy() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void AddEEModel(G4Vee2hadrons *)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
#define DBL_MAX
Definition: templates.hh:83
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX)
void SetCrossSecFactor(G4double fac)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121