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

#include <G4SynchrotronRadiation.hh>

+ Inheritance diagram for G4SynchrotronRadiation:

Public Member Functions

 G4SynchrotronRadiation (const G4String &pName="SynRad", G4ProcessType type=fElectromagnetic)
 
virtual ~G4SynchrotronRadiation ()
 
G4SynchrotronRadiationoperator= (const G4SynchrotronRadiation &right)=delete
 
 G4SynchrotronRadiation (const G4SynchrotronRadiation &)=delete
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step) override
 
G4double GetPhotonEnergy (const G4Track &trackData, const G4Step &stepData)
 
G4double GetRandomEnergySR (G4double, G4double, G4double)
 
G4double InvSynFracInt (G4double x)
 
G4double Chebyshev (G4double a, G4double b, const G4double c[], G4int n, G4double x)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void ProcessDescription (std::ostream &) const override
 
void DumpInfo () const override
 
void SetAngularGenerator (G4VEmAngularDistribution *p)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 50 of file G4SynchrotronRadiation.hh.

Constructor & Destructor Documentation

◆ G4SynchrotronRadiation() [1/2]

G4SynchrotronRadiation::G4SynchrotronRadiation ( const G4String pName = "SynRad",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 54 of file G4SynchrotronRadiation.cc.

56 : G4VDiscreteProcess(processName, type)
57 , theGamma(G4Gamma::Gamma())
58{
59 G4TransportationManager* transportMgr =
61
62 fFieldPropagator = transportMgr->GetPropagatorInField();
63
64 secID = G4PhysicsModelCatalog::GetModelID("model_SynRad");
66 verboseLevel = 1;
67 FirstTime = true;
68 FirstTime1 = true;
69 genAngle = nullptr;
71 theManager = G4LossTableManager::Instance();
72 theManager->Register(this);
73}
@ fSynchrotronRadiation
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
static G4int GetModelID(const G4int modelIndex)
void SetAngularGenerator(G4VEmAngularDistribution *p)
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4SynchrotronRadiation()

G4SynchrotronRadiation::~G4SynchrotronRadiation ( )
virtual

Definition at line 77 of file G4SynchrotronRadiation.cc.

78{
79 delete genAngle;
80 theManager->DeRegister(this);
81}
void DeRegister(G4VEnergyLossProcess *p)

◆ G4SynchrotronRadiation() [2/2]

G4SynchrotronRadiation::G4SynchrotronRadiation ( const G4SynchrotronRadiation )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4SynchrotronRadiation::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 419 of file G4SynchrotronRadiation.cc.

420{
421 if(0 < verboseLevel && &part == G4Electron::Electron())
423 // same for all particles, print only for one (electron)
424}
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
void ProcessDescription(std::ostream &) const override

◆ Chebyshev()

G4double G4SynchrotronRadiation::Chebyshev ( G4double  a,
G4double  b,
const G4double  c[],
G4int  n,
G4double  x 
)
inline

Definition at line 99 of file G4SynchrotronRadiation.hh.

102{
103 G4double y;
104 G4double y2 = 2.0 * (y = (2.0 * x - a - b) / (b - a)); // Change of variable.
105 G4double d = 0., dd = 0.;
106 for(G4int j = n - 1; j >= 1; --j) // Clenshaw's recurrence.
107 {
108 G4double sv = d;
109 d = y2 * d - dd + c[j];
110 dd = sv;
111 }
112 return y * d - dd + 0.5 * c[0];
113}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85

Referenced by InvSynFracInt().

◆ DumpInfo()

void G4SynchrotronRadiation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4SynchrotronRadiation.hh.

◆ GetMeanFreePath()

G4double G4SynchrotronRadiation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 103 of file G4SynchrotronRadiation.cc.

106{
107 // gives the MeanFreePath in Geant4 internal units
108 G4double MeanFreePath = DBL_MAX;
109
110 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
111
113
114 G4double gamma =
115 aDynamicParticle->GetTotalEnergy() / aDynamicParticle->GetMass();
116
117 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
118
119 if(gamma < 1.0e3 || 0.0 == particleCharge)
120 {
121 MeanFreePath = DBL_MAX;
122 }
123 else
124 {
125 G4ThreeVector FieldValue;
126 const G4Field* pField = nullptr;
127 G4bool fieldExertsForce = false;
128
129 G4FieldManager* fieldMgr =
130 fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
131
132 if(fieldMgr != nullptr)
133 {
134 // If the field manager has no field, there is no field !
135 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
136 }
137
138 if(fieldExertsForce)
139 {
140 pField = fieldMgr->GetDetectorField();
141 G4ThreeVector globPosition = trackData.GetPosition();
142
143 G4double globPosVec[4], FieldValueVec[6];
144
145 globPosVec[0] = globPosition.x();
146 globPosVec[1] = globPosition.y();
147 globPosVec[2] = globPosition.z();
148 globPosVec[3] = trackData.GetGlobalTime();
149
150 pField->GetFieldValue(globPosVec, FieldValueVec);
151
152 FieldValue =
153 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
154
155 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
156 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
157 G4double perpB = unitMcrossB.mag();
158
159 static const G4double fLambdaConst =
160 std::sqrt(3.0) * eplus / (2.5 * fine_structure_const * c_light);
161
162 if(perpB > 0.0)
163 {
164 MeanFreePath = fLambdaConst *
165 aDynamicParticle->GetDefinition()->GetPDGMass() /
166 (perpB * particleCharge * particleCharge);
167 }
168 if(verboseLevel > 0 && FirstTime)
169 {
170 G4cout << "G4SynchrotronRadiation::GetMeanFreePath "
171 << " for particle "
172 << aDynamicParticle->GetDefinition()->GetParticleName() << ":"
173 << '\n'
174 << " MeanFreePath = " << G4BestUnit(MeanFreePath, "Length")
175 << G4endl;
176 if(verboseLevel > 1)
177 {
178 G4ThreeVector pvec = aDynamicParticle->GetMomentum();
179 G4double Btot = FieldValue.getR();
180 G4double ptot = pvec.getR();
181 G4double rho = ptot / (MeV * c_light * Btot);
182 // full bending radius
183 G4double Theta = unitMomentum.theta(FieldValue);
184 // angle between particle and field
185 G4cout << " B = " << Btot / tesla << " Tesla"
186 << " perpB = " << perpB / tesla << " Tesla"
187 << " Theta = " << Theta
188 << " std::sin(Theta)=" << std::sin(Theta) << '\n'
189 << " ptot = " << G4BestUnit(ptot, "Energy")
190 << " rho = " << G4BestUnit(rho, "Length") << G4endl;
191 }
192 FirstTime = false;
193 }
194 }
195 }
196 return MeanFreePath;
197}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
double z() const
double theta() const
double x() const
double getR() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
#define DBL_MAX
Definition: templates.hh:62

◆ GetPhotonEnergy()

G4double G4SynchrotronRadiation::GetPhotonEnergy ( const G4Track trackData,
const G4Step stepData 
)

◆ GetRandomEnergySR()

G4double G4SynchrotronRadiation::GetRandomEnergySR ( G4double  gamma,
G4double  perpB,
G4double  mass_c2 
)

Definition at line 391 of file G4SynchrotronRadiation.cc.

394{
395 static const G4double fEnergyConst =
396 1.5 * c_light * c_light * eplus * hbar_Planck;
397 G4double Ecr = fEnergyConst * gamma * gamma * perpB / mass_c2;
398
399 if(verboseLevel > 0 && FirstTime1)
400 {
401 // mean and rms of photon energy
402 G4double Emean = 8. / (15. * std::sqrt(3.)) * Ecr;
403 G4double E_rms = std::sqrt(211. / 675.) * Ecr;
404 G4long prec = G4cout.precision();
405 G4cout << "G4SynchrotronRadiation::GetRandomEnergySR :" << '\n'
406 << std::setprecision(4) << " Ecr = " << G4BestUnit(Ecr, "Energy")
407 << '\n'
408 << " Emean = " << G4BestUnit(Emean, "Energy") << '\n'
409 << " E_rms = " << G4BestUnit(E_rms, "Energy") << G4endl;
410 FirstTime1 = false;
411 G4cout.precision(prec);
412 }
413
414 G4double energySR = Ecr * InvSynFracInt(G4UniformRand());
415 return energySR;
416}
long G4long
Definition: G4Types.hh:87
#define G4UniformRand()
Definition: Randomize.hh:52
G4double InvSynFracInt(G4double x)

Referenced by PostStepDoIt().

◆ InvSynFracInt()

G4double G4SynchrotronRadiation::InvSynFracInt ( G4double  x)

Definition at line 298 of file G4SynchrotronRadiation.cc.

300{
301 // from 0 to 0.7
302 static constexpr G4double aa1 = 0;
303 static constexpr G4double aa2 = 0.7;
304 static constexpr G4int ncheb1 = 27;
305 static constexpr G4double cheb1[ncheb1] = {
306 1.22371665676046468821, 0.108956475422163837267,
307 0.0383328524358594396134, 0.00759138369340257753721,
308 0.00205712048644963340914, 0.000497810783280019308661,
309 0.000130743691810302187818, 0.0000338168760220395409734,
310 8.97049680900520817728e-6, 2.38685472794452241466e-6,
311 6.41923109149104165049e-7, 1.73549898982749277843e-7,
312 4.72145949240790029153e-8, 1.29039866111999149636e-8,
313 3.5422080787089834182e-9, 9.7594757336403784905e-10,
314 2.6979510184976065731e-10, 7.480422622550977077e-11,
315 2.079598176402699913e-11, 5.79533622220841193e-12,
316 1.61856011449276096e-12, 4.529450993473807e-13,
317 1.2698603951096606e-13, 3.566117394511206e-14,
318 1.00301587494091e-14, 2.82515346447219e-15,
319 7.9680747949792e-16
320 };
321 // from 0.7 to 0.9132260271183847
322 static constexpr G4double aa3 = 0.9132260271183847;
323 static constexpr G4int ncheb2 = 27;
324 static constexpr G4double cheb2[ncheb2] = {
325 1.1139496701107756, 0.3523967429328067, 0.0713849171926623,
326 0.01475818043595387, 0.003381255637322462, 0.0008228057599452224,
327 0.00020785506681254216, 0.00005390169253706556, 0.000014250571923902464,
328 3.823880733161044e-6, 1.0381966089136036e-6, 2.8457557457837253e-7,
329 7.86223332179956e-8, 2.1866609342508474e-8, 6.116186259857143e-9,
330 1.7191233618437565e-9, 4.852755117740807e-10, 1.3749966961763457e-10,
331 3.908961987062447e-11, 1.1146253766895824e-11, 3.1868887323415814e-12,
332 9.134319791300977e-13, 2.6211077371181566e-13, 7.588643377757906e-14,
333 2.1528376972619e-14, 6.030906040404772e-15, 1.9549163926819867e-15
334 };
335 // Chebyshev with exp/log scale
336 // a = -Log[1 - SynFracInt[1]]; b = -Log[1 - SynFracInt[7]];
337 static constexpr G4double aa4 = 2.4444485538746025480;
338 static constexpr G4double aa5 = 9.3830728608909477079;
339 static constexpr G4int ncheb3 = 28;
340 static constexpr G4double cheb3[ncheb3] = {
341 1.2292683840435586977, 0.160353449247864455879,
342 -0.0353559911947559448721, 0.00776901561223573936985,
343 -0.00165886451971685133259, 0.000335719118906954279467,
344 -0.0000617184951079161143187, 9.23534039743246708256e-6,
345 -6.06747198795168022842e-7, -3.07934045961999778094e-7,
346 1.98818772614682367781e-7, -8.13909971567720135413e-8,
347 2.84298174969641838618e-8, -9.12829766621316063548e-9,
348 2.77713868004820551077e-9, -8.13032767247834023165e-10,
349 2.31128525568385247392e-10, -6.41796873254200220876e-11,
350 1.74815310473323361543e-11, -4.68653536933392363045e-12,
351 1.24016595805520752748e-12, -3.24839432979935522159e-13,
352 8.44601465226513952994e-14, -2.18647276044246803998e-14,
353 5.65407548745690689978e-15, -1.46553625917463067508e-15,
354 3.82059606377570462276e-16, -1.00457896653436912508e-16
355 };
356 static constexpr G4double aa6 = 33.122936966163038145;
357 static constexpr G4int ncheb4 = 27;
358 static constexpr G4double cheb4[ncheb4] = {
359 1.69342658227676741765, 0.0742766400841232319225,
360 -0.019337880608635717358, 0.00516065527473364110491,
361 -0.00139342012990307729473, 0.000378549864052022522193,
362 -0.000103167085583785340215, 0.0000281543441271412178337,
363 -7.68409742018258198651e-6, 2.09543221890204537392e-6,
364 -5.70493140367526282946e-7, 1.54961164548564906446e-7,
365 -4.19665599629607704794e-8, 1.13239680054166507038e-8,
366 -3.04223563379021441863e-9, 8.13073745977562957997e-10,
367 -2.15969415476814981374e-10, 5.69472105972525594811e-11,
368 -1.48844799572430829499e-11, 3.84901514438304484973e-12,
369 -9.82222575944247161834e-13, 2.46468329208292208183e-13,
370 -6.04953826265982691612e-14, 1.44055805710671611984e-14,
371 -3.28200813577388740722e-15, 6.96566359173765367675e-16,
372 -1.294122794852896275e-16
373 };
374
375 if(x < aa2)
376 return x * x * x * Chebyshev(aa1, aa2, cheb1, ncheb1, x);
377 else if(x < aa3)
378 return Chebyshev(aa2, aa3, cheb2, ncheb2, x);
379 else if(x < 1 - 0.0000841363)
380 {
381 G4double y = -G4Log(1 - x);
382 return y * Chebyshev(aa4, aa5, cheb3, ncheb3, y);
383 }
384 else
385 {
386 G4double y = -G4Log(1 - x);
387 return y * Chebyshev(aa5, aa6, cheb4, ncheb4, y);
388 }
389}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double Chebyshev(G4double a, G4double b, const G4double c[], G4int n, G4double x)

Referenced by GetRandomEnergySR().

◆ IsApplicable()

G4bool G4SynchrotronRadiation::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 94 of file G4SynchrotronRadiation.cc.

96{
97 return (particle.GetPDGCharge() != 0.0 && !particle.IsShortLived());
98}

◆ operator=()

G4SynchrotronRadiation & G4SynchrotronRadiation::operator= ( const G4SynchrotronRadiation right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4SynchrotronRadiation::PostStepDoIt ( const G4Track track,
const G4Step Step 
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 200 of file G4SynchrotronRadiation.cc.

203{
204 aParticleChange.Initialize(trackData);
205
206 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
207
208 G4double gamma = aDynamicParticle->GetTotalEnergy() /
209 (aDynamicParticle->GetDefinition()->GetPDGMass());
210
211 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
212 if(gamma <= 1.0e3 || 0.0 == particleCharge)
213 {
214 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
215 }
216
217 G4ThreeVector FieldValue;
218 const G4Field* pField = nullptr;
219
220 G4bool fieldExertsForce = false;
221 G4FieldManager* fieldMgr =
222 fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
223
224 if(fieldMgr != nullptr)
225 {
226 // If the field manager has no field, there is no field !
227 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
228 }
229
230 if(fieldExertsForce)
231 {
232 pField = fieldMgr->GetDetectorField();
233 G4ThreeVector globPosition = trackData.GetPosition();
234 G4double globPosVec[4], FieldValueVec[6];
235 globPosVec[0] = globPosition.x();
236 globPosVec[1] = globPosition.y();
237 globPosVec[2] = globPosition.z();
238 globPosVec[3] = trackData.GetGlobalTime();
239
240 pField->GetFieldValue(globPosVec, FieldValueVec);
241 FieldValue =
242 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
243
244 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
245 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
246 G4double perpB = unitMcrossB.mag();
247 if(perpB > 0.0)
248 {
249 // M-C of synchrotron photon energy
250 G4double energyOfSR = GetRandomEnergySR(
251 gamma, perpB, aDynamicParticle->GetDefinition()->GetPDGMass());
252
253 // check against insufficient energy
254 if(energyOfSR <= 0.0)
255 {
256 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
257 }
258 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
259 G4ThreeVector gammaDirection =
260 genAngle->SampleDirection(aDynamicParticle, energyOfSR, 1, nullptr);
261
262 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
263 gammaPolarization = gammaPolarization.unit();
264
265 // create G4DynamicParticle object for the SR photon
266 auto aGamma =
267 new G4DynamicParticle(theGamma, gammaDirection, energyOfSR);
268 aGamma->SetPolarization(gammaPolarization.x(), gammaPolarization.y(),
269 gammaPolarization.z());
270
272
273 // Update the incident particle
274 G4double newKinEnergy = kineticEnergy - energyOfSR;
275
276 if(newKinEnergy > 0.)
277 {
278 aParticleChange.ProposeEnergy(newKinEnergy);
279 }
280 else
281 {
283 }
284
285 // Create the G4Track
286 G4Track* aSecondaryTrack = new G4Track(aGamma, trackData.GetGlobalTime(), trackData.GetPosition());
287 aSecondaryTrack->SetTouchableHandle(stepData.GetPostStepPoint()->GetTouchableHandle());
288 aSecondaryTrack->SetParentID(trackData.GetTrackID());
289 aSecondaryTrack->SetCreatorModelID(secID);
290 aParticleChange.AddSecondary(aSecondaryTrack);
291
292 }
293 }
294 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
295}
Hep3Vector unit() const
G4double GetKineticEnergy() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
G4double GetRandomEnergySR(G4double, G4double, G4double)
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

◆ ProcessDescription()

void G4SynchrotronRadiation::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 427 of file G4SynchrotronRadiation.cc.

428{
429 out << GetProcessName()
430 << ": Incoherent Synchrotron Radiation\n"
431 "Good description for long magnets at all energies.\n";
432}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

Referenced by BuildPhysicsTable(), and DumpInfo().

◆ SetAngularGenerator()

void G4SynchrotronRadiation::SetAngularGenerator ( G4VEmAngularDistribution p)

Definition at line 85 of file G4SynchrotronRadiation.cc.

86{
87 if(p != genAngle)
88 {
89 delete genAngle;
90 genAngle = p;
91 }
92}

Referenced by G4SynchrotronRadiation().


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