Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BaierKatkov.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#ifndef G4BaierKatkov_h
28#define G4BaierKatkov_h 1
29
30#include "globals.hh"
32#include <vector>
33#include "G4ThreeVector.hh"
34
36
37/** \file G4BaierKatkov.hh
38* \brief Definition of the G4BaierKatkov class
39* This class is designed for the calculation of radiation probability, radiation point
40* and the parameters of the photon produced as well as spectrum accumulation using
41* the Baier-Katkov integral:
42* V. N. Baier, V. M. Katkov, and V. M. Strakhovenko,
43* Electromagnetic Processes at High Energies in Oriented Single Crystals
44* (World Scientific, Singapore, 1998).
45*/
46
48{
49public:
50 // default constructor
52
53 // destructor
54 ~G4BaierKatkov() = default;
55
56 /**
57 You may call DoRadiation at each step of your trajectory
58 CAUTION: please ensure that your steps are physically small enough for calculation
59 of the radiation type you are interested in
60 CAUTION: do ResetRadIntegral() before the start of a new trajectory
61
62 1) change some model defaults if necessary
63 (SetSinglePhotonRadiationProbabilityLimit,
64 SetNSmallTrajectorySteps, SetSpectrumEnergyRange)
65 2) call DoRadiation at each step of your trajectory
66 3) if DoRadiation returns TRUE, this means that a photon is produced (not added
67 as a secondary yet) and its parameters are calculated.
68 4) You may generate a new photon using GeneratePhoton either with
69 the parameters calculated in DoRadiation or your own parameters.
70 CAUTION: By now GeneratePhoton works only for a FastSim model
71 5) Use GetPhotonEnergyInSpectrum() and GetTotalSpectrum() to return calculated
72 total spectrum (all the photons altogether)
73 Caution: is not normalized on the event number
74 6) Get the charged particle parameters in the radiation point:
75 GetParticleNewTotalEnergy(),
76 GetParticleNewAngleX(), GetParticleNewAngleY(),
77 GetNewGlobalTime(),
78 GetParticleNewCoordinateXYZ()
79 */
80
81 ///get functions
82
83 /// get maximal radiation probability to preserve single photon radiation
85 {return fSinglePhotonRadiationProbabilityLimit;}
86
87 ///CAUTION! : use the get functions below ONLY AFTER the call of DoRadiation
88 /// and ONLY IF IT RETURNS true
89
90 ///total probability of radiation: needs calculation of DoRadiation first
91 G4double GetTotalRadiationProbability(){return fTotalRadiationProbability;}
92 ///get new parameters of the particle
93 ///(the parameters at the point of radiation emission)
94 ///needs calculation of DoRadiation first
95 G4double GetParticleNewTotalEnergy(){return fNewParticleEnergy;}
96 G4double GetParticleNewAngleX(){return fNewParticleAngleX;}
97 G4double GetParticleNewAngleY(){return fNewParticleAngleY;}
98 G4double GetNewGlobalTime(){return fNewGlobalTime;}
99 const G4ThreeVector& GetParticleNewCoordinateXYZ(){return fNewParticleCoordinateXYZ;}
100
101 ///get photon energies (x-value in spectrum)
102 const std::vector<G4double>& GetPhotonEnergyInSpectrum()
103 {return fPhotonEnergyInSpectrum;}
104
105 ///get fTotalSpectrum after finishing the trajectory part with DoRadiation
106 const std::vector<G4double>& GetTotalSpectrum(){return fTotalSpectrum;}
107
108 ///set functions
109
110 ///set maximal radiation probability to preserve single photon radiation
112 {fSinglePhotonRadiationProbabilityLimit = wmax;}
113
114 ///number of steps in a trajectory small piece before
115 ///the next call of the radiation integral
116 void SetNSmallTrajectorySteps(G4double nSmallTrajectorySteps)
117 {fNSmallTrajectorySteps = nSmallTrajectorySteps;}
118
119 ///reinitialize intermediate integrals fFa, fSs, fSc, fSsx, fSsy, fScx, fScy;
120 ///reset radiation integral internal variables to defaults;
121 ///reset the trajectory and radiation probability along the trajectory
122 void ResetRadIntegral();
123
124 ///setting the number of photons in sampling of Baier-Katkov Integral
125 ///(MC integration by photon energy and angles <=> photon momentum)
126 void SetSamplingPhotonsNumber(G4int nPhotons){fNMCPhotons = nPhotons;}
127
128 ///setting the number of radiation angles 1/gamma, defining the width of
129 ///the angular distribution of photon sampling in the Baier-Katkov Integral
130 void SetRadiationAngleFactor(G4double radiationAngleFactor)
131 {fRadiationAngleFactor = radiationAngleFactor;}
132
133 ///CAUTION, the bins width is logarithmic
134 ///Do not worry if the maximal energy > particle energy.
135 ///This elements of spectrum with non-physical energies
136 ///will not be processed (they will be 0).
138 G4double emax,
139 G4int numberOfBins);
140 /// SetSpectrumEnergyRange also calls ResetRadIntegral()
141
143 fMaxPhotonEnergy,
144 fNBinsSpectrum);}
146 emax,
147 fNBinsSpectrum);}
148 void SetNBinsSpectrum(G4int nbin){SetSpectrumEnergyRange(fMinPhotonEnergy,
149 fMaxPhotonEnergy,
150 nbin);}
151
152 /// Increase the statistic of virtual photons in a certain energy region
153 /// CAUTION! : don't do it before SetSpectrumEnergyRange or SetMinPhotonEnergy
155 G4int timesPhotonStatistics);
156
157 /// Virtual collimator masks the selection of photon angles in fTotalSpectrum
158 /// Virtual collimator doesn't influence on Geant4 simulations.
159 void SetVirtualCollimator(G4double virtualCollimatorAngularDiameter)
160 {fVirtualCollimatorAngularDiameter=virtualCollimatorAngularDiameter;}
161
162 /// add the new elements of the trajectory, calculate radiation in a crystal
163 /// see complete description in G4BaierKatkov::DoRadiation
164 /// calls RadIntegral and all the necessary functions
165 /// sets the parameters of a photon produced (if any)
166 /// using SetPhotonProductionParameters()
167 /// returns true in the case of photon generation, false if not
168 G4bool DoRadiation(G4double etotal, G4double mass,
169 G4double angleX, G4double angleY,
170 G4double angleScatteringX, G4double angleScatteringY,
171 G4double step, G4double globalTime,
172 G4ThreeVector coordinateXYZ,
173 G4bool flagEndTrajectory=false);
174
175 /// generates secondary photon belonging to fastStep with variables
176 /// photon energy, momentum direction, coordinates and global time
177 /// CALCULATED IN DoRadiation => USE IT ONLY AFTER DoRadiation returns true
178 void GeneratePhoton(G4FastStep &fastStep);
179
180private:
181
182 ///set functions
183
184 ///function setting the photon sampling parameters in the Baier-Katkov integral;
185 ///only the maximal energy is set, while fMinPhotonEnergy is used as a minimal energy;
186 ///the angles set the angular distribution (the tails are infinite)
187 void SetPhotonSamplingParameters(G4double ekin,
188 G4double minPhotonAngleX, G4double maxPhotonAngleX,
189 G4double minPhotonAngleY, G4double maxPhotonAngleY);
190
191 ///main functions:
192
193 ///generation of the photons in sampling of Baier-Katkov Integral
194 ///(MC integration by photon energy and angles <=> by photon momentum)
195 void GeneratePhotonSampling();
196
197 ///Baier-Katkov method: calculation of integral, spectrum, full probability;
198 ///returns the total radiation probability;
199 ///calculates the radiation spectrum on this trajectory piece
200 G4double RadIntegral(G4double etotal, G4double mass,
201 std::vector<G4double> &vectorParticleAnglesX,
202 std::vector<G4double> &vectorParticleAnglesY,
203 std::vector<G4double> &vectorScatteringAnglesX,
204 std::vector<G4double> &vectorScatteringAnglesY,
205 std::vector<G4double> &vectorSteps,
206 G4int imin);
207
208 ///set photon production parameters (returns false if no photon produced)
209 ///accumulates fTotalSpectrum
210 ///CAUTION: it is an accessory function of DoRadiation, do not use it separately
211 G4bool SetPhotonProductionParameters(G4double etotal, G4double mass);
212
213 G4int FindVectorIndex(std::vector<G4double> &myvector, G4double value);
214
215 G4double fTotalRadiationProbability = 0.;
216 G4double fSinglePhotonRadiationProbabilityLimit=0.25;//Maximal radiation
217 //probability to preserve single photon radiation
218
219 //number of steps in a trajectory piece before the next call of the radiation integral
220 G4int fNSmallTrajectorySteps=10000;
221 ///trajectory element No (the first element of the array feeded in RadIntegral)
222 G4int fImin0 = 0;
223 ///Monte Carlo statistics of photon sampling in Baier-Katkov with 1 trajectory
224 G4int fNMCPhotons =150;
225 ///the number of bins in photon spectrum
226 G4int fNBinsSpectrum = 110;
227 G4double fMinPhotonEnergy = 0.1*CLHEP::MeV;//min energy in spectrum output
228 G4double fMaxPhotonEnergy = 1*CLHEP::GeV; //max energy in spectrum output
229 G4double fLogEmaxdEmin = 1.;// = log(fMaxPhotonEnergy/fMinPhotonEnergy),
230 // 1/normalizing coefficient in
231 // 1/E distribution between
232 // fMinPhotonEnergy and fMaxPhotonEnergy
233 // is used only for spectrum output, not for simulations
234 //(we take bremsstrahlung for photon sampling)
235 G4double fLogEdEmin = 1.; // = log(E/fMinPhotonEnergy), the same as fLogEmaxdEmin
236 // but with the particle energy as the maximal limit
237
238 G4double fVirtualCollimatorAngularDiameter=1.;//default, infinite angle
239 std::vector<G4bool> fInsideVirtualCollimator;
240
241 ///data of the phootn energy range with additional statistics
242 std::vector<G4double> fLogAddRangeEmindEmin;//=G4Log(emin/fMinPhotonEnergy)
243 std::vector<G4double> fLogAddRangeEmaxdEmin;//=G4Log(emax/fMinPhotonEnergy)
244 std::vector<G4int> fTimesPhotonStatistics;
245
246 ///number of trajectories
247 //(at each of the Baier-Katkov Integral is calculated for the same photons)
248 G4int fItrajectories = 0;
249
250 G4double fEph0=0; //energy of the photon produced
251 G4ThreeVector PhMomentumDirection; //momentum direction of the photon produced
252
253 ///Radiation integral variables
254 G4double fMeanPhotonAngleX =0.; //average angle of radiated photon direction
255 //in sampling, x-plane
256 G4double fParamPhotonAngleX=1.e-3*CLHEP::rad; //a parameter radiated photon
257 //sampling distribution, x-plane
258 G4double fMeanPhotonAngleY =0.; //average angle of radiated photon direction
259 //in sampling, y-plane
260 G4double fParamPhotonAngleY=1.e-3*CLHEP::rad; //a parameter radiated photon
261 //sampling distribution, y-plane
262 G4double fRadiationAngleFactor = 1.; // number of radiation angles 1/gamma:
263 // more fRadiationAngleFactor =>
264 // higher fParamPhotonAngleX and Y
265
266 ///new particle parameters (the parameters at the point of radiation emission)
267 G4double fNewParticleEnergy=0;
268 G4double fNewParticleAngleX=0;
269 G4double fNewParticleAngleY=0;
270 G4double fNewGlobalTime=0;
271 G4ThreeVector fNewParticleCoordinateXYZ;
272
273 ///sampling of the energy and the angles of a photon emission
274 ///(integration variables, Monte Carlo integration)
275 std::vector<G4double> fPhotonEnergyInIntegral;
276 std::vector<G4double> fPhotonAngleInIntegralX;
277 std::vector<G4double> fPhotonAngleInIntegralY;
278 std::vector<G4double> fPhotonAngleNormCoef;
279 ///spectrum bin index for each photon
280 std::vector<G4double> fIBinsSpectrum;
281 ///the vector of the discrete CDF of the radiation of sampling photons
282 std::vector<G4double> fPhotonProductionCDF;
283
284 ///vectors of the trajectory
285 std::vector<G4double> fParticleAnglesX;
286 std::vector<G4double> fParticleAnglesY;
287 std::vector<G4double> fScatteringAnglesX;
288 std::vector<G4double> fScatteringAnglesY;
289 std::vector<G4double> fSteps;
290 std::vector<G4double> fGlobalTimes;
291 std::vector<G4ThreeVector> fParticleCoordinatesXYZ;
292
293 ///intermediate integrals (different for each photon energy value)!!!
294 std::vector<G4double> fFa;//phase
295 std::vector<G4double> fSs;
296 std::vector<G4double> fSc;
297 std::vector<G4double> fSsx;
298 std::vector<G4double> fSsy;
299 std::vector<G4double> fScx;
300 std::vector<G4double> fScy;
301
302 ///output
303 std::vector<G4double> fPhotonEnergyInSpectrum; //energy values in spectrum
304
305 std::vector<G4int> fNPhotonsPerBin; //number of photons per spectrum bin
306 //(accumulating during total run)
307
308 std::vector<G4double> fSpectrum; //spectrum normalized by the total
309 //radiation probability of one particle at one call of RadIntegral
310
311 std::vector<std::vector<G4double>> fAccumSpectrum; //accumulate Spectrum during
312 //the part of a trajectory
313
314 std::vector<G4double> fAccumTotalSpectrum; //spectrum normalized by the total
315 //radiation probability summed
316 //for all the particles (is not divided
317 //of one particle number fNPhotonsPerBin)
318
319 std::vector<G4double> fTotalSpectrum; //spectrum normalized by
320 //the total radiation probability summed
321 //for all the particles
322 //(is divided by the photon number fNPhotonsPerBin)
323 //multiplied by the number of trajectories
324 //(fItrajectories)
325
326 std::vector<G4double> fImax0; //trajectory element numbers at the end of each
327 //small piece; G4double just for security of some operations
328 ///total radiation probability along this trajectory
329 std::vector<G4double> fTotalRadiationProbabilityAlongTrajectory;
330};
331
332#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetRadiationAngleFactor(G4double radiationAngleFactor)
void ResetRadIntegral()
const std::vector< G4double > & GetPhotonEnergyInSpectrum()
get photon energies (x-value in spectrum)
void GeneratePhoton(G4FastStep &fastStep)
void SetSpectrumEnergyRange(G4double emin, G4double emax, G4int numberOfBins)
void AddStatisticsInPhotonEnergyRegion(G4double emin, G4double emax, G4int timesPhotonStatistics)
G4double GetParticleNewAngleX()
G4double GetSinglePhotonRadiationProbabilityLimit()
get functions
G4double GetParticleNewTotalEnergy()
void SetSamplingPhotonsNumber(G4int nPhotons)
void SetVirtualCollimator(G4double virtualCollimatorAngularDiameter)
void SetSinglePhotonRadiationProbabilityLimit(G4double wmax)
set functions
void SetMaxPhotonEnergy(G4double emax)
~G4BaierKatkov()=default
void SetMinPhotonEnergy(G4double emin)
SetSpectrumEnergyRange also calls ResetRadIntegral()
void SetNBinsSpectrum(G4int nbin)
G4bool DoRadiation(G4double etotal, G4double mass, G4double angleX, G4double angleY, G4double angleScatteringX, G4double angleScatteringY, G4double step, G4double globalTime, G4ThreeVector coordinateXYZ, G4bool flagEndTrajectory=false)
G4double GetTotalRadiationProbability()
total probability of radiation: needs calculation of DoRadiation first
G4double GetParticleNewAngleY()
const std::vector< G4double > & GetTotalSpectrum()
get fTotalSpectrum after finishing the trajectory part with DoRadiation
G4double GetNewGlobalTime()
const G4ThreeVector & GetParticleNewCoordinateXYZ()
void SetNSmallTrajectorySteps(G4double nSmallTrajectorySteps)