source: trunk/source/processes/electromagnetic/highenergy/src/G4eeToHadronsMultiModel.cc@ 969

Last change on this file since 969 was 961, checked in by garnier, 17 years ago

update processes

File size: 7.5 KB
Line 
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,v 1.6 2008/07/11 17:49:11 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4eeToHadronsMultiModel
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 02.08.2004
39//
40// Modifications:
41// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
42// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
43//
44
45//
46// -------------------------------------------------------------------
47//
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
51#include "G4eeToHadronsMultiModel.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
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
62using namespace std;
63
64G4eeToHadronsMultiModel::G4eeToHadronsMultiModel(G4int ver, const G4String& name)
65 : G4VEmModel(name),
66 csFactor(1.0),
67 nModels(0),
68 verbose(ver),
69 isInitialised(false)
70{}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
74G4eeToHadronsMultiModel::~G4eeToHadronsMultiModel()
75{
76 G4int n = models.size();
77 if(n>0) {
78 for(G4int i=0; i<n; i++) {
79 delete models[i];
80 }
81 }
82 delete cross;
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
87void G4eeToHadronsMultiModel::Initialise(const G4ParticleDefinition*,
88 const G4DataVector&)
89{
90 if(!isInitialised) {
91 isInitialised = true;
92
93 thKineticEnergy = DBL_MAX;
94 maxKineticEnergy = 1.2*GeV;
95
96 cross = new G4eeCrossSections();
97
98 G4eeToTwoPiModel* m2pi = new G4eeToTwoPiModel(cross);
99 m2pi->SetHighEnergy(maxKineticEnergy);
100 AddEEModel(m2pi);
101
102 G4eeTo3PiModel* m3pi1 = new G4eeTo3PiModel(cross);
103 m3pi1->SetHighEnergy(0.95*GeV);
104 AddEEModel(m3pi1);
105
106 G4eeTo3PiModel* m3pi2 = new G4eeTo3PiModel(cross);
107 m3pi2->SetLowEnergy(0.95*GeV);
108 m3pi2->SetHighEnergy(maxKineticEnergy);
109 AddEEModel(m3pi2);
110
111 G4ee2KChargedModel* m2kc = new G4ee2KChargedModel(cross);
112 m2kc->SetHighEnergy(maxKineticEnergy);
113 AddEEModel(m2kc);
114
115 G4ee2KNeutralModel* m2kn = new G4ee2KNeutralModel(cross);
116 m2kn->SetHighEnergy(maxKineticEnergy);
117 AddEEModel(m2kn);
118
119 G4eeToPGammaModel* mpg1 = new G4eeToPGammaModel(cross,"pi0");
120 mpg1->SetLowEnergy(0.7*GeV);
121 mpg1->SetHighEnergy(maxKineticEnergy);
122 AddEEModel(mpg1);
123
124 G4eeToPGammaModel* mpg2 = new G4eeToPGammaModel(cross,"eta");
125 mpg2->SetLowEnergy(0.7*GeV);
126 mpg2->SetHighEnergy(maxKineticEnergy);
127 AddEEModel(mpg2);
128
129 nModels = models.size();
130
131 if(pParticleChange) {
132 fParticleChange =
133 reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
134 } else {
135 fParticleChange = new G4ParticleChangeForGamma();
136 }
137 }
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
142void G4eeToHadronsMultiModel::AddEEModel(G4Vee2hadrons* mod)
143{
144 G4eeToHadronsModel* model = new G4eeToHadronsModel(mod, verbose);
145 model->SetLowEnergyLimit(LowEnergyLimit());
146 model->SetHighEnergyLimit(HighEnergyLimit());
147 models.push_back(model);
148 G4double elow = mod->ThresholdEnergy();
149 ekinMin.push_back(elow);
150 if(thKineticEnergy > elow) thKineticEnergy = elow;
151 ekinMax.push_back(mod->HighEnergy());
152 ekinPeak.push_back(mod->PeakEnergy());
153 cumSum.push_back(0.0);
154}
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
157
158G4double G4eeToHadronsMultiModel::CrossSectionPerVolume(
159 const G4Material* mat,
160 const G4ParticleDefinition* p,
161 G4double kineticEnergy,
162 G4double, G4double)
163{
164 return mat->GetElectronDensity()*
165 ComputeCrossSectionPerElectron(p, kineticEnergy);
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169
170G4double G4eeToHadronsMultiModel::ComputeCrossSectionPerAtom(
171 const G4ParticleDefinition* p,
172 G4double kineticEnergy,
173 G4double Z, G4double,
174 G4double, G4double)
175{
176 return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
177}
178
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181
182void G4eeToHadronsMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
183 const G4MaterialCutsCouple* couple,
184 const G4DynamicParticle* dp,
185 G4double, G4double)
186{
187 G4double kinEnergy = dp->GetKineticEnergy();
188 if (kinEnergy > thKineticEnergy) {
189 G4double q = cumSum[nModels-1]*G4UniformRand();
190 for(G4int i=0; i<nModels; i++) {
191 if(q <= cumSum[i]) {
192 (models[i])->SampleSecondaries(newp, couple,dp);
193 if(newp->size() > 0) fParticleChange->ProposeTrackStatus(fStopAndKill);
194 break;
195 }
196 }
197 }
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
201
202void G4eeToHadronsMultiModel::PrintInfo()
203{
204 if(verbose > 0) {
205 G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
206 - 2.0*electron_mass_c2;
207 G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
208 - 2.0*electron_mass_c2;
209 G4cout << " e+ annihilation into hadrons active from "
210 << e1/GeV << " GeV to " << e2/GeV << " GeV"
211 << G4endl;
212 }
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
216
217void G4eeToHadronsMultiModel::SetCrossSecFactor(G4double fac)
218{
219 if(fac > 1.0) {
220 csFactor = fac;
221 if(verbose > 0)
222 G4cout << "### G4eeToHadronsMultiModel: The cross section for G4eeToHadronsMultiModel "
223 << " is increased by the Factor= " << csFactor << G4endl;
224 }
225}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.