Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BaierKatkov Class Reference

#include <G4BaierKatkov.hh>

Public Member Functions

 G4BaierKatkov ()
 
 ~G4BaierKatkov ()=default
 
G4double GetSinglePhotonRadiationProbabilityLimit ()
 get functions
 
G4double GetTotalRadiationProbability ()
 total probability of radiation: needs calculation of DoRadiation first
 
G4double GetParticleNewTotalEnergy ()
 
G4double GetParticleNewAngleX ()
 
G4double GetParticleNewAngleY ()
 
G4double GetNewGlobalTime ()
 
const G4ThreeVectorGetParticleNewCoordinateXYZ ()
 
const std::vector< G4double > & GetPhotonEnergyInSpectrum ()
 get photon energies (x-value in spectrum)
 
const std::vector< G4double > & GetTotalSpectrum ()
 get fTotalSpectrum after finishing the trajectory part with DoRadiation
 
void SetSinglePhotonRadiationProbabilityLimit (G4double wmax)
 set functions
 
void SetNSmallTrajectorySteps (G4int nSmallTrajectorySteps)
 
void ResetRadIntegral ()
 
void SetSamplingPhotonsNumber (G4int nPhotons)
 
void SetRadiationAngleFactor (G4double radiationAngleFactor)
 
void SetSpectrumEnergyRange (G4double emin, G4double emax, G4int numberOfBins)
 
void SetMinPhotonEnergy (G4double emin)
 SetSpectrumEnergyRange also calls ResetRadIntegral()
 
void SetMaxPhotonEnergy (G4double emax)
 
void SetNBinsSpectrum (G4int nbin)
 
void AddStatisticsInPhotonEnergyRegion (G4double emin, G4double emax, G4int timesPhotonStatistics)
 
void SetVirtualCollimator (G4double virtualCollimatorAngularDiameter)
 
G4bool DoRadiation (G4double etotal, G4double mass, G4double angleX, G4double angleY, G4double angleScatteringX, G4double angleScatteringY, G4double step, G4double globalTime, G4ThreeVector coordinateXYZ, G4bool flagEndTrajectory=false)
 
void GeneratePhoton (G4FastStep &fastStep)
 

Detailed Description

Definition at line 51 of file G4BaierKatkov.hh.

Constructor & Destructor Documentation

◆ G4BaierKatkov()

G4BaierKatkov::G4BaierKatkov ( )

Definition at line 39 of file G4BaierKatkov.cc.

40{
41 //sets the default spectrum energy range of integration and
42 //calls ResetRadIntegral()
43 SetSpectrumEnergyRange(0.1*MeV,1.*GeV,110);
44
45 //Do not worry if the maximal energy > particle energy
46 //this elements of spectrum with non-physical energies
47 //will not be processed (they will be 0)
48
49 G4cout << " "<< G4endl;
50 G4cout << "G4BaierKatkov model is activated."<< G4endl;
51 G4cout << " "<< G4endl;
52}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetSpectrumEnergyRange(G4double emin, G4double emax, G4int numberOfBins)

◆ ~G4BaierKatkov()

G4BaierKatkov::~G4BaierKatkov ( )
default

Member Function Documentation

◆ AddStatisticsInPhotonEnergyRegion()

void G4BaierKatkov::AddStatisticsInPhotonEnergyRegion ( G4double emin,
G4double emax,
G4int timesPhotonStatistics )

Increase the statistic of virtual photons in a certain energy region CAUTION! : don't do it before SetSpectrumEnergyRange or SetMinPhotonEnergy

Definition at line 152 of file G4BaierKatkov.cc.

155{
156
157 if(timesPhotonStatistics<=1)
158 {
159 G4cout << "G4BaierKatkov model, "
160 "function AddStatisticsInPhotonEnergyRegion("
161 << emin/CLHEP::MeV << " MeV, "
162 << emax/CLHEP::MeV << " MeV, "
163 << timesPhotonStatistics << ")" << G4endl;
164 G4cout << "Warning: the statistics factor cannot be <=1." << G4endl;
165 G4cout << "The statistics was not added." << G4endl;
166 G4cout << " "<< G4endl;
167 }
168 else if(fMinPhotonEnergy>emin)
169 {
170 G4cout << "G4BaierKatkov model, "
171 "function AddStatisticsInPhotonEnergyRegion("
172 << emin/CLHEP::MeV << " MeV, "
173 << emax/CLHEP::MeV << " MeV, "
174 << timesPhotonStatistics << ")" << G4endl;
175 G4cout << "Warning: the minimal energy inserted is less then "
176 "the minimal energy cut of the spectrum: "
177 << fMinPhotonEnergy/CLHEP::MeV << " MeV." << G4endl;
178 G4cout << "The statistics was not added." << G4endl;
179 G4cout << " "<< G4endl;
180 }
181 else if(emax-emin<DBL_EPSILON)
182 {
183 G4cout << "G4BaierKatkov model, "
184 "function AddStatisticsInPhotonEnergyRegion("
185 << emin/CLHEP::MeV << " MeV, "
186 << emax/CLHEP::MeV << " MeV, "
187 << timesPhotonStatistics << ")" << G4endl;
188 G4cout << "Warning: the maximal energy <= the minimal energy." << G4endl;
189 G4cout << "The statistics was not added." << G4endl;
190 G4cout << " "<< G4endl;
191 }
192 else
193 {
194 G4bool setrange = true;
195 G4double logAddRangeEmindEmin = G4Log(emin/fMinPhotonEnergy);
196 G4double logAddRangeEmaxdEmin = G4Log(emax/fMinPhotonEnergy);
197
198 G4int nAddRange = (G4int)fTimesPhotonStatistics.size();
199 for (G4int j=0;j<nAddRange;j++)
200 {
201 if((logAddRangeEmindEmin>=fLogAddRangeEmindEmin[j]&&
202 logAddRangeEmindEmin< fLogAddRangeEmaxdEmin[j])||
203 (logAddRangeEmaxdEmin> fLogAddRangeEmindEmin[j]&&
204 logAddRangeEmaxdEmin<=fLogAddRangeEmaxdEmin[j])||
205 (logAddRangeEmindEmin<=fLogAddRangeEmindEmin[j]&&
206 logAddRangeEmaxdEmin>=fLogAddRangeEmaxdEmin[j]))
207 {
208 G4cout << "G4BaierKatkov model, "
209 "function AddStatisticsInPhotonEnergyRegion("
210 << emin/CLHEP::MeV << " MeV, "
211 << emax/CLHEP::MeV << " MeV, "
212 << timesPhotonStatistics << ")" << G4endl;
213 G4cout << "Warning: the energy range intersects another "
214 "added energy range." << G4endl;
215 G4cout << "The statistics was not added." << G4endl;
216 G4cout << " "<< G4endl;
217 setrange = false;
218 break;
219 }
220 }
221 if (setrange)
222 {
223 fLogAddRangeEmindEmin.push_back(logAddRangeEmindEmin);
224 fLogAddRangeEmaxdEmin.push_back(logAddRangeEmaxdEmin);
225 fTimesPhotonStatistics.push_back(timesPhotonStatistics);
226
227 G4cout << "G4BaierKatkov model: increasing the statistics of photon sampling "
228 "in Baier-Katkov with a factor of "
229 << timesPhotonStatistics << G4endl;
230 G4cout << "in the energy spectrum range: ("
231 << emin/CLHEP::MeV << " MeV, "
232 << emax/CLHEP::MeV << " MeV)" << G4endl;
233 }
234 }
235}
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define DBL_EPSILON
Definition templates.hh:66

◆ DoRadiation()

G4bool G4BaierKatkov::DoRadiation ( G4double etotal,
G4double mass,
G4double angleX,
G4double angleY,
G4double angleScatteringX,
G4double angleScatteringY,
G4double step,
G4double globalTime,
G4ThreeVector coordinateXYZ,
G4bool flagEndTrajectory = false )

add the new elements of the trajectory, calculate radiation in a crystal see complete description in G4BaierKatkov::DoRadiation calls RadIntegral and all the necessary functions sets the parameters of a photon produced (if any) using SetPhotonProductionParameters() returns true in the case of photon generation, false if not

DoRadiation description

The real trajectory is divided onto parts. Every part is accumulated until the radiation cannot be considered as single photon, i.e. the radiation probability exceeds some threshold fSinglePhotonRadiationProbabilityLimit. Typically fSinglePhotonRadiationProbabilityLimit should be between 0.01 and 0.1. Then the photon generation as well as its parameters are simulated in SetPhotonProductionParameters()

Finally if radiation took place, this function sets the new particle parameters (the parameters at the point of radiation emission) fNewParticleEnergy; fNewParticleAngleX; fNewParticleAngleY; NewStep; fNewGlobalTime; fNewParticleCoordinateXYZ;

In order to control if the trajectory part reaches this limit, one needs to divide it into smaller pieces (consisting of fNSmallTrajectorySteps elements, typically several hundred) and to calculate the total radiation probability along the trajectory part after each small piece accomplished. If it exceeds the fSinglePhotonRadiationProbabilityLimit, the trajectory part is over and the radiation can be generated.

In order to calculate the radiation probability the Baier-Katkov integral is called after each small piece of the trajectory. It is summarized with the integral along the previous small pieces for this part of the trajectory.

To speed-up the simulations, the photon angles are generated only once at the start of the trajectory part.

Definition at line 652 of file G4BaierKatkov.cc.

658{
659 /**
660 DoRadiation description
661
662 The real trajectory is divided onto parts.
663 Every part is accumulated until the radiation cannot be considered as single photon,
664 i.e. the radiation probability exceeds some threshold
665 fSinglePhotonRadiationProbabilityLimit.
666 Typically fSinglePhotonRadiationProbabilityLimit should be between 0.01 and 0.1.
667 Then the photon generation as well as its parameters are simulated
668 in SetPhotonProductionParameters()
669
670 Finally if radiation took place, this function sets the new particle parameters
671 (the parameters at the point of radiation emission)
672 fNewParticleEnergy; fNewParticleAngleX;
673 fNewParticleAngleY; NewStep; fNewGlobalTime; fNewParticleCoordinateXYZ;
674
675 In order to control if the trajectory part reaches this limit,
676 one needs to divide it into smaller pieces (consisting of
677 fNSmallTrajectorySteps elements, typically several hundred) and to calculate the total
678 radiation probability along the trajectory part after each small piece accomplished.
679 If it exceeds the fSinglePhotonRadiationProbabilityLimit,
680 the trajectory part is over and the radiation can be generated.
681
682 In order to calculate the radiation probability the Baier-Katkov integral is called
683 after each small piece of the trajectory. It is summarized with the integral along
684 the previous small pieces for this part of the trajectory.
685
686 To speed-up the simulations, the photon angles are generated only once
687 at the start of the trajectory part.
688 */
689
690 //flag if a photon was produced
691 G4bool flagPhotonProduced = false;
692
693 //adding the next trajectory element to the vectors
694 fParticleAnglesX.push_back(angleX);
695 fParticleAnglesY.push_back(angleY);
696 fScatteringAnglesX.push_back(angleScatteringX);
697 fScatteringAnglesY.push_back(angleScatteringY);
698 fSteps.push_back(step);
699 fGlobalTimes.push_back(globalTime);
700 fParticleCoordinatesXYZ.push_back(coordinateXYZ);
701
702 G4double imax = fSteps.size();
703 if((imax==fImin0+fNSmallTrajectorySteps)||flagEndTrajectory)
704 {
705 //set the angular limits at the start of the trajectory part
706 if(fImin0==0)
707 {
708 //radiation within the angle = +-fRadiationAngleFactor/gamma
709 G4double radiationAngleLimit=fRadiationAngleFactor*mass/etotal;
710
711 SetPhotonSamplingParameters(etotal-mass,
712 *std::min_element(fParticleAnglesX.begin(),
713 fParticleAnglesX.end())-radiationAngleLimit,
714 *std::max_element(fParticleAnglesX.begin(),
715 fParticleAnglesX.end())+radiationAngleLimit,
716 *std::min_element(fParticleAnglesY.begin(),
717 fParticleAnglesY.end())-radiationAngleLimit,
718 *std::max_element(fParticleAnglesY.begin(),
719 fParticleAnglesY.end())+radiationAngleLimit);
720
721 //calculation of angles of photon emission
722 //(these angles are integration variables, Monte Carlo integration)
723 GeneratePhotonSampling();
724 }
725
726 //calculate Baier-Katkov integral after this
727 //small piece of trajectory (fNSmallTrajectorySteps elements)
728 fTotalRadiationProbability = RadIntegral(etotal,mass,
729 fParticleAnglesX,fParticleAnglesY,
730 fScatteringAnglesX,fScatteringAnglesY,
731 fSteps, fImin0);
732
733 //setting the last element of this small trajectory piece
734 fImin0 = imax;
735 fImax0.push_back(imax*1.);
736
737 //if the radiation probability is high enough or if we are finishing
738 //our trajectory => to simulate the photon emission
739 if(fTotalRadiationProbability>fSinglePhotonRadiationProbabilityLimit||
740 flagEndTrajectory)
741 {
742 fItrajectories += 1; //count this trajectory !!!correction 19.07.2023
743
744 flagPhotonProduced = SetPhotonProductionParameters(etotal,mass);
745
746 // correction 19.07.2023 fItrajectories += 1; //count this trajectory
747
748 //reinitialize intermediate integrals fFa, fSs, fSc, fSsx, fSsy, fScx, fScy;
749 //reset radiation integral internal variables to defaults;
750 //reset the trajectory and radiation probability along the trajectory
752 }
753 }
754
755 return flagPhotonProduced;
756}
void ResetRadIntegral()

◆ GeneratePhoton()

void G4BaierKatkov::GeneratePhoton ( G4FastStep & fastStep)

generates secondary photon belonging to fastStep with variables photon energy, momentum direction, coordinates and global time CALCULATED IN DoRadiation => USE IT ONLY AFTER DoRadiation returns true

Definition at line 760 of file G4BaierKatkov.cc.

761{
762 const G4DynamicParticle theGamma = G4DynamicParticle(G4Gamma::Gamma(),
763 PhMomentumDirection,fEph0);
764 //generation of a secondary photon
765 fastStep.CreateSecondaryTrack(theGamma,fNewParticleCoordinateXYZ,fNewGlobalTime,true);
766}
G4Track * CreateSecondaryTrack(const G4DynamicParticle &, G4ThreeVector, G4ThreeVector, G4double, G4bool localCoordinates=true)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81

◆ GetNewGlobalTime()

G4double G4BaierKatkov::GetNewGlobalTime ( )
inline

Definition at line 102 of file G4BaierKatkov.hh.

102{return fNewGlobalTime;}

◆ GetParticleNewAngleX()

G4double G4BaierKatkov::GetParticleNewAngleX ( )
inline

Definition at line 100 of file G4BaierKatkov.hh.

100{return fNewParticleAngleX;}

◆ GetParticleNewAngleY()

G4double G4BaierKatkov::GetParticleNewAngleY ( )
inline

Definition at line 101 of file G4BaierKatkov.hh.

101{return fNewParticleAngleY;}

◆ GetParticleNewCoordinateXYZ()

const G4ThreeVector & G4BaierKatkov::GetParticleNewCoordinateXYZ ( )
inline

Definition at line 103 of file G4BaierKatkov.hh.

103{return fNewParticleCoordinateXYZ;}

◆ GetParticleNewTotalEnergy()

G4double G4BaierKatkov::GetParticleNewTotalEnergy ( )
inline

get new parameters of the particle (the parameters at the point of radiation emission) needs calculation of DoRadiation first

Definition at line 99 of file G4BaierKatkov.hh.

99{return fNewParticleEnergy;}

◆ GetPhotonEnergyInSpectrum()

const std::vector< G4double > & G4BaierKatkov::GetPhotonEnergyInSpectrum ( )
inline

get photon energies (x-value in spectrum)

Definition at line 106 of file G4BaierKatkov.hh.

107 {return fPhotonEnergyInSpectrum;}

◆ GetSinglePhotonRadiationProbabilityLimit()

G4double G4BaierKatkov::GetSinglePhotonRadiationProbabilityLimit ( )
inline

get functions

You may call DoRadiation at each step of your trajectory CAUTION: please ensure that your steps are physically small enough for calculation of the radiation type you are interested in CAUTION: do ResetRadIntegral() before the start of a new trajectory

1) change some model defaults if necessary (SetSinglePhotonRadiationProbabilityLimit, SetNSmallTrajectorySteps, SetSpectrumEnergyRange) 2) call DoRadiation at each step of your trajectory 3) if DoRadiation returns TRUE, this means that a photon is produced (not added as a secondary yet) and its parameters are calculated. 4) You may generate a new photon using GeneratePhoton either with the parameters calculated in DoRadiation or your own parameters. CAUTION: By now GeneratePhoton works only for a FastSim model 5) Use GetPhotonEnergyInSpectrum() and GetTotalSpectrum() to return calculated total spectrum (all the photons altogether) Caution: is not normalized on the event number 6) Get the charged particle parameters in the radiation point: GetParticleNewTotalEnergy(), GetParticleNewAngleX(), GetParticleNewAngleY(), GetNewGlobalTime(), GetParticleNewCoordinateXYZ() get maximal radiation probability to preserve single photon radiation

Definition at line 88 of file G4BaierKatkov.hh.

89 {return fSinglePhotonRadiationProbabilityLimit;}

◆ GetTotalRadiationProbability()

G4double G4BaierKatkov::GetTotalRadiationProbability ( )
inline

total probability of radiation: needs calculation of DoRadiation first

CAUTION! : use the get functions below ONLY AFTER the call of DoRadiation and ONLY IF IT RETURNS true

Definition at line 95 of file G4BaierKatkov.hh.

95{return fTotalRadiationProbability;}

◆ GetTotalSpectrum()

const std::vector< G4double > & G4BaierKatkov::GetTotalSpectrum ( )
inline

get fTotalSpectrum after finishing the trajectory part with DoRadiation

Definition at line 110 of file G4BaierKatkov.hh.

110{return fTotalSpectrum;}

◆ ResetRadIntegral()

void G4BaierKatkov::ResetRadIntegral ( )

reinitialize intermediate integrals fFa, fSs, fSc, fSsx, fSsy, fScx, fScy; reset radiation integral internal variables to defaults; reset the trajectory and radiation probability along the trajectory

Definition at line 56 of file G4BaierKatkov.cc.

57{
58 fAccumSpectrum.clear();
59
60 //reinitialize intermediate integrals with zeros
61 fFa.resize(fNMCPhotons);
62 fSs.resize(fNMCPhotons);
63 fSc.resize(fNMCPhotons);
64 fSsx.resize(fNMCPhotons);
65 fSsy.resize(fNMCPhotons);
66 fScx.resize(fNMCPhotons);
67 fScy.resize(fNMCPhotons);
68 std::fill(fFa.begin(), fFa.end(), 0.);
69 std::fill(fSs.begin(), fSs.end(), 0.);
70 std::fill(fSc.begin(), fSc.end(), 0.);
71 std::fill(fSsx.begin(), fSsx.end(), 0.);
72 std::fill(fSsy.begin(), fSsy.end(), 0.);
73 std::fill(fScx.begin(), fScx.end(), 0.);
74 std::fill(fScy.begin(), fScy.end(), 0.);
75
76 //Reset radiation integral internal variables to defaults
77 fMeanPhotonAngleX =0.; //average angle of
78 //radiated photon direction in sampling, x-plane
79 fParamPhotonAngleX=1.e-3*rad; //a parameter of
80 //radiated photon sampling distribution, x-plane
81 fMeanPhotonAngleY =0.; //average angle of
82 //radiated photon direction in sampling, y-plane
83 fParamPhotonAngleY=1.e-3*rad; //a parameter of
84 //radiated photon sampling distribution, y-plane
85
86 fImin0 = 0;//set the first vector element to 0
87
88 //reset the trajectory
89 fParticleAnglesX.clear();
90 fParticleAnglesY.clear();
91 fScatteringAnglesX.clear();
92 fScatteringAnglesY.clear();
93 fSteps.clear();
94 fGlobalTimes.clear();
95 fParticleCoordinatesXYZ.clear();
96
97 //resets the vector of element numbers at the trajectory start
98 fImax0.clear();
99 //sets 0 element of the vector of element numbers
100 fImax0.push_back(0.);
101
102 //resets the radiation probability
103 fTotalRadiationProbabilityAlongTrajectory.clear();
104 //sets the radiation probability at the trajectory start
105 fTotalRadiationProbabilityAlongTrajectory.push_back(0.);
106}

Referenced by DoRadiation(), and SetSpectrumEnergyRange().

◆ SetMaxPhotonEnergy()

void G4BaierKatkov::SetMaxPhotonEnergy ( G4double emax)
inline

Definition at line 149 of file G4BaierKatkov.hh.

149 {SetSpectrumEnergyRange(fMinPhotonEnergy,
150 emax,
151 fNBinsSpectrum);}

◆ SetMinPhotonEnergy()

void G4BaierKatkov::SetMinPhotonEnergy ( G4double emin)
inline

SetSpectrumEnergyRange also calls ResetRadIntegral()

Definition at line 146 of file G4BaierKatkov.hh.

147 fMaxPhotonEnergy,
148 fNBinsSpectrum);}

◆ SetNBinsSpectrum()

void G4BaierKatkov::SetNBinsSpectrum ( G4int nbin)
inline

Definition at line 152 of file G4BaierKatkov.hh.

152 {SetSpectrumEnergyRange(fMinPhotonEnergy,
153 fMaxPhotonEnergy,
154 nbin);}

◆ SetNSmallTrajectorySteps()

void G4BaierKatkov::SetNSmallTrajectorySteps ( G4int nSmallTrajectorySteps)
inline

number of steps in a trajectory small piece before the next call of the radiation integral

Definition at line 120 of file G4BaierKatkov.hh.

121 {fNSmallTrajectorySteps = nSmallTrajectorySteps;}

◆ SetRadiationAngleFactor()

void G4BaierKatkov::SetRadiationAngleFactor ( G4double radiationAngleFactor)
inline

setting the number of radiation angles 1/gamma, defining the width of the angular distribution of photon sampling in the Baier-Katkov Integral

Definition at line 134 of file G4BaierKatkov.hh.

135 {fRadiationAngleFactor = radiationAngleFactor;}

◆ SetSamplingPhotonsNumber()

void G4BaierKatkov::SetSamplingPhotonsNumber ( G4int nPhotons)
inline

setting the number of photons in sampling of Baier-Katkov Integral (MC integration by photon energy and angles <=> photon momentum)

Definition at line 130 of file G4BaierKatkov.hh.

130{fNMCPhotons = nPhotons;}

◆ SetSinglePhotonRadiationProbabilityLimit()

void G4BaierKatkov::SetSinglePhotonRadiationProbabilityLimit ( G4double wmax)
inline

set functions

set maximal radiation probability to preserve single photon radiation

Definition at line 115 of file G4BaierKatkov.hh.

116 {fSinglePhotonRadiationProbabilityLimit = wmax;}

◆ SetSpectrumEnergyRange()

void G4BaierKatkov::SetSpectrumEnergyRange ( G4double emin,
G4double emax,
G4int numberOfBins )

CAUTION, the bins width is logarithmic Do not worry if the maximal energy > particle energy. This elements of spectrum with non-physical energies will not be processed (they will be 0).

Definition at line 110 of file G4BaierKatkov.cc.

113{
114 fMinPhotonEnergy = emin;
115 fMaxPhotonEnergy = emax;
116 fNBinsSpectrum = numberOfBins;
117
118 fLogEmaxdEmin = G4Log(fMaxPhotonEnergy/fMinPhotonEnergy);
119
120 //in initializing fNPhotonsPerBin
121 fNPhotonsPerBin.resize(fNBinsSpectrum);
122 std::fill(fNPhotonsPerBin.begin(), fNPhotonsPerBin.end(), 0);
123
124 //initializing the Spectrum
125 fSpectrum.resize(fNBinsSpectrum);
126 std::fill(fSpectrum.begin(), fSpectrum.end(), 0);
127
128 //initializing the fAccumTotalSpectrum
129 fAccumTotalSpectrum.resize(fNBinsSpectrum);
130 std::fill(fAccumTotalSpectrum.begin(), fAccumTotalSpectrum.end(), 0);
131
132 //initializing the fTotalSpectrum
133 fTotalSpectrum.resize(fNBinsSpectrum);
134 std::fill(fTotalSpectrum.begin(), fTotalSpectrum.end(), 0);
135
136 fPhotonEnergyInSpectrum.clear();
137 for (G4int j=0;j<fNBinsSpectrum;j++)
138 {
139 //bin position (mean between 2 bin limits)
140 fPhotonEnergyInSpectrum.push_back(fMinPhotonEnergy*
141 (std::exp(fLogEmaxdEmin*j/fNBinsSpectrum)+
142 std::exp(fLogEmaxdEmin*(j+1)/fNBinsSpectrum))/2.);
143 }
144
145 fItrajectories = 0;
146
148}

Referenced by G4BaierKatkov(), SetMaxPhotonEnergy(), SetMinPhotonEnergy(), and SetNBinsSpectrum().

◆ SetVirtualCollimator()

void G4BaierKatkov::SetVirtualCollimator ( G4double virtualCollimatorAngularDiameter)
inline

Virtual collimator masks the selection of photon angles in fTotalSpectrum Virtual collimator doesn't influence on Geant4 simulations.

Definition at line 163 of file G4BaierKatkov.hh.

164 {fVirtualCollimatorAngularDiameter=virtualCollimatorAngularDiameter;}

The documentation for this class was generated from the following files: