48const G4double G4SynchrotronRadiationInMat::fIntegralProbabilityOfSR[200] = {
49 1.000000e+00, 9.428859e-01, 9.094095e-01, 8.813971e-01, 8.565154e-01,
50 8.337008e-01, 8.124961e-01, 7.925217e-01, 7.735517e-01, 7.554561e-01,
51 7.381233e-01, 7.214521e-01, 7.053634e-01, 6.898006e-01, 6.747219e-01,
52 6.600922e-01, 6.458793e-01, 6.320533e-01, 6.185872e-01, 6.054579e-01,
53 5.926459e-01, 5.801347e-01, 5.679103e-01, 5.559604e-01, 5.442736e-01,
54 5.328395e-01, 5.216482e-01, 5.106904e-01, 4.999575e-01, 4.894415e-01,
55 4.791351e-01, 4.690316e-01, 4.591249e-01, 4.494094e-01, 4.398800e-01,
56 4.305320e-01, 4.213608e-01, 4.123623e-01, 4.035325e-01, 3.948676e-01,
57 3.863639e-01, 3.780179e-01, 3.698262e-01, 3.617858e-01, 3.538933e-01,
58 3.461460e-01, 3.385411e-01, 3.310757e-01, 3.237474e-01, 3.165536e-01,
59 3.094921e-01, 3.025605e-01, 2.957566e-01, 2.890784e-01, 2.825237e-01,
60 2.760907e-01, 2.697773e-01, 2.635817e-01, 2.575020e-01, 2.515365e-01,
61 2.456834e-01, 2.399409e-01, 2.343074e-01, 2.287812e-01, 2.233607e-01,
62 2.180442e-01, 2.128303e-01, 2.077174e-01, 2.027040e-01, 1.977885e-01,
63 1.929696e-01, 1.882457e-01, 1.836155e-01, 1.790775e-01, 1.746305e-01,
64 1.702730e-01, 1.660036e-01, 1.618212e-01, 1.577243e-01, 1.537117e-01,
65 1.497822e-01, 1.459344e-01, 1.421671e-01, 1.384791e-01, 1.348691e-01,
66 1.313360e-01, 1.278785e-01, 1.244956e-01, 1.211859e-01, 1.179483e-01,
67 1.147818e-01, 1.116850e-01, 1.086570e-01, 1.056966e-01, 1.028026e-01,
68 9.997405e-02, 9.720975e-02, 9.450865e-02, 9.186969e-02, 8.929179e-02,
69 8.677391e-02, 8.431501e-02, 8.191406e-02, 7.957003e-02, 7.728192e-02,
70 7.504872e-02, 7.286944e-02, 7.074311e-02, 6.866874e-02, 6.664538e-02,
71 6.467208e-02, 6.274790e-02, 6.087191e-02, 5.904317e-02, 5.726079e-02,
72 5.552387e-02, 5.383150e-02, 5.218282e-02, 5.057695e-02, 4.901302e-02,
73 4.749020e-02, 4.600763e-02, 4.456450e-02, 4.315997e-02, 4.179325e-02,
74 4.046353e-02, 3.917002e-02, 3.791195e-02, 3.668855e-02, 3.549906e-02,
75 3.434274e-02, 3.321884e-02, 3.212665e-02, 3.106544e-02, 3.003452e-02,
76 2.903319e-02, 2.806076e-02, 2.711656e-02, 2.619993e-02, 2.531021e-02,
77 2.444677e-02, 2.360897e-02, 2.279620e-02, 2.200783e-02, 2.124327e-02,
78 2.050194e-02, 1.978324e-02, 1.908662e-02, 1.841151e-02, 1.775735e-02,
79 1.712363e-02, 1.650979e-02, 1.591533e-02, 1.533973e-02, 1.478250e-02,
80 1.424314e-02, 1.372117e-02, 1.321613e-02, 1.272755e-02, 1.225498e-02,
81 1.179798e-02, 1.135611e-02, 1.092896e-02, 1.051609e-02, 1.011712e-02,
82 9.731635e-03, 9.359254e-03, 8.999595e-03, 8.652287e-03, 8.316967e-03,
83 7.993280e-03, 7.680879e-03, 7.379426e-03, 7.088591e-03, 6.808051e-03,
84 6.537491e-03, 6.276605e-03, 6.025092e-03, 5.782661e-03, 5.549027e-03,
85 5.323912e-03, 5.107045e-03, 4.898164e-03, 4.697011e-03, 4.503336e-03,
86 4.316896e-03, 4.137454e-03, 3.964780e-03, 3.798649e-03, 3.638843e-03,
87 3.485150e-03, 3.337364e-03, 3.195284e-03, 3.058715e-03, 2.927469e-03,
88 2.801361e-03, 2.680213e-03, 2.563852e-03, 2.452110e-03, 2.344824e-03
99 , LowestKineticEnergy(10. * keV)
102 , fVerboseLevel(verboseLevel)
110 CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
111 PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
112 fPsiGamma = fEta = fOrderAngleK = 0.0;
149 if(KineticEnergy < LowestKineticEnergy || gamma < 1.0e3)
154 const G4Field* pField =
nullptr;
157 G4bool fieldExertsForce =
false;
159 if((particleCharge != 0.0))
164 if(fieldMgr !=
nullptr)
175 G4double globPosVec[4], FieldValueVec[6];
177 globPosVec[0] = globPosition.
x();
178 globPosVec[1] = globPosition.
y();
179 globPosVec[2] = globPosition.
z();
185 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
194 MeanFreePath = fLambdaConst * beta / perpB;
201 if(fVerboseLevel > 0)
203 G4cout <<
"G4SynchrotronRadiationInMat::MeanFreePath = " << MeanFreePath / m
228 const G4Field* pField =
nullptr;
231 G4bool fieldExertsForce =
false;
233 if((particleCharge != 0.0))
236 if(fieldMgr !=
nullptr)
246 G4double globPosVec[4], FieldValueVec[6];
247 globPosVec[0] = globPosition.
x();
248 globPosVec[1] = globPosition.
y();
249 globPosVec[2] = globPosition.
z();
254 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
264 if(fVerboseLevel > 0)
266 G4cout <<
"SR photon energy = " << energyOfSR / keV <<
" keV" <<
G4endl;
270 if(energyOfSR <= 0.0)
279 G4double cosTheta, sinTheta, fcos, beta;
284 fcos = (1 + cosTheta * cosTheta) * 0.5;
289 beta = std::sqrt(1. - 1. / (gamma * gamma));
291 cosTheta = (cosTheta + beta) / (1. + beta * cosTheta);
298 sinTheta = std::sqrt(1. - cosTheta * cosTheta);
302 G4double dirx = sinTheta * std::cos(Phi);
303 G4double diry = sinTheta * std::sin(Phi);
307 gammaDirection.
rotateUz(particleDirection);
311 gammaPolarization = gammaPolarization.
unit();
316 aGamma->SetPolarization(gammaPolarization.
x(), gammaPolarization.
y(),
317 gammaPolarization.
z());
322 G4double newKinEnergy = kineticEnergy - energyOfSR;
324 if(newKinEnergy > 0.)
374 const G4Field* pField =
nullptr;
377 G4bool fieldExertsForce =
false;
379 if((particleCharge != 0.0))
382 if(fieldMgr !=
nullptr)
392 G4double globPosVec[3], FieldValueVec[3];
394 globPosVec[0] = globPosition.
x();
395 globPosVec[1] = globPosition.
y();
396 globPosVec[2] = globPosition.
z();
400 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
409 for(i = 0; i < 200; ++i)
411 if(random >= fIntegralProbabilityOfSR[i])
414 energyOfSR = 0.0001 * i * i * fEnergyConst * gamma * gamma * perpB;
417 if(energyOfSR <= 0.0)
435 static constexpr G4int iMax = 200;
440 for(i = 0; i < iMax; ++i)
442 if(random >= fIntegralProbabilityOfSR[i])
464 G4double result, hypCos2, hypCos = std::cosh(t);
466 hypCos2 = hypCos * hypCos;
467 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos);
491 result = integral.Laguerre(
494 result *= 3. / 5. / pi;
503 G4double result, hypCos = std::cosh(t);
505 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos);
529 result = integral.Laguerre(
532 result *= 9. * std::sqrt(3.) * ksi / 8. / pi;
540 G4double result, hypCos = std::cosh(t);
543 std::cosh(fOrderAngleK * t) * std::exp(t - fEta * hypCos);
563 result = integral.Laguerre(
573 G4double result, funK, funK2, gpsi2 = gpsi * gpsi;
576 fEta = 0.5 * fKsi * (1. + gpsi2) * std::sqrt(1. + gpsi2);
578 fOrderAngleK = 1. / 3.;
582 result = gpsi2 * funK2 / (1. + gpsi2);
584 fOrderAngleK = 2. / 3.;
589 result *= (1. + gpsi2) * fKsi;
G4double condition(const G4ErrorSymMatrix &m)
CLHEP::Hep3Vector G4ThreeVector
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetPDGCharge() const
static G4int GetModelID(const G4int modelIndex)
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
const G4TouchableHandle & GetTouchableHandle() const
G4StepPoint * GetPostStepPoint() const
G4double GetEnergyProbSR(G4double)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &Step) override
static G4double GetEnergyConst()
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
static G4double GetLambdaConst()
G4double GetAngleNumberAtGammaKsi(G4double)
G4SynchrotronRadiationInMat(const G4String &processName="SynchrotronRadiation", G4ProcessType type=fElectromagnetic)
G4double GetRandomEnergySR(G4double, G4double)
G4double GetProbSpectrumSRforEnergy(G4double)
G4double GetIntProbSR(G4double)
G4double GetIntegrandForAngleK(G4double)
G4double GetPhotonEnergy(const G4Track &trackData, const G4Step &stepData)
G4bool IsApplicable(const G4ParticleDefinition &) override
G4double GetAngleK(G4double)
~G4SynchrotronRadiationInMat()
G4double GetProbSpectrumSRforInt(G4double)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
void SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
void SetProcessSubType(G4int)