Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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),
67 csFactor(1.0),
68 nModels(0),
69 verbose(ver),
70 isInitialised(false)
71{
72 thKineticEnergy = DBL_MAX;
73 maxKineticEnergy = 4.521*GeV; //crresponding to 10TeV in lab
74 fParticleChange = nullptr;
75 cross = nullptr;
76 delta = 1.0*MeV; //for bin width
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{
83 delete cross;
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
89 const G4DataVector& cuts)
90{
91 if(!isInitialised) {
92 isInitialised = true;
93
94 //G4cout<<"###Initialise in HadronMultiModel###"<<G4endl;
95
96 cross = new G4eeCrossSections();
97
98 G4eeToTwoPiModel* m2pi =
99 new G4eeToTwoPiModel(cross,maxKineticEnergy,delta);
100 AddEEModel(m2pi,cuts);
101
102 G4eeTo3PiModel* m3pi =
103 new G4eeTo3PiModel(cross,maxKineticEnergy,delta);
104 AddEEModel(m3pi,cuts);
105
106 G4ee2KChargedModel* m2kc =
107 new G4ee2KChargedModel(cross,maxKineticEnergy,delta);
108 AddEEModel(m2kc,cuts);
109
110 G4ee2KNeutralModel* m2kn =
111 new G4ee2KNeutralModel(cross,maxKineticEnergy,delta);
112 AddEEModel(m2kn,cuts);
113
114 G4eeToPGammaModel* mpg1 =
115 new G4eeToPGammaModel(cross,"pi0",maxKineticEnergy,delta);
116 AddEEModel(mpg1,cuts);
117
118 G4eeToPGammaModel* mpg2 =
119 new G4eeToPGammaModel(cross,"eta",maxKineticEnergy,delta);
120 AddEEModel(mpg2,cuts);
121
122 nModels = models.size();
123
124 fParticleChange = GetParticleChangeForGamma();
125 }
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
129
130
131
132void G4eeToHadronsMultiModel::AddEEModel(G4Vee2hadrons* mod,
133 const G4DataVector& cuts)
134{
135 G4eeToHadronsModel* model = new G4eeToHadronsModel(mod, verbose);
136 models.push_back(model);
137 G4double elow = mod->LowEnergy();
138 ekinMin.push_back(elow);
139 if(thKineticEnergy > elow) { thKineticEnergy = elow; }
140 ekinMax.push_back(mod->HighEnergy());
141 ekinPeak.push_back(mod->PeakEnergy());
142 cumSum.push_back(0.0);
143
145 model->Initialise(positron,cuts);
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
151 const G4Material* mat,
152 const G4ParticleDefinition* p,
153 G4double kineticEnergy,
155{
156 return mat->GetElectronDensity()*
157 ComputeCrossSectionPerElectron(p, kineticEnergy);
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
163 const G4ParticleDefinition* p,
164 G4double kineticEnergy,
167{
168 return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
169}
170
171
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173
175 std::vector<G4DynamicParticle*>* newp,
176 const G4MaterialCutsCouple* couple,
177 const G4DynamicParticle* dp,
179{
180 G4double kinEnergy = dp->GetKineticEnergy();
181 G4double energy = LabToCM(kinEnergy);
182 if (energy > thKineticEnergy) {
183 G4double q = cumSum[nModels-1]*G4UniformRand();
184 for(G4int i=0; i<nModels; i++) {
185 if(q <= cumSum[i]) {
186 (models[i])->SampleSecondaries(newp, couple,dp);
187 if(newp->size() > 0) {
188 fParticleChange->ProposeTrackStatus(fStopAndKill);
189 }
190 break;
191 }
192 }
193 }
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197
198void G4eeToHadronsMultiModel::ModelDescription(std::ostream& outFile) const
199{
200 if(verbose > 0) {
201 G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
202 - 2.0*electron_mass_c2;
203 G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
204 - 2.0*electron_mass_c2;
205 outFile << " e+ annihilation into hadrons active from "
206 << e1/GeV << " GeV to " << e2/GeV << " GeV" << G4endl;
207 }
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
213{
214 if(fac > 1.0) {
215 csFactor = fac;
216 if(verbose > 0) {
217 G4cout << "### G4eeToHadronsMultiModel: The cross section for "
218 << "G4eeToHadronsMultiModel is increased by "
219 << csFactor << " times" << G4endl;
220 }
221 }
222}
223
224//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#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:215
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
void ProposeTrackStatus(G4TrackStatus status)
G4double LowEnergy() const
virtual G4double PeakEnergy() const =0
G4double HighEnergy() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void ModelDescription(std::ostream &outFile) const override
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SetCrossSecFactor(G4double fac)
virtual 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")
#define DBL_MAX
Definition: templates.hh:62