Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eeToTwoGammaModel.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: G4eeToTwoGammaModel
33//
34// Author: Vladimir Ivanchenko on base of Michel Maire code
35//
36// Creation date: 02.08.2004
37//
38// Modifications:
39// 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
40// 18-04-05 Compute CrossSectionPerVolume (V.Ivanchenko)
41// 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
42// 20-10-06 Add theGamma as a member (V.Ivanchenko)
43//
44//
45// Class Description:
46//
47// Implementation of e+ annihilation into 2 gamma on fly
48// Annihilation at rest is sampled by G4VPositronAtRestModel
49//
50// -------------------------------------------------------------------
51//
52
53#ifndef G4eeToTwoGammaModel_h
54#define G4eeToTwoGammaModel_h 1
55
56#include "G4VEmModel.hh"
57#include <vector>
58
60
62{
63
64public:
65
66 explicit G4eeToTwoGammaModel(const G4ParticleDefinition* p = nullptr,
67 const G4String& nam = "eplus2gg");
68
70
71 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
72
74
77 G4double kinEnergy,
78 G4double Z,
79 G4double A = 0.,
80 G4double cutEnergy = 0.,
81 G4double maxEnergy = DBL_MAX) override;
82
85 G4double kineticEnergy,
86 G4double cutEnergy = 0.0,
87 G4double maxEnergy = DBL_MAX) override;
88
89 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
91 const G4DynamicParticle*,
92 G4double tmin,
93 G4double maxEnergy) override;
94
95 // hide assignment operator
98
99private:
100
101 G4double pi_rcl2;
102 const G4ParticleDefinition* theGamma;
103 G4ParticleChangeForGamma* fParticleChange;
104};
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
108#endif
double G4double
Definition G4Types.hh:83
const G4double A[17]
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
G4eeToTwoGammaModel & operator=(const G4eeToTwoGammaModel &right)=delete
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4eeToTwoGammaModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="eplus2gg")
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
G4eeToTwoGammaModel(const G4eeToTwoGammaModel &)=delete
~G4eeToTwoGammaModel() override
#define DBL_MAX
Definition templates.hh:62