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

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 7.3 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.9 2010/10/26 14:15:40 vnivanch Exp $
27// GEANT4 tag $Name: emhighenergy-V09-03-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  thKineticEnergy  = DBL_MAX;
72  maxKineticEnergy = 1.2*GeV;
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77G4eeToHadronsMultiModel::~G4eeToHadronsMultiModel()
78{
79  G4int n = models.size();
80  if(n>0) {
81    for(G4int i=0; i<n; i++) {
82      delete models[i];
83    }
84  }
85  delete cross;
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
90void G4eeToHadronsMultiModel::Initialise(const G4ParticleDefinition*, 
91                                         const G4DataVector&)
92{
93  if(!isInitialised) {
94    isInitialised = true;
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    fParticleChange = GetParticleChangeForGamma();
132  }
133}
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137void G4eeToHadronsMultiModel::AddEEModel(G4Vee2hadrons* mod)
138{
139  G4eeToHadronsModel* model = new G4eeToHadronsModel(mod, verbose);
140  model->SetLowEnergyLimit(LowEnergyLimit());
141  model->SetHighEnergyLimit(HighEnergyLimit());
142  models.push_back(model);
143  G4double elow = mod->ThresholdEnergy();
144  ekinMin.push_back(elow);
145  if(thKineticEnergy > elow) thKineticEnergy = elow;
146  ekinMax.push_back(mod->HighEnergy());
147  ekinPeak.push_back(mod->PeakEnergy());
148  cumSum.push_back(0.0);
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
152
153G4double G4eeToHadronsMultiModel::CrossSectionPerVolume(
154                                      const G4Material* mat,
155                                      const G4ParticleDefinition* p,
156                                      G4double kineticEnergy,
157                                      G4double, G4double)
158{
159  return mat->GetElectronDensity()*
160    ComputeCrossSectionPerElectron(p, kineticEnergy);
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164
165G4double G4eeToHadronsMultiModel::ComputeCrossSectionPerAtom(
166                                      const G4ParticleDefinition* p,
167                                      G4double kineticEnergy,
168                                      G4double Z, G4double,
169                                      G4double, G4double)
170{
171  return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
172}
173
174
175//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
176
177void G4eeToHadronsMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
178                                                const G4MaterialCutsCouple* couple,
179                                                const G4DynamicParticle* dp,
180                                                G4double, G4double)
181{
182  G4double kinEnergy = dp->GetKineticEnergy();
183  if (kinEnergy > 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) fParticleChange->ProposeTrackStatus(fStopAndKill);
189        break;
190      }
191    }
192  }
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
197void G4eeToHadronsMultiModel::PrintInfo()
198{
199  if(verbose > 0) {
200    G4double e1 = 0.5*thKineticEnergy*thKineticEnergy/electron_mass_c2
201      - 2.0*electron_mass_c2; 
202    G4double e2 = 0.5*maxKineticEnergy*maxKineticEnergy/electron_mass_c2
203      - 2.0*electron_mass_c2; 
204    G4cout << "      e+ annihilation into hadrons active from "
205           << e1/GeV << " GeV to " << e2/GeV << " GeV"
206           << G4endl;
207  }
208}
209
210//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211
212void G4eeToHadronsMultiModel::SetCrossSecFactor(G4double fac)
213{
214  if(fac > 1.0) {
215    csFactor = fac;
216    if(verbose > 0)
217      G4cout << "### G4eeToHadronsMultiModel: The cross section for G4eeToHadronsMultiModel "
218             << " is increased by the Factor= " << csFactor << G4endl;
219  }
220}
221
222//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.