![]() |
|
|||
File indexing completed on 2025-09-18 09:14:06
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // Author: Alexei Sytov 0027 // Co-author: Gianfranco Paterno (modifications & testing) 0028 // On the base of the CRYSTALRAD realization of the Baier-Katkov integral: 0029 // A. I. Sytov, V. V. Tikhomirov, and L. Bandiera PRAB 22, 064601 (2019) 0030 0031 #ifndef G4BaierKatkov_h 0032 #define G4BaierKatkov_h 1 0033 0034 #include "globals.hh" 0035 #include <CLHEP/Units/SystemOfUnits.h> 0036 #include <vector> 0037 #include "G4ThreeVector.hh" 0038 0039 #include "G4VFastSimulationModel.hh" 0040 0041 /** \file G4BaierKatkov.hh 0042 * \brief Definition of the G4BaierKatkov class 0043 * This class is designed for the calculation of radiation probability, radiation point 0044 * and the parameters of the photon produced as well as spectrum accumulation using 0045 * the Baier-Katkov integral: 0046 * V. N. Baier, V. M. Katkov, and V. M. Strakhovenko, 0047 * Electromagnetic Processes at High Energies in Oriented Single Crystals 0048 * (World Scientific, Singapore, 1998). 0049 */ 0050 0051 class G4BaierKatkov 0052 { 0053 public: 0054 // default constructor 0055 G4BaierKatkov(); 0056 0057 // destructor 0058 ~G4BaierKatkov() = default; 0059 0060 /** 0061 You may call DoRadiation at each step of your trajectory 0062 CAUTION: please ensure that your steps are physically small enough for calculation 0063 of the radiation type you are interested in 0064 CAUTION: do ResetRadIntegral() before the start of a new trajectory 0065 0066 1) change some model defaults if necessary 0067 (SetSinglePhotonRadiationProbabilityLimit, 0068 SetNSmallTrajectorySteps, SetSpectrumEnergyRange) 0069 2) call DoRadiation at each step of your trajectory 0070 3) if DoRadiation returns TRUE, this means that a photon is produced (not added 0071 as a secondary yet) and its parameters are calculated. 0072 4) You may generate a new photon using GeneratePhoton either with 0073 the parameters calculated in DoRadiation or your own parameters. 0074 CAUTION: By now GeneratePhoton works only for a FastSim model 0075 5) Use GetPhotonEnergyInSpectrum() and GetTotalSpectrum() to return calculated 0076 total spectrum (all the photons altogether) 0077 Caution: is not normalized on the event number 0078 6) Get the charged particle parameters in the radiation point: 0079 GetParticleNewTotalEnergy(), 0080 GetParticleNewAngleX(), GetParticleNewAngleY(), 0081 GetNewGlobalTime(), 0082 GetParticleNewCoordinateXYZ() 0083 */ 0084 0085 ///get functions 0086 0087 /// get maximal radiation probability to preserve single photon radiation 0088 G4double GetSinglePhotonRadiationProbabilityLimit() 0089 {return fSinglePhotonRadiationProbabilityLimit;} 0090 0091 ///CAUTION! : use the get functions below ONLY AFTER the call of DoRadiation 0092 /// and ONLY IF IT RETURNS true 0093 0094 ///total probability of radiation: needs calculation of DoRadiation first 0095 G4double GetTotalRadiationProbability(){return fTotalRadiationProbability;} 0096 ///get new parameters of the particle 0097 ///(the parameters at the point of radiation emission) 0098 ///needs calculation of DoRadiation first 0099 G4double GetParticleNewTotalEnergy(){return fNewParticleEnergy;} 0100 G4double GetParticleNewAngleX(){return fNewParticleAngleX;} 0101 G4double GetParticleNewAngleY(){return fNewParticleAngleY;} 0102 G4double GetNewGlobalTime(){return fNewGlobalTime;} 0103 const G4ThreeVector& GetParticleNewCoordinateXYZ(){return fNewParticleCoordinateXYZ;} 0104 0105 ///get photon energies (x-value in spectrum) 0106 const std::vector<G4double>& GetPhotonEnergyInSpectrum() 0107 {return fPhotonEnergyInSpectrum;} 0108 0109 ///get fTotalSpectrum after finishing the trajectory part with DoRadiation 0110 const std::vector<G4double>& GetTotalSpectrum(){return fTotalSpectrum;} 0111 0112 ///set functions 0113 0114 ///set maximal radiation probability to preserve single photon radiation 0115 void SetSinglePhotonRadiationProbabilityLimit(G4double wmax) 0116 {fSinglePhotonRadiationProbabilityLimit = wmax;} 0117 0118 ///number of steps in a trajectory small piece before 0119 ///the next call of the radiation integral 0120 void SetNSmallTrajectorySteps(G4int nSmallTrajectorySteps) 0121 {fNSmallTrajectorySteps = nSmallTrajectorySteps;} 0122 0123 ///reinitialize intermediate integrals fFa, fSs, fSc, fSsx, fSsy, fScx, fScy; 0124 ///reset radiation integral internal variables to defaults; 0125 ///reset the trajectory and radiation probability along the trajectory 0126 void ResetRadIntegral(); 0127 0128 ///setting the number of photons in sampling of Baier-Katkov Integral 0129 ///(MC integration by photon energy and angles <=> photon momentum) 0130 void SetSamplingPhotonsNumber(G4int nPhotons){fNMCPhotons = nPhotons;} 0131 0132 ///setting the number of radiation angles 1/gamma, defining the width of 0133 ///the angular distribution of photon sampling in the Baier-Katkov Integral 0134 void SetRadiationAngleFactor(G4double radiationAngleFactor) 0135 {fRadiationAngleFactor = radiationAngleFactor;} 0136 0137 ///CAUTION, the bins width is logarithmic 0138 ///Do not worry if the maximal energy > particle energy. 0139 ///This elements of spectrum with non-physical energies 0140 ///will not be processed (they will be 0). 0141 void SetSpectrumEnergyRange(G4double emin, 0142 G4double emax, 0143 G4int numberOfBins); 0144 /// SetSpectrumEnergyRange also calls ResetRadIntegral() 0145 0146 void SetMinPhotonEnergy(G4double emin){SetSpectrumEnergyRange(emin, 0147 fMaxPhotonEnergy, 0148 fNBinsSpectrum);} 0149 void SetMaxPhotonEnergy(G4double emax){SetSpectrumEnergyRange(fMinPhotonEnergy, 0150 emax, 0151 fNBinsSpectrum);} 0152 void SetNBinsSpectrum(G4int nbin){SetSpectrumEnergyRange(fMinPhotonEnergy, 0153 fMaxPhotonEnergy, 0154 nbin);} 0155 0156 /// Increase the statistic of virtual photons in a certain energy region 0157 /// CAUTION! : don't do it before SetSpectrumEnergyRange or SetMinPhotonEnergy 0158 void AddStatisticsInPhotonEnergyRegion(G4double emin, G4double emax, 0159 G4int timesPhotonStatistics); 0160 0161 /// Virtual collimator masks the selection of photon angles in fTotalSpectrum 0162 /// Virtual collimator doesn't influence on Geant4 simulations. 0163 void SetVirtualCollimator(G4double virtualCollimatorAngularDiameter) 0164 {fVirtualCollimatorAngularDiameter=virtualCollimatorAngularDiameter;} 0165 0166 /// add the new elements of the trajectory, calculate radiation in a crystal 0167 /// see complete description in G4BaierKatkov::DoRadiation 0168 /// calls RadIntegral and all the necessary functions 0169 /// sets the parameters of a photon produced (if any) 0170 /// using SetPhotonProductionParameters() 0171 /// returns true in the case of photon generation, false if not 0172 G4bool DoRadiation(G4double etotal, G4double mass, 0173 G4double angleX, G4double angleY, 0174 G4double angleScatteringX, G4double angleScatteringY, 0175 G4double step, G4double globalTime, 0176 G4ThreeVector coordinateXYZ, 0177 G4bool flagEndTrajectory=false); 0178 0179 /// generates secondary photon belonging to fastStep with variables 0180 /// photon energy, momentum direction, coordinates and global time 0181 /// CALCULATED IN DoRadiation => USE IT ONLY AFTER DoRadiation returns true 0182 void GeneratePhoton(G4FastStep &fastStep); 0183 0184 private: 0185 0186 ///set functions 0187 0188 ///function setting the photon sampling parameters in the Baier-Katkov integral; 0189 ///only the maximal energy is set, while fMinPhotonEnergy is used as a minimal energy; 0190 ///the angles set the angular distribution (the tails are infinite) 0191 void SetPhotonSamplingParameters(G4double ekin, 0192 G4double minPhotonAngleX, G4double maxPhotonAngleX, 0193 G4double minPhotonAngleY, G4double maxPhotonAngleY); 0194 0195 ///main functions: 0196 0197 ///generation of the photons in sampling of Baier-Katkov Integral 0198 ///(MC integration by photon energy and angles <=> by photon momentum) 0199 void GeneratePhotonSampling(); 0200 0201 ///Baier-Katkov method: calculation of integral, spectrum, full probability; 0202 ///returns the total radiation probability; 0203 ///calculates the radiation spectrum on this trajectory piece 0204 G4double RadIntegral(G4double etotal, G4double mass, 0205 std::vector<G4double> &vectorParticleAnglesX, 0206 std::vector<G4double> &vectorParticleAnglesY, 0207 std::vector<G4double> &vectorScatteringAnglesX, 0208 std::vector<G4double> &vectorScatteringAnglesY, 0209 std::vector<G4double> &vectorSteps, 0210 G4int imin); 0211 0212 ///set photon production parameters (returns false if no photon produced) 0213 ///accumulates fTotalSpectrum 0214 ///CAUTION: it is an accessory function of DoRadiation, do not use it separately 0215 G4bool SetPhotonProductionParameters(G4double etotal, G4double mass); 0216 0217 G4int FindVectorIndex(std::vector<G4double> &myvector, G4double value); 0218 0219 G4double fTotalRadiationProbability = 0.; 0220 G4double fSinglePhotonRadiationProbabilityLimit=0.25;//Maximal radiation 0221 //probability to preserve single photon radiation 0222 0223 //number of steps in a trajectory piece before the next call of the radiation integral 0224 G4int fNSmallTrajectorySteps=10000; 0225 ///trajectory element No (the first element of the array feeded in RadIntegral) 0226 G4int fImin0 = 0; 0227 ///Monte Carlo statistics of photon sampling in Baier-Katkov with 1 trajectory 0228 G4int fNMCPhotons =150; 0229 ///the number of bins in photon spectrum 0230 G4int fNBinsSpectrum = 110; 0231 G4double fMinPhotonEnergy = 0.1*CLHEP::MeV;//min energy in spectrum output 0232 G4double fMaxPhotonEnergy = 1*CLHEP::GeV; //max energy in spectrum output 0233 G4double fLogEmaxdEmin = 1.;// = log(fMaxPhotonEnergy/fMinPhotonEnergy), 0234 // 1/normalizing coefficient in 0235 // 1/E distribution between 0236 // fMinPhotonEnergy and fMaxPhotonEnergy 0237 // is used only for spectrum output, not for simulations 0238 //(we take bremsstrahlung for photon sampling) 0239 G4double fLogEdEmin = 1.; // = log(E/fMinPhotonEnergy), the same as fLogEmaxdEmin 0240 // but with the particle energy as the maximal limit 0241 0242 G4double fVirtualCollimatorAngularDiameter=1.;//default, infinite angle 0243 std::vector<G4bool> fInsideVirtualCollimator; 0244 0245 ///data of the phootn energy range with additional statistics 0246 std::vector<G4double> fLogAddRangeEmindEmin;//=G4Log(emin/fMinPhotonEnergy) 0247 std::vector<G4double> fLogAddRangeEmaxdEmin;//=G4Log(emax/fMinPhotonEnergy) 0248 std::vector<G4int> fTimesPhotonStatistics; 0249 0250 ///number of trajectories 0251 //(at each of the Baier-Katkov Integral is calculated for the same photons) 0252 G4int fItrajectories = 0; 0253 0254 G4double fEph0=0; //energy of the photon produced 0255 G4ThreeVector PhMomentumDirection; //momentum direction of the photon produced 0256 0257 ///Radiation integral variables 0258 G4double fMeanPhotonAngleX =0.; //average angle of radiated photon direction 0259 //in sampling, x-plane 0260 G4double fParamPhotonAngleX=1.e-3*CLHEP::rad; //a parameter radiated photon 0261 //sampling distribution, x-plane 0262 G4double fMeanPhotonAngleY =0.; //average angle of radiated photon direction 0263 //in sampling, y-plane 0264 G4double fParamPhotonAngleY=1.e-3*CLHEP::rad; //a parameter radiated photon 0265 //sampling distribution, y-plane 0266 G4double fRadiationAngleFactor = 4.; // number of radiation angles 1/gamma: 0267 // more fRadiationAngleFactor => 0268 // higher fParamPhotonAngleX and Y 0269 0270 ///new particle parameters (the parameters at the point of radiation emission) 0271 G4double fNewParticleEnergy=0; 0272 G4double fNewParticleAngleX=0; 0273 G4double fNewParticleAngleY=0; 0274 G4double fNewGlobalTime=0; 0275 G4ThreeVector fNewParticleCoordinateXYZ; 0276 0277 ///sampling of the energy and the angles of a photon emission 0278 ///(integration variables, Monte Carlo integration) 0279 std::vector<G4double> fPhotonEnergyInIntegral; 0280 std::vector<G4double> fPhotonAngleInIntegralX; 0281 std::vector<G4double> fPhotonAngleInIntegralY; 0282 std::vector<G4double> fPhotonAngleNormCoef; 0283 ///spectrum bin index for each photon 0284 std::vector<G4double> fIBinsSpectrum; 0285 ///the vector of the discrete CDF of the radiation of sampling photons 0286 std::vector<G4double> fPhotonProductionCDF; 0287 0288 ///vectors of the trajectory 0289 std::vector<G4double> fParticleAnglesX; 0290 std::vector<G4double> fParticleAnglesY; 0291 std::vector<G4double> fScatteringAnglesX; 0292 std::vector<G4double> fScatteringAnglesY; 0293 std::vector<G4double> fSteps; 0294 std::vector<G4double> fGlobalTimes; 0295 std::vector<G4ThreeVector> fParticleCoordinatesXYZ; 0296 0297 ///intermediate integrals (different for each photon energy value)!!! 0298 std::vector<G4double> fFa;//phase 0299 std::vector<G4double> fSs; 0300 std::vector<G4double> fSc; 0301 std::vector<G4double> fSsx; 0302 std::vector<G4double> fSsy; 0303 std::vector<G4double> fScx; 0304 std::vector<G4double> fScy; 0305 0306 ///output 0307 std::vector<G4double> fPhotonEnergyInSpectrum; //energy values in spectrum 0308 0309 std::vector<G4int> fNPhotonsPerBin; //number of photons per spectrum bin 0310 //(accumulating during total run) 0311 0312 std::vector<G4double> fSpectrum; //spectrum normalized by the total 0313 //radiation probability of one particle at one call of RadIntegral 0314 0315 std::vector<std::vector<G4double>> fAccumSpectrum; //accumulate Spectrum during 0316 //the part of a trajectory 0317 0318 std::vector<G4double> fAccumTotalSpectrum; //spectrum normalized by the total 0319 //radiation probability summed 0320 //for all the particles (is not divided 0321 //of one particle number fNPhotonsPerBin) 0322 0323 std::vector<G4double> fTotalSpectrum; //spectrum normalized by 0324 //the total radiation probability summed 0325 //for all the particles 0326 //(is divided by the photon number fNPhotonsPerBin) 0327 //multiplied by the number of trajectories 0328 //(fItrajectories) 0329 0330 std::vector<G4double> fImax0; //trajectory element numbers at the end of each 0331 //small piece; G4double just for security of some operations 0332 ///total radiation probability along this trajectory 0333 std::vector<G4double> fTotalRadiationProbabilityAlongTrajectory; 0334 }; 0335 0336 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |