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

#include <G4SynchrotronRadiationInMat.hh>

+ Inheritance diagram for G4SynchrotronRadiationInMat:

Public Member Functions

 G4SynchrotronRadiationInMat (const G4String &processName="SynchrotronRadiation", G4ProcessType type=fElectromagnetic)
 
 ~G4SynchrotronRadiationInMat ()
 
G4SynchrotronRadiationInMatoperator= (const G4SynchrotronRadiationInMat &right)=delete
 
 G4SynchrotronRadiationInMat (const G4SynchrotronRadiationInMat &)=delete
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step) override
 
G4double GetPhotonEnergy (const G4Track &trackData, const G4Step &stepData)
 
G4double GetRandomEnergySR (G4double, G4double)
 
G4double GetProbSpectrumSRforInt (G4double)
 
G4double GetIntProbSR (G4double)
 
G4double GetProbSpectrumSRforEnergy (G4double)
 
G4double GetEnergyProbSR (G4double)
 
G4double GetIntegrandForAngleK (G4double)
 
G4double GetAngleK (G4double)
 
G4double GetAngleNumberAtGammaKsi (G4double)
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void SetRootNumber (G4int rn)
 
void SetVerboseLevel (G4int v)
 
void SetKsi (G4double ksi)
 
void SetEta (G4double eta)
 
void SetPsiGamma (G4double psg)
 
void SetOrderAngleK (G4double ord)
 
- 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 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
 
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 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 &)
 

Static Public Member Functions

static G4double GetLambdaConst ()
 
static G4double GetEnergyConst ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VDiscreteProcess
- 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 55 of file G4SynchrotronRadiationInMat.hh.

Constructor & Destructor Documentation

◆ G4SynchrotronRadiationInMat() [1/2]

G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat ( const G4String & processName = "SynchrotronRadiation",
G4ProcessType type = fElectromagnetic )
explicit

Definition at line 93 of file G4SynchrotronRadiationInMat.cc.

95 : G4VDiscreteProcess(processName, type)
96 , theGamma(G4Gamma::Gamma())
97 , theElectron(G4Electron::Electron())
98 , thePositron(G4Positron::Positron())
99 , LowestKineticEnergy(10. * keV)
100 , fAlpha(0.0)
101 , fRootNumber(80)
102 , fVerboseLevel(verboseLevel)
103{
104 G4TransportationManager* transportMgr =
106
107 fFieldPropagator = transportMgr->GetPropagatorInField();
108 secID = G4PhysicsModelCatalog::GetModelID("model_SynchrotronRadiation");
110 CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
111 PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
112 fPsiGamma = fEta = fOrderAngleK = 0.0;
113}
@ fSynchrotronRadiation
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4int GetModelID(const G4int modelIndex)
static G4Positron * Positron()
Definition G4Positron.cc:90
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4int verboseLevel
void SetProcessSubType(G4int)

Referenced by GetAngleK(), GetEnergyProbSR(), and GetIntProbSR().

◆ ~G4SynchrotronRadiationInMat()

G4SynchrotronRadiationInMat::~G4SynchrotronRadiationInMat ( )
default

◆ G4SynchrotronRadiationInMat() [2/2]

G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat ( const G4SynchrotronRadiationInMat & )
delete

Member Function Documentation

◆ GetAngleK()

G4double G4SynchrotronRadiationInMat::GetAngleK ( G4double eta)

Definition at line 550 of file G4SynchrotronRadiationInMat.cc.

551{
552 fEta = eta; // should be > 0. !
553 G4int n;
554 G4double result, a;
555
556 a = fAlpha; // always = 0.
557 n = fRootNumber; // around default = 80
558
561 integral;
562
563 result = integral.Laguerre(
565
566 return result;
567}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double Laguerre(T &typeT, F f, G4double alpha, G4int n)
G4SynchrotronRadiationInMat(const G4String &processName="SynchrotronRadiation", G4ProcessType type=fElectromagnetic)

Referenced by GetAngleNumberAtGammaKsi().

◆ GetAngleNumberAtGammaKsi()

G4double G4SynchrotronRadiationInMat::GetAngleNumberAtGammaKsi ( G4double gpsi)

Definition at line 571 of file G4SynchrotronRadiationInMat.cc.

572{
573 G4double result, funK, funK2, gpsi2 = gpsi * gpsi;
574
575 fPsiGamma = gpsi;
576 fEta = 0.5 * fKsi * (1. + gpsi2) * std::sqrt(1. + gpsi2);
577
578 fOrderAngleK = 1. / 3.;
579 funK = GetAngleK(fEta);
580 funK2 = funK * funK;
581
582 result = gpsi2 * funK2 / (1. + gpsi2);
583
584 fOrderAngleK = 2. / 3.;
585 funK = GetAngleK(fEta);
586 funK2 = funK * funK;
587
588 result += funK2;
589 result *= (1. + gpsi2) * fKsi;
590
591 return result;
592}

◆ GetEnergyConst()

G4double G4SynchrotronRadiationInMat::GetEnergyConst ( )
static

Definition at line 128 of file G4SynchrotronRadiationInMat.cc.

128{ return fEnergyConst; }

◆ GetEnergyProbSR()

G4double G4SynchrotronRadiationInMat::GetEnergyProbSR ( G4double ksi)

Definition at line 514 of file G4SynchrotronRadiationInMat.cc.

515{
516 if(ksi <= 0.)
517 return 1.0;
518 fKsi = ksi; // should be > 0. !
519 G4int n;
520 G4double result, a;
521
522 a = fAlpha; // always = 0.
523 n = fRootNumber; // around default = 80
524
527 integral;
528
529 result = integral.Laguerre(
531
532 result *= 9. * std::sqrt(3.) * ksi / 8. / pi;
533
534 return result;
535}

◆ GetIntegrandForAngleK()

G4double G4SynchrotronRadiationInMat::GetIntegrandForAngleK ( G4double t)

Definition at line 538 of file G4SynchrotronRadiationInMat.cc.

539{
540 G4double result, hypCos = std::cosh(t);
541
542 result =
543 std::cosh(fOrderAngleK * t) * std::exp(t - fEta * hypCos); // fEta > 0. !
544 result /= hypCos;
545 return result;
546}

Referenced by GetAngleK().

◆ GetIntProbSR()

G4double G4SynchrotronRadiationInMat::GetIntProbSR ( G4double ksi)

Definition at line 476 of file G4SynchrotronRadiationInMat.cc.

477{
478 if(ksi <= 0.)
479 return 1.0;
480 fKsi = ksi; // should be > 0. !
481 G4int n;
482 G4double result, a;
483
484 a = fAlpha; // always = 0.
485 n = fRootNumber; // around default = 80
486
489 integral;
490
491 result = integral.Laguerre(
493
494 result *= 3. / 5. / pi;
495
496 return result;
497}

◆ GetLambdaConst()

G4double G4SynchrotronRadiationInMat::GetLambdaConst ( )
static

Definition at line 126 of file G4SynchrotronRadiationInMat.cc.

126{ return fLambdaConst; }

◆ GetMeanFreePath()

G4double G4SynchrotronRadiationInMat::GetMeanFreePath ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 132 of file G4SynchrotronRadiationInMat.cc.

134{
135 // gives the MeanFreePath in GEANT4 internal units
136 G4double MeanFreePath;
137
138 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
139
141
142 G4double gamma =
143 aDynamicParticle->GetTotalEnergy() / aDynamicParticle->GetMass();
144
145 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
146
147 G4double KineticEnergy = aDynamicParticle->GetKineticEnergy();
148
149 if(KineticEnergy < LowestKineticEnergy || gamma < 1.0e3)
150 MeanFreePath = DBL_MAX;
151 else
152 {
153 G4ThreeVector FieldValue;
154 const G4Field* pField = nullptr;
155
156 G4FieldManager* fieldMgr = nullptr;
157 G4bool fieldExertsForce = false;
158
159 if((particleCharge != 0.0))
160 {
161 fieldMgr =
162 fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
163
164 if(fieldMgr != nullptr)
165 {
166 // If the field manager has no field, there is no field !
167 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
168 }
169 }
170 if(fieldExertsForce)
171 {
172 pField = fieldMgr->GetDetectorField();
173 G4ThreeVector globPosition = trackData.GetPosition();
174
175 G4double globPosVec[4], FieldValueVec[6];
176
177 globPosVec[0] = globPosition.x();
178 globPosVec[1] = globPosition.y();
179 globPosVec[2] = globPosition.z();
180 globPosVec[3] = trackData.GetGlobalTime();
181
182 pField->GetFieldValue(globPosVec, FieldValueVec);
183
184 FieldValue =
185 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
186
187 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
188 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
189 G4double perpB = unitMcrossB.mag();
190 G4double beta = aDynamicParticle->GetTotalMomentum() /
191 (aDynamicParticle->GetTotalEnergy());
192
193 if(perpB > 0.0)
194 MeanFreePath = fLambdaConst * beta / perpB;
195 else
196 MeanFreePath = DBL_MAX;
197 }
198 else
199 MeanFreePath = DBL_MAX;
200 }
201 if(fVerboseLevel > 0)
202 {
203 G4cout << "G4SynchrotronRadiationInMat::MeanFreePath = " << MeanFreePath / m
204 << " m" << G4endl;
205 }
206 return MeanFreePath;
207}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double mag() const
G4double GetMass() const
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
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
#define DBL_MAX
Definition templates.hh:62

◆ GetPhotonEnergy()

G4double G4SynchrotronRadiationInMat::GetPhotonEnergy ( const G4Track & trackData,
const G4Step & stepData )

Definition at line 360 of file G4SynchrotronRadiationInMat.cc.

362{
363 G4int i;
364 G4double energyOfSR = -1.0;
365
366 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
367
368 G4double gamma =
369 aDynamicParticle->GetTotalEnergy() / (aDynamicParticle->GetMass());
370
371 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
372
373 G4ThreeVector FieldValue;
374 const G4Field* pField = nullptr;
375
376 G4FieldManager* fieldMgr = nullptr;
377 G4bool fieldExertsForce = false;
378
379 if((particleCharge != 0.0))
380 {
381 fieldMgr = fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
382 if(fieldMgr != nullptr)
383 {
384 // If the field manager has no field, there is no field !
385 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
386 }
387 }
388 if(fieldExertsForce)
389 {
390 pField = fieldMgr->GetDetectorField();
391 G4ThreeVector globPosition = trackData.GetPosition();
392 G4double globPosVec[3], FieldValueVec[3];
393
394 globPosVec[0] = globPosition.x();
395 globPosVec[1] = globPosition.y();
396 globPosVec[2] = globPosition.z();
397
398 pField->GetFieldValue(globPosVec, FieldValueVec);
399 FieldValue =
400 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
401
402 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
403 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
404 G4double perpB = unitMcrossB.mag();
405 if(perpB > 0.0)
406 {
407 // M-C of synchrotron photon energy
408 G4double random = G4UniformRand();
409 for(i = 0; i < 200; ++i)
410 {
411 if(random >= fIntegralProbabilityOfSR[i])
412 break;
413 }
414 energyOfSR = 0.0001 * i * i * fEnergyConst * gamma * gamma * perpB;
415
416 // check against insufficient energy
417 if(energyOfSR <= 0.0)
418 {
419 return -1.0;
420 }
421 }
422 else
423 {
424 return -1.0;
425 }
426 }
427 return energyOfSR;
428}
#define G4UniformRand()
Definition Randomize.hh:52
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
const G4DynamicParticle * GetDynamicParticle() const

◆ GetProbSpectrumSRforEnergy()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforEnergy ( G4double t)

Definition at line 501 of file G4SynchrotronRadiationInMat.cc.

502{
503 G4double result, hypCos = std::cosh(t);
504
505 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos); // fKsi > 0. !
506 result /= hypCos;
507 return result;
508}

Referenced by GetEnergyProbSR().

◆ GetProbSpectrumSRforInt()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforInt ( G4double t)

Definition at line 462 of file G4SynchrotronRadiationInMat.cc.

463{
464 G4double result, hypCos2, hypCos = std::cosh(t);
465
466 hypCos2 = hypCos * hypCos;
467 result = std::cosh(5. * t / 3.) * std::exp(t - fKsi * hypCos); // fKsi > 0. !
468 result /= hypCos2;
469 return result;
470}

Referenced by GetIntProbSR().

◆ GetRandomEnergySR()

G4double G4SynchrotronRadiationInMat::GetRandomEnergySR ( G4double gamma,
G4double perpB )

Definition at line 431 of file G4SynchrotronRadiationInMat.cc.

433{
434 G4int i;
435 static constexpr G4int iMax = 200;
436 G4double energySR, random, position;
437
438 random = G4UniformRand();
439
440 for(i = 0; i < iMax; ++i)
441 {
442 if(random >= fIntegralProbabilityOfSR[i])
443 break;
444 }
445 if(i <= 0)
447 else if(i >= iMax)
448 position = G4double(iMax);
449 else
450 position = i + G4UniformRand();
451
452 energySR =
453 0.0001 * position * position * fEnergyConst * gamma * gamma * perpB;
454
455 if(energySR < 0.)
456 energySR = 0.;
457
458 return energySR;
459}

Referenced by PostStepDoIt().

◆ IsApplicable()

G4bool G4SynchrotronRadiationInMat::IsApplicable ( const G4ParticleDefinition & particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 119 of file G4SynchrotronRadiationInMat.cc.

121{
122 return ((&particle == (const G4ParticleDefinition*) theElectron) ||
123 (&particle == (const G4ParticleDefinition*) thePositron));
124}

◆ operator=()

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 210 of file G4SynchrotronRadiationInMat.cc.

213{
214 aParticleChange.Initialize(trackData);
215
216 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
217
218 G4double gamma =
219 aDynamicParticle->GetTotalEnergy() / (aDynamicParticle->GetMass());
220
221 if(gamma <= 1.0e3)
222 {
223 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
224 }
225 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
226
227 G4ThreeVector FieldValue;
228 const G4Field* pField = nullptr;
229
230 G4FieldManager* fieldMgr = nullptr;
231 G4bool fieldExertsForce = false;
232
233 if((particleCharge != 0.0))
234 {
235 fieldMgr = fFieldPropagator->FindAndSetFieldManager(trackData.GetVolume());
236 if(fieldMgr != nullptr)
237 {
238 // If the field manager has no field, there is no field !
239 fieldExertsForce = (fieldMgr->GetDetectorField() != nullptr);
240 }
241 }
242 if(fieldExertsForce)
243 {
244 pField = fieldMgr->GetDetectorField();
245 G4ThreeVector globPosition = trackData.GetPosition();
246 G4double globPosVec[4], FieldValueVec[6];
247 globPosVec[0] = globPosition.x();
248 globPosVec[1] = globPosition.y();
249 globPosVec[2] = globPosition.z();
250 globPosVec[3] = trackData.GetGlobalTime();
251
252 pField->GetFieldValue(globPosVec, FieldValueVec);
253 FieldValue =
254 G4ThreeVector(FieldValueVec[0], FieldValueVec[1], FieldValueVec[2]);
255
256 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
257 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
258 G4double perpB = unitMcrossB.mag();
259 if(perpB > 0.0)
260 {
261 // M-C of synchrotron photon energy
262 G4double energyOfSR = GetRandomEnergySR(gamma, perpB);
263
264 if(fVerboseLevel > 0)
265 {
266 G4cout << "SR photon energy = " << energyOfSR / keV << " keV" << G4endl;
267 }
268
269 // check against insufficient energy
270 if(energyOfSR <= 0.0)
271 {
272 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
273 }
274 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
275 G4ParticleMomentum particleDirection =
276 aDynamicParticle->GetMomentumDirection();
277
278 // M-C of its direction, simplified dipole busted approach
279 G4double cosTheta, sinTheta, fcos, beta;
280
281 do
282 {
283 cosTheta = 1. - 2. * G4UniformRand();
284 fcos = (1 + cosTheta * cosTheta) * 0.5;
285 }
286 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
287 while(fcos < G4UniformRand());
288
289 beta = std::sqrt(1. - 1. / (gamma * gamma));
290
291 cosTheta = (cosTheta + beta) / (1. + beta * cosTheta);
292
293 if(cosTheta > 1.)
294 cosTheta = 1.;
295 if(cosTheta < -1.)
296 cosTheta = -1.;
297
298 sinTheta = std::sqrt(1. - cosTheta * cosTheta);
299
300 G4double Phi = twopi * G4UniformRand();
301
302 G4double dirx = sinTheta * std::cos(Phi);
303 G4double diry = sinTheta * std::sin(Phi);
304 G4double dirz = cosTheta;
305
306 G4ThreeVector gammaDirection(dirx, diry, dirz);
307 gammaDirection.rotateUz(particleDirection);
308
309 // polarization of new gamma
310 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
311 gammaPolarization = gammaPolarization.unit();
312
313 // create G4DynamicParticle object for the SR photon
314 auto aGamma =
315 new G4DynamicParticle(G4Gamma::Gamma(), gammaDirection, energyOfSR);
316 aGamma->SetPolarization(gammaPolarization.x(), gammaPolarization.y(),
317 gammaPolarization.z());
318
320
321 // Update the incident particle
322 G4double newKinEnergy = kineticEnergy - energyOfSR;
323
324 if(newKinEnergy > 0.)
325 {
326 aParticleChange.ProposeMomentumDirection(particleDirection);
327 aParticleChange.ProposeEnergy(newKinEnergy);
329 }
330 else
331 {
334 G4double charge = aDynamicParticle->GetDefinition()->GetPDGCharge();
335 if(charge < 0.)
336 {
338 }
339 else
340 {
342 }
343 }
344
345 // Create the G4Track
346 G4Track* aSecondaryTrack = new G4Track(aGamma, trackData.GetGlobalTime(), trackData.GetPosition());
347 aSecondaryTrack->SetTouchableHandle(stepData.GetPostStepPoint()->GetTouchableHandle());
348 aSecondaryTrack->SetParentID(trackData.GetTrackID());
349 aSecondaryTrack->SetCreatorModelID(secID);
350 aParticleChange.AddSecondary(aSecondaryTrack);
351 }
352 else
353 {
354 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
355 }
356 }
357 return G4VDiscreteProcess::PostStepDoIt(trackData, stepData);
358}
@ fStopAndKill
@ fStopButAlive
Hep3Vector unit() const
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetRandomEnergySR(G4double, G4double)
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetCreatorModelID(const G4int id)
void SetParentID(const G4int aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange

◆ SetEta()

void G4SynchrotronRadiationInMat::SetEta ( G4double eta)
inline

Definition at line 96 of file G4SynchrotronRadiationInMat.hh.

96{ fEta = eta; };

◆ SetKsi()

void G4SynchrotronRadiationInMat::SetKsi ( G4double ksi)
inline

Definition at line 95 of file G4SynchrotronRadiationInMat.hh.

95{ fKsi = ksi; };

◆ SetOrderAngleK()

void G4SynchrotronRadiationInMat::SetOrderAngleK ( G4double ord)
inline

Definition at line 98 of file G4SynchrotronRadiationInMat.hh.

99 {
100 fOrderAngleK = ord;
101 }; // should be 1/3 or 2/3

◆ SetPsiGamma()

void G4SynchrotronRadiationInMat::SetPsiGamma ( G4double psg)
inline

Definition at line 97 of file G4SynchrotronRadiationInMat.hh.

97{ fPsiGamma = psg; };

◆ SetRootNumber()

void G4SynchrotronRadiationInMat::SetRootNumber ( G4int rn)
inline

Definition at line 93 of file G4SynchrotronRadiationInMat.hh.

93{ fRootNumber = rn; };

◆ SetVerboseLevel()

void G4SynchrotronRadiationInMat::SetVerboseLevel ( G4int v)
inline

Definition at line 94 of file G4SynchrotronRadiationInMat.hh.

94{ fVerboseLevel = v; };

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