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