Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eeToHadronsMultiModel.hh
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.hh 97391 2016-06-02 10:08:45Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eeToHadronsMultiModel
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 18.05.2005
38 //
39 // Modifications:
40 //
41 
42 //
43 // Class Description: vector of e+e- -> hadrons models
44 //
45 
46 // -------------------------------------------------------------------
47 //
48 
49 #ifndef G4eeToHadronsMultiModel_h
50 #define G4eeToHadronsMultiModel_h 1
51 
52 #include "G4VEmModel.hh"
53 #include "G4eeToHadronsModel.hh"
55 #include "G4TrackStatus.hh"
56 #include "Randomize.hh"
59 #include <vector>
60 
61 class G4eeCrossSections;
62 class G4Vee2hadrons;
63 
65 {
66 
67 public:
68 
69  explicit G4eeToHadronsMultiModel(G4int ver=0,
70  const G4String& nam = "eeToHadrons");
71 
72  virtual ~G4eeToHadronsMultiModel();
73 
74  virtual void Initialise(const G4ParticleDefinition*,
75  const G4DataVector&) override;
76 
78  const G4ParticleDefinition*,
79  G4double kineticEnergy,
80  G4double cutEnergy,
81  G4double maxEnergy) override;
82 
84  const G4ParticleDefinition*,
85  G4double kineticEnergy,
87  G4double cutEnergy = 0.0,
88  G4double maxEnergy = DBL_MAX) override;
89 
90  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
91  const G4MaterialCutsCouple*,
92  const G4DynamicParticle*,
93  G4double tmin = 0.0,
94  G4double maxEnergy = DBL_MAX) override;
95 
96  virtual void PrintInfo();
97 
98  // Set the factor to artificially increase the crossSection (default 1)
100 
102  const G4ParticleDefinition*,
103  G4double kineticEnergy,
104  G4double cutEnergy = 0.0,
105  G4double maxEnergy = DBL_MAX);
106 
107 private:
108 
109  void AddEEModel(G4Vee2hadrons*, const G4DataVector&);
110 
111  //change incident e+ kinetic energy into CM total energy(sum of e+ and e-)
112  inline G4double LabToCM(G4double);
113 
114  // hide assignment operator
115  G4eeToHadronsMultiModel & operator=(const G4eeToHadronsMultiModel &right);
117 
118  G4eeCrossSections* cross;
119  G4ParticleChangeForGamma* fParticleChange;
120  G4double delta;
121 
122  std::vector<G4eeToHadronsModel*> models;
123 
124  std::vector<G4double> ekinMin;
125  std::vector<G4double> ekinPeak;
126  std::vector<G4double> ekinMax;
127  std::vector<G4double> cumSum;
128 
129  G4double thKineticEnergy;
130  G4double maxKineticEnergy;
131  G4double csFactor;
132 
133  G4int nModels;
134  G4int verbose;
135  G4bool isInitialised;
136 };
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
140 //change incident e+ kinetic energy into CM total energy(sum of e+ and e-)
141 inline G4double G4eeToHadronsMultiModel::LabToCM(G4double kinE_lab)
142 {
143  G4double totE_CM = 0.0;
145  G4double totE_lab = kinE_lab + mass;
146  totE_CM = std::sqrt(2*mass*(mass+totE_lab));
147 
148  return totE_CM;
149 }
150 
152  const G4ParticleDefinition*,
153  G4double kineticEnergy,
155 {
156  G4double res = 0.0;
157 
158  G4double energy = LabToCM(kineticEnergy);
159 
160  if (energy > thKineticEnergy) {
161  for(G4int i=0; i<nModels; i++) {
162  if(energy >= ekinMin[i] && energy <= ekinMax[i]){
163  res += (models[i])->ComputeCrossSectionPerElectron(0,energy);
164  }
165  cumSum[i] = res;
166  }
167  }
168  return res*csFactor;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
173 #endif
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
int G4int
Definition: G4Types.hh:78
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
static constexpr double electron_mass_c2
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double maxEnergy=DBL_MAX) override
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4eeToHadronsMultiModel(G4int ver=0, const G4String &nam="eeToHadrons")
G4double energy(const ThreeVector &p, const G4double m)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static const G4double fac
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void SetCrossSecFactor(G4double fac)