Geant4 11.4.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 SetRoundVirtualCollimator (G4double virtualCollimatorAngularRadius, G4double virtualCollimatorAngularCenterX=0., G4double virtualCollimatorAngularCenterY=0.)
void SetEllipticVirtualCollimator (G4double virtualCollimatorAngularRadiusX, G4double virtualCollimatorAngularRadiusY, G4double virtualCollimatorAngularCenterX=0., G4double virtualCollimatorAngularCenterY=0.)
void SetRectangularVirtualCollimator (G4double virtualCollimatorAngularHalfWidthX, G4double virtualCollimatorAngularHalfWidthY, G4double virtualCollimatorAngularCenterX=0., G4double virtualCollimatorAngularCenterY=0.)
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 37 of file G4BaierKatkov.cc.

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

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

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

768{
769 const G4DynamicParticle theGamma = G4DynamicParticle(G4Gamma::Gamma(),
770 PhMomentumDirection,fEph0);
771 //generation of a secondary photon
772 fastStep.CreateSecondaryTrack(theGamma,fNewParticleCoordinateXYZ,fNewGlobalTime,true);
773}
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 54 of file G4BaierKatkov.cc.

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

Referenced by DoRadiation(), and SetSpectrumEnergyRange().

◆ SetEllipticVirtualCollimator()

void G4BaierKatkov::SetEllipticVirtualCollimator ( G4double virtualCollimatorAngularRadiusX,
G4double virtualCollimatorAngularRadiusY,
G4double virtualCollimatorAngularCenterX = 0.,
G4double virtualCollimatorAngularCenterY = 0. )
inline

Definition at line 174 of file G4BaierKatkov.hh.

178 {fVirtualCollimatorAngularHalfWidthX2 = virtualCollimatorAngularRadiusX*
179 virtualCollimatorAngularRadiusX;
180 fVirtualCollimatorAngularHalfWidthY2 = virtualCollimatorAngularRadiusY*
181 virtualCollimatorAngularRadiusY;
182 fVirtualCollimatorAngularCenterX = virtualCollimatorAngularCenterX;
183 fVirtualCollimatorAngularCenterY = virtualCollimatorAngularCenterY;
184 fVirtualCollimatorTypeID = 1;}

◆ 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;}

◆ SetRectangularVirtualCollimator()

void G4BaierKatkov::SetRectangularVirtualCollimator ( G4double virtualCollimatorAngularHalfWidthX,
G4double virtualCollimatorAngularHalfWidthY,
G4double virtualCollimatorAngularCenterX = 0.,
G4double virtualCollimatorAngularCenterY = 0. )
inline

Definition at line 186 of file G4BaierKatkov.hh.

190 {fVirtualCollimatorAngularHalfWidthX = virtualCollimatorAngularHalfWidthX;
191 fVirtualCollimatorAngularHalfWidthY = virtualCollimatorAngularHalfWidthY;
192 fVirtualCollimatorAngularCenterX = virtualCollimatorAngularCenterX;
193 fVirtualCollimatorAngularCenterY = virtualCollimatorAngularCenterY;
194 fVirtualCollimatorTypeID = 2;}

◆ SetRoundVirtualCollimator()

void G4BaierKatkov::SetRoundVirtualCollimator ( G4double virtualCollimatorAngularRadius,
G4double virtualCollimatorAngularCenterX = 0.,
G4double virtualCollimatorAngularCenterY = 0. )
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.

166 {fVirtualCollimatorAngularHalfWidthX2 = virtualCollimatorAngularRadius*
167 virtualCollimatorAngularRadius;
168 fVirtualCollimatorAngularHalfWidthY2 = fVirtualCollimatorAngularHalfWidthX2;
169
170 fVirtualCollimatorAngularCenterX = virtualCollimatorAngularCenterX;
171 fVirtualCollimatorAngularCenterY = virtualCollimatorAngularCenterY;
172 fVirtualCollimatorTypeID = 1;}

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

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

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


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