Geant4 11.2.2
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 (G4double 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 47 of file G4BaierKatkov.hh.

Constructor & Destructor Documentation

◆ G4BaierKatkov()

G4BaierKatkov::G4BaierKatkov ( )

Definition at line 35 of file G4BaierKatkov.cc.

36{
37 //sets the default spectrum energy range of integration and
38 //calls ResetRadIntegral()
39 SetSpectrumEnergyRange(0.1*MeV,1.*GeV,110);
40
41 //Do not worry if the maximal energy > particle energy
42 //this elements of spectrum with non-physical energies
43 //will not be processed (they will be 0)
44
45 G4cout << " "<< G4endl;
46 G4cout << "G4BaierKatkov model is activated."<< G4endl;
47 G4cout << " "<< G4endl;
48}
#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 148 of file G4BaierKatkov.cc.

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

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

Referenced by G4ChannelingFastSimModel::DoIt().

◆ 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 756 of file G4BaierKatkov.cc.

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

Referenced by G4ChannelingFastSimModel::DoIt().

◆ GetNewGlobalTime()

G4double G4BaierKatkov::GetNewGlobalTime ( )
inline

Definition at line 98 of file G4BaierKatkov.hh.

98{return fNewGlobalTime;}

Referenced by G4ChannelingFastSimModel::DoIt().

◆ GetParticleNewAngleX()

G4double G4BaierKatkov::GetParticleNewAngleX ( )
inline

Definition at line 96 of file G4BaierKatkov.hh.

96{return fNewParticleAngleX;}

Referenced by G4ChannelingFastSimModel::DoIt().

◆ GetParticleNewAngleY()

G4double G4BaierKatkov::GetParticleNewAngleY ( )
inline

Definition at line 97 of file G4BaierKatkov.hh.

97{return fNewParticleAngleY;}

Referenced by G4ChannelingFastSimModel::DoIt().

◆ GetParticleNewCoordinateXYZ()

const G4ThreeVector & G4BaierKatkov::GetParticleNewCoordinateXYZ ( )
inline

Definition at line 99 of file G4BaierKatkov.hh.

99{return fNewParticleCoordinateXYZ;}

Referenced by G4ChannelingFastSimModel::DoIt().

◆ 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 95 of file G4BaierKatkov.hh.

95{return fNewParticleEnergy;}

Referenced by G4ChannelingFastSimModel::DoIt().

◆ GetPhotonEnergyInSpectrum()

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

get photon energies (x-value in spectrum)

Definition at line 102 of file G4BaierKatkov.hh.

103 {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 84 of file G4BaierKatkov.hh.

85 {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 91 of file G4BaierKatkov.hh.

91{return fTotalRadiationProbability;}

◆ GetTotalSpectrum()

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

get fTotalSpectrum after finishing the trajectory part with DoRadiation

Definition at line 106 of file G4BaierKatkov.hh.

106{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 52 of file G4BaierKatkov.cc.

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

Referenced by G4ChannelingFastSimModel::DoIt(), DoRadiation(), and SetSpectrumEnergyRange().

◆ SetMaxPhotonEnergy()

void G4BaierKatkov::SetMaxPhotonEnergy ( G4double emax)
inline

Definition at line 145 of file G4BaierKatkov.hh.

145 {SetSpectrumEnergyRange(fMinPhotonEnergy,
146 emax,
147 fNBinsSpectrum);}

◆ SetMinPhotonEnergy()

void G4BaierKatkov::SetMinPhotonEnergy ( G4double emin)
inline

SetSpectrumEnergyRange also calls ResetRadIntegral()

Definition at line 142 of file G4BaierKatkov.hh.

143 fMaxPhotonEnergy,
144 fNBinsSpectrum);}

◆ SetNBinsSpectrum()

void G4BaierKatkov::SetNBinsSpectrum ( G4int nbin)
inline

Definition at line 148 of file G4BaierKatkov.hh.

148 {SetSpectrumEnergyRange(fMinPhotonEnergy,
149 fMaxPhotonEnergy,
150 nbin);}

◆ SetNSmallTrajectorySteps()

void G4BaierKatkov::SetNSmallTrajectorySteps ( G4double nSmallTrajectorySteps)
inline

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

Definition at line 116 of file G4BaierKatkov.hh.

117 {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 130 of file G4BaierKatkov.hh.

131 {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 126 of file G4BaierKatkov.hh.

126{fNMCPhotons = nPhotons;}

◆ SetSinglePhotonRadiationProbabilityLimit()

void G4BaierKatkov::SetSinglePhotonRadiationProbabilityLimit ( G4double wmax)
inline

set functions

set maximal radiation probability to preserve single photon radiation

Definition at line 111 of file G4BaierKatkov.hh.

112 {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 106 of file G4BaierKatkov.cc.

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

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 159 of file G4BaierKatkov.hh.

160 {fVirtualCollimatorAngularDiameter=virtualCollimatorAngularDiameter;}

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