Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CoherentPairProductionPhysics.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
28#ifndef G4CoherentPairProductionPhysics_h
29#define G4CoherentPairProductionPhysics_h 1
30
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35
37{
38public:
40 "Coherent Pair Production Physics");
41
43
44 void ConstructParticle() override;
45 void ConstructProcess() override;
46
47 ///activate incoherent scattering
48 ///(standard gamma conversion should be switched off in physics list)
49 void ActivateIncoherentScattering(){fIncoherentScattering = true;}
50
51 ///set functions
52
53 ///set name of G4ChannelingFastSimModel from which an auto input should be performed
54 void SetNameChannelingModel(const G4String& nameChannelingModel)
55 {fNameChannelingModel=nameChannelingModel;}
56
57 ///set name of G4Region to where the G4ChannelingFastSimModel is active
58 void SetNameG4Region(const G4String& nameG4Region)
59 {fNameRegion=nameG4Region;}
60
61 void SetLowEnergyLimit(G4double energy){fLowEnergyLimit=energy;}
62 void SetHighAngleLimit(G4double angle) {fHighAngleLimit=angle;}
63 void SetPPKineticEnergyCut(G4double kineticEnergyCut) {fPPKineticEnergyCut=kineticEnergyCut;}
64
65 /// set the number of pairs in sampling of Baier-Katkov Integral
66 /// (MC integration by e+- energy and angles <=> e+- momentum)
67 void SetSamplingPairsNumber(G4int nPairs){fNMCPairs = nPairs;}
68
69 /// set the number of particle angles 1/gamma in pair production
70 /// defining the width of the angular distribution of pair sampling
71 /// in the Baier-Katkov Integral
72 void SetChargeParticleAngleFactor(G4double chargeParticleAngleFactor)
73 {fChargeParticleAngleFactor = chargeParticleAngleFactor;}
74
75 /// set number of trajectory steps of a single particle (e- or e+)
76 void SetNTrajectorySteps(G4int nTrajectorySteps)
77 {fNTrajectorySteps = nTrajectorySteps;}
78
79private:
80
81 ///flag of simulation of incoherent scattering
82 G4bool fIncoherentScattering = false;
83
84 ///name of G4ChannelingFastSimModel from which an auto input should be performed
85 G4String fNameChannelingModel = "ChannelingModel";
86
87 ///name of G4Region to where the G4ChannelingFastSimModel is active
88 G4String fNameRegion = "Crystal";
89
90 G4double fLowEnergyLimit = 1*CLHEP::GeV;
91 G4double fHighAngleLimit = 50*CLHEP::mrad;
92
93 ///minimal kinetic energy of a charged particle produced
94 G4double fPPKineticEnergyCut = 1*CLHEP::MeV;
95
96 ///Monte Carlo statistics of e+- pair sampling in Baier-Katkov for 1 photon
97 G4int fNMCPairs = 150;
98
99 G4double fChargeParticleAngleFactor = 4; // number of particle angles 1/gamma:
100 // more fChargeParticleAngleFactor => higher paramParticleAngle
101
102 ///number of trajectory steps of a single particle (e- or e+)
103 G4int fNTrajectorySteps=250;
104
105};
106
107#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4CoherentPairProductionPhysics(const G4String &name="Coherent Pair Production Physics")
void SetNameG4Region(const G4String &nameG4Region)
set name of G4Region to where the G4ChannelingFastSimModel is active
void SetChargeParticleAngleFactor(G4double chargeParticleAngleFactor)
void SetPPKineticEnergyCut(G4double kineticEnergyCut)
void SetNameChannelingModel(const G4String &nameChannelingModel)
set functions
void SetNTrajectorySteps(G4int nTrajectorySteps)
set number of trajectory steps of a single particle (e- or e+)
G4VPhysicsConstructor(const G4String &="")