Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VEmAngularDistribution.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// File name: G4VEmAngularDistribution
32//
33// Author: V. Ivanchenko using design of existing
34// interface G4VBremAngularDistribution
35//
36// Creation date: 13 October 2010
37//
38// Modifications:
39//
40// Class Description:
41//
42// Abstract base class for polar angle sampling
43//
44// Class Description: End
45
46// -------------------------------------------------------------------
47//
48
49#ifndef G4VEmAngularDistribution_h
50#define G4VEmAngularDistribution_h 1
51
52#include "globals.hh"
53#include "G4ThreeVector.hh"
54#include "G4DynamicParticle.hh"
55
56class G4Material;
57
59{
60public:
61
62 explicit G4VEmAngularDistribution(const G4String& name);
63
65
66 // Sample direction in global coordinate system,
67 // this means for zero scattering angle this direction is the same
68 // as the direction of primary
70 G4double finalTotalEnergy,
71 G4int Z,
72 const G4Material*) = 0;
73
74 // Sample direction in global coordinate system for given electronic shell,
75 // this means for zero scattering angle this direction is the same
76 // as the direction of primary
78 const G4DynamicParticle* dp,
79 G4double finalTotalEnergy,
80 G4int Z,
81 G4int shellID,
82 const G4Material*);
83
84 virtual void SamplePairDirections(const G4DynamicParticle* dp,
85 G4double elecKinEnergy,
86 G4double posiKinEnergy,
87 G4ThreeVector& dirElectron,
88 G4ThreeVector& dirPositron,
89 G4int Z = 0,
90 const G4Material* mat = nullptr);
91
92 virtual void PrintGeneratorInformation() const;
93
94 inline const G4String& GetName() const;
95
96 // hide assignment operator
98 operator=(const G4VEmAngularDistribution &right) = delete;
100
101protected:
102
105
106private:
107
108 G4String fName;
109};
110
112{
113 return fName;
114}
115
116#endif
117
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4String & GetName() const
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution(const G4String &name)
G4VEmAngularDistribution(const G4VEmAngularDistribution &)=delete
G4VEmAngularDistribution & operator=(const G4VEmAngularDistribution &right)=delete
virtual void PrintGeneratorInformation() const
virtual void SamplePairDirections(const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
virtual ~G4VEmAngularDistribution()
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0