Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CoherentPairProduction.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// Author: Alexei Sytov
27// Co-author: Gianfranco Paterno (testing)
28// Using the key points of G4BaierKatkov and developments of V.V. Tikhomirov,
29// partially described in L. Bandiera et al. Eur. Phys. J. C 82, 699 (2022)
30
31#ifndef G4CoherentPairProduction_h
32#define G4CoherentPairProduction_h 1
33
34#include "G4VDiscreteProcess.hh"
35
36#include <vector>
40
42#include "G4LogicalVolume.hh"
43#include "G4ParticleTable.hh"
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46
48{
49public:
50 G4CoherentPairProduction(const G4String& processName = "cpp",
52
54
55 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
56
58 {
59 return(aPD.GetParticleName() == "gamma");
60 }
61
62 // print documentation in html format
63 void ProcessDescription(std::ostream&) const override;
64
65 ///special functions
66 void Input(const G4Material* crystal,
67 const G4String &lattice)
68 {Input(crystal,lattice,"");}
69
70 void Input(const G4Material* crystal,
71 const G4String &lattice,
72 const G4String &filePath);
73
74 // an option to use crystal data already created outside this class
75 void Input(const G4ChannelingFastSimCrystalData* crystalData);
76
77 ///activate incoherent scattering
78 ///(standard gamma conversion should be switched off in physics list)
79 void ActivateIncoherentScattering(){fIncoherentScattering = true;}
80
82
83 ///get cuts
84 // minimal energy for non-zero cross section
85 G4double ModelMinPrimaryEnergy() { return fLowEnergyLimit;}
86 G4double GetHighAngleLimit() {return fHighAngleLimit;}
87 G4double GetPPKineticEnergyCut() {return fPPKineticEnergyCut;}
88
89 /// get the number of pairs in sampling of Baier-Katkov Integral
90 /// (MC integration by e+- energy and angles <=> e+- momentum)
91 G4int GetSamplingPairsNumber(){return fNMCPairs;}
92
93 /// get the number of particle angles 1/gamma in pair production
94 /// defining the width of the angular distribution of pair sampling
95 /// in the Baier-Katkov Integral
96 G4double GetChargeParticleAngleFactor(){return fChargeParticleAngleFactor;}
97
98 /// get number of trajectory steps of a single particle (e- or e+)
99 G4double GetNTrajectorySteps(){return fNTrajectorySteps;}
100
101 /// get effective radiation length
102 /// (due to coherent process of pair production)
103 /// simulated for the current photon
104 G4double GetEffectiveLrad(){return fEffectiveLrad;}
105
106 ///get the name of G4Region in which the model is applicable
107 G4String GetG4RegionName() {return fG4RegionName;}
108
109 ///set cuts
110 void SetLowEnergyLimit(G4double energy){fLowEnergyLimit=energy;}
111 void SetHighAngleLimit(G4double angle) {fHighAngleLimit=angle;}
112 void SetPPKineticEnergyCut(G4double kineticEnergyCut) {fPPKineticEnergyCut=kineticEnergyCut;}
113
114 /// set the number of pairs in sampling of Baier-Katkov Integral
115 /// (MC integration by e+- energy and angles <=> e+- momentum)
116 void SetSamplingPairsNumber(G4int nPairs){fNMCPairs = nPairs;}
117
118 /// set the number of particle angles 1/gamma in pair production
119 /// defining the width of the angular distribution of pair sampling
120 /// in the Baier-Katkov Integral
121 void SetChargeParticleAngleFactor(G4double chargeParticleAngleFactor)
122 {fChargeParticleAngleFactor = chargeParticleAngleFactor;}
123
124 /// set number of trajectory steps of a single particle (e- or e+)
125 void SetNTrajectorySteps(G4int nTrajectorySteps)
126 {fNTrajectorySteps = nTrajectorySteps;}
127
128 ///set the name of G4Region in which the model is applicable
129 void SetG4RegionName(const G4String& nameG4Region){fG4RegionName=nameG4Region;}
130
131 G4double GetMeanFreePath(const G4Track& aTrack,
132 G4double,
133 G4ForceCondition* condition) override;
134
135private:
136
137 G4int FindVectorIndex(std::vector<G4double> &myvector, G4double value);
138
139 G4ChannelingFastSimCrystalData* fCrystalData{nullptr};
140
141 //collection of etotal
142 std::vector <CLHEP::Hep2Vector> fullVectorEtotal;
143
144 //collection of x
145 std::vector <CLHEP::Hep2Vector> fullVectorX;
146
147 //collection of y
148 std::vector <CLHEP::Hep2Vector> fullVectorY;
149
150 //collection of tx
151 std::vector <CLHEP::Hep2Vector> fullVectorTX;
152
153 //collection of tx
154 std::vector <CLHEP::Hep2Vector> fullVectorTY;
155
156 //the vector of the discrete CDF of the production of sampling e+e- pairs
157 //(in reality per distance along the photon direction)
158 std::vector <G4double> fPairProductionCDFdz;
159
160 G4double fLowEnergyLimit = 1*CLHEP::GeV;
161 G4double fHighAngleLimit = 50*CLHEP::mrad;
162
163 ///minimal kinetic energy of a charged particle produced
164 G4double fPPKineticEnergyCut = 1*CLHEP::MeV;
165
166 ///Monte Carlo statistics of e+- pair sampling in Baier-Katkov for 1 photon
167 G4int fNMCPairs = 150;
168
169 G4double fChargeParticleAngleFactor = 4; // number of particle angles 1/gamma:
170 // more fChargeParticleAngleFactor => higher paramParticleAngle
171
172 ///number of trajectory steps of a single particle (e- or e+)
173 G4int fNTrajectorySteps=250;
174
175 ///effective radiation length (due to coherent process of pair production)
176 G4double fEffectiveLrad = 0.;
177
178 ///the name of G4Region in which the model is applicable
179 G4String fG4RegionName = "Crystal";
180
181 ///charged particle mass
182 const G4double fMass = CLHEP::electron_mass_c2;
183
184 ///flag of simulation of incoherent scattering
185 G4bool fIncoherentScattering = false;
186
187};
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190
191#endif
192
Definition of the G4ChannelingFastSimCrystalData class The class inherits G4VChannelingFastSimCrystal...
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4CoherentPairProduction(const G4String &processName="cpp", G4ProcessType aType=fElectromagnetic)
G4double GetNTrajectorySteps()
get number of trajectory steps of a single particle (e- or e+)
G4String GetG4RegionName()
get the name of G4Region in which the model is applicable
void SetPPKineticEnergyCut(G4double kineticEnergyCut)
void SetLowEnergyLimit(G4double energy)
set cuts
void SetG4RegionName(const G4String &nameG4Region)
set the name of G4Region in which the model is applicable
G4bool IsApplicable(const G4ParticleDefinition &aPD) override
void Input(const G4Material *crystal, const G4String &lattice)
special functions
~G4CoherentPairProduction()=default
void SetNTrajectorySteps(G4int nTrajectorySteps)
set number of trajectory steps of a single particle (e- or e+)
G4ChannelingFastSimCrystalData * GetCrystalData()
void SetChargeParticleAngleFactor(G4double chargeParticleAngleFactor)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *condition) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void ProcessDescription(std::ostream &) const override
const G4String & GetParticleName() const
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)