Geant4-11
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4eeToHadronsMultiModel
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 18.05.2005
37//
38// Modifications:
39//
40
41//
42// Class Description: vector of e+e- -> hadrons models
43//
44
45// -------------------------------------------------------------------
46//
47
48#ifndef G4eeToHadronsMultiModel_h
49#define G4eeToHadronsMultiModel_h 1
50
51#include "G4VEmModel.hh"
52#include "G4eeToHadronsModel.hh"
54#include "G4TrackStatus.hh"
55#include "Randomize.hh"
58#include <vector>
59
61class G4Vee2hadrons;
62
64{
65
66public:
67
68 explicit G4eeToHadronsMultiModel(G4int ver=0,
69 const G4String& nam = "eeToHadrons");
70
71 ~G4eeToHadronsMultiModel() override;
72
73 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
74
77 G4double kineticEnergy,
78 G4double cutEnergy,
79 G4double maxEnergy) override;
80
83 G4double kineticEnergy,
85 G4double cutEnergy = 0.0,
86 G4double maxEnergy = DBL_MAX) override;
87
88 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
90 const G4DynamicParticle*,
91 G4double tmin = 0.0,
92 G4double maxEnergy = DBL_MAX) override;
93
94 void ModelDescription(std::ostream& outFile) const override;
95
96 // Set the factor to artificially increase the crossSection (default 1)
98
100 G4double kineticEnergy,
101 G4double cutEnergy = 0.0,
102 G4double maxEnergy = DBL_MAX);
103
104 // hide assignment operator
107
108private:
109
111
112 //change incident e+ kinetic energy into CM total energy(sum of e+ and e-)
113 inline G4double LabToCM(G4double);
114
117
118 std::vector<G4eeToHadronsModel*> models;
119
120 std::vector<G4double> ekinMin;
121 std::vector<G4double> ekinPeak;
122 std::vector<G4double> ekinMax;
123 std::vector<G4double> cumSum;
124
129
133};
134
135//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
136
137//change incident e+ kinetic energy into CM total energy(sum of e+ and e-)
139{
140 G4double totE_CM = 0.0;
142 G4double totE_lab = kinE_lab + mass;
143 totE_CM = std::sqrt(2*mass*(mass+totE_lab));
144
145 return totE_CM;
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
149
150#endif
static const G4double fac
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void ModelDescription(std::ostream &outFile) const override
G4eeToHadronsMultiModel & operator=(const G4eeToHadronsMultiModel &right)=delete
std::vector< G4eeToHadronsModel * > models
void AddEEModel(G4Vee2hadrons *, const G4DataVector &)
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)
G4eeToHadronsMultiModel(const G4eeToHadronsMultiModel &)=delete
std::vector< G4double > ekinPeak
std::vector< G4double > cumSum
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy) override
std::vector< G4double > ekinMin
std::vector< G4double > ekinMax
G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4eeToHadronsMultiModel(G4int ver=0, const G4String &nam="eeToHadrons")
G4ParticleChangeForGamma * fParticleChange
static constexpr double electron_mass_c2
#define DBL_MAX
Definition: templates.hh:62