Geant4 9.6.0
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)
 
virtual ~G4SynchrotronRadiationInMat ()
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)
 
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 &)
 
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 &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Static Public Member Functions

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

Additional Inherited Members

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

Detailed Description

Definition at line 68 of file G4SynchrotronRadiationInMat.hh.

Constructor & Destructor Documentation

◆ G4SynchrotronRadiationInMat()

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

Definition at line 123 of file G4SynchrotronRadiationInMat.cc.

124 :G4VDiscreteProcess (processName, type),
125 LowestKineticEnergy (10.*keV),
126 HighestKineticEnergy (100.*TeV),
127 TotBin(200),
128 theGamma (G4Gamma::Gamma() ),
129 theElectron ( G4Electron::Electron() ),
130 thePositron ( G4Positron::Positron() ),
131 GammaCutInKineticEnergy(0),
132 ElectronCutInKineticEnergy(0),
133 PositronCutInKineticEnergy(0),
134 ParticleCutInKineticEnergy(0),
135 fAlpha(0.0), fRootNumber(80),
136 fVerboseLevel( verboseLevel )
137{
139
140 fFieldPropagator = transportMgr->GetPropagatorInField();
142 CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
143 PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
144 fPsiGamma = fEta = fOrderAngleK = 0.0;
145}
@ fSynchrotronRadiation
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4SynchrotronRadiationInMat()

G4SynchrotronRadiationInMat::~G4SynchrotronRadiationInMat ( )
virtual

Definition at line 152 of file G4SynchrotronRadiationInMat.cc.

153{}

Member Function Documentation

◆ GetAngleK()

G4double G4SynchrotronRadiationInMat::GetAngleK ( G4double  eta)

Definition at line 636 of file G4SynchrotronRadiationInMat.cc.

637{
638 fEta = eta; // should be > 0. !
639 G4int n;
640 G4double result, a;
641
642 a = fAlpha; // always = 0.
643 n = fRootNumber; // around default = 80
644
646
647 result = integral.Laguerre(this,
649
650 return result;
651}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66

Referenced by GetAngleNumberAtGammaKsi().

◆ GetAngleNumberAtGammaKsi()

G4double G4SynchrotronRadiationInMat::GetAngleNumberAtGammaKsi ( G4double  gpsi)

Definition at line 657 of file G4SynchrotronRadiationInMat.cc.

658{
659 G4double result, funK, funK2, gpsi2 = gpsi*gpsi;
660
661 fPsiGamma = gpsi;
662 fEta = 0.5*fKsi*(1 + gpsi2)*std::sqrt(1 + gpsi2);
663
664 fOrderAngleK = 1./3.;
665 funK = GetAngleK(fEta);
666 funK2 = funK*funK;
667
668 result = gpsi2*funK2/(1 + gpsi2);
669
670 fOrderAngleK = 2./3.;
671 funK = GetAngleK(fEta);
672 funK2 = funK*funK;
673
674 result += funK2;
675 result *= (1 + gpsi2)*fKsi;
676
677 return result;
678}

◆ GetEnergyConst()

static G4double G4SynchrotronRadiationInMat::GetEnergyConst ( )
inlinestatic

Definition at line 110 of file G4SynchrotronRadiationInMat.hh.

110{ return fEnergyConst; };

◆ GetEnergyProbSR()

G4double G4SynchrotronRadiationInMat::GetEnergyProbSR ( G4double  ksi)

Definition at line 599 of file G4SynchrotronRadiationInMat.cc.

600{
601 if (ksi <= 0.) return 1.0;
602 fKsi = ksi; // should be > 0. !
603 G4int n;
604 G4double result, a;
605
606 a = fAlpha; // always = 0.
607 n = fRootNumber; // around default = 80
608
610
611 result = integral.Laguerre(this,
613
614 result *= 9.*std::sqrt(3.)*ksi/8./pi;
615
616 return result;
617}
const G4double pi

◆ GetIntegrandForAngleK()

G4double G4SynchrotronRadiationInMat::GetIntegrandForAngleK ( G4double  t)

Definition at line 623 of file G4SynchrotronRadiationInMat.cc.

624{
625 G4double result, hypCos=std::cosh(t);
626
627 result = std::cosh(fOrderAngleK*t)*std::exp(t - fEta*hypCos); // fEta > 0. !
628 result /= hypCos;
629 return result;
630}

Referenced by GetAngleK().

◆ GetIntProbSR()

G4double G4SynchrotronRadiationInMat::GetIntProbSR ( G4double  ksi)

Definition at line 560 of file G4SynchrotronRadiationInMat.cc.

561{
562 if (ksi <= 0.) return 1.0;
563 fKsi = ksi; // should be > 0. !
564 G4int n;
565 G4double result, a;
566
567 a = fAlpha; // always = 0.
568 n = fRootNumber; // around default = 80
569
571
572 result = integral.Laguerre(this,
574
575 result *= 3./5./pi;
576
577 return result;
578}

◆ GetLambdaConst()

static G4double G4SynchrotronRadiationInMat::GetLambdaConst ( )
inlinestatic

Definition at line 109 of file G4SynchrotronRadiationInMat.hh.

109{ return fLambdaConst; };

◆ GetMeanFreePath()

G4double G4SynchrotronRadiationInMat::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 174 of file G4SynchrotronRadiationInMat.cc.

177{
178 // gives the MeanFreePath in GEANT4 internal units
179 G4double MeanFreePath;
180
181 const G4DynamicParticle* aDynamicParticle = trackData.GetDynamicParticle();
182 // G4Material* aMaterial = trackData.GetMaterial();
183
184 //G4bool isOutRange ;
185
187
188 G4double gamma = aDynamicParticle->GetTotalEnergy()/
189 aDynamicParticle->GetMass();
190
191 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
192
193 G4double KineticEnergy = aDynamicParticle->GetKineticEnergy();
194
195 if ( KineticEnergy < LowestKineticEnergy || gamma < 1.0e3 ) MeanFreePath = DBL_MAX;
196 else
197 {
198
199 G4ThreeVector FieldValue;
200 const G4Field* pField = 0;
201
202 G4FieldManager* fieldMgr=0;
203 G4bool fieldExertsForce = false;
204
205 if( (particleCharge != 0.0) )
206 {
207 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
208
209 if ( fieldMgr != 0 )
210 {
211 // If the field manager has no field, there is no field !
212
213 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
214 }
215 }
216 if ( fieldExertsForce )
217 {
218 pField = fieldMgr->GetDetectorField() ;
219 G4ThreeVector globPosition = trackData.GetPosition();
220
221 G4double globPosVec[4], FieldValueVec[6];
222
223 globPosVec[0] = globPosition.x();
224 globPosVec[1] = globPosition.y();
225 globPosVec[2] = globPosition.z();
226 globPosVec[3] = trackData.GetGlobalTime();
227
228 pField->GetFieldValue( globPosVec, FieldValueVec );
229
230 FieldValue = G4ThreeVector( FieldValueVec[0],
231 FieldValueVec[1],
232 FieldValueVec[2] );
233
234
235
236 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
237 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum) ;
238 G4double perpB = unitMcrossB.mag() ;
239 G4double beta = aDynamicParticle->GetTotalMomentum()/
240 (aDynamicParticle->GetTotalEnergy() );
241
242 if( perpB > 0.0 ) MeanFreePath = fLambdaConst*beta/perpB;
243 else MeanFreePath = DBL_MAX;
244 }
245 else MeanFreePath = DBL_MAX;
246 }
247 if(fVerboseLevel > 0)
248 {
249 G4cout<<"G4SynchrotronRadiationInMat::MeanFreePath = "<<MeanFreePath/m<<" m"<<G4endl;
250 }
251 return MeanFreePath;
252}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
CLHEP::Hep3Vector G4ThreeVector
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 double Point[4], double *fieldArr) const =0
G4double GetPDGCharge() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
#define DBL_MAX
Definition: templates.hh:83

◆ GetPhotonEnergy()

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

Definition at line 421 of file G4SynchrotronRadiationInMat.cc.

424{
425 G4int i ;
426 G4double energyOfSR = -1.0 ;
427 //G4Material* aMaterial=trackData.GetMaterial() ;
428
429 const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
430
431 G4double gamma = aDynamicParticle->GetTotalEnergy()/
432 (aDynamicParticle->GetMass() ) ;
433
434 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
435
436 G4ThreeVector FieldValue;
437 const G4Field* pField = 0 ;
438
439 G4FieldManager* fieldMgr=0;
440 G4bool fieldExertsForce = false;
441
442 if( (particleCharge != 0.0) )
443 {
444 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
445 if ( fieldMgr != 0 )
446 {
447 // If the field manager has no field, there is no field !
448
449 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
450 }
451 }
452 if ( fieldExertsForce )
453 {
454 pField = fieldMgr->GetDetectorField();
455 G4ThreeVector globPosition = trackData.GetPosition();
456 G4double globPosVec[3], FieldValueVec[3];
457
458 globPosVec[0] = globPosition.x();
459 globPosVec[1] = globPosition.y();
460 globPosVec[2] = globPosition.z();
461
462 pField->GetFieldValue( globPosVec, FieldValueVec );
463 FieldValue = G4ThreeVector( FieldValueVec[0],
464 FieldValueVec[1],
465 FieldValueVec[2] );
466
467 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
468 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum) ;
469 G4double perpB = unitMcrossB.mag();
470 if( perpB > 0.0 )
471 {
472 // M-C of synchrotron photon energy
473
474 G4double random = G4UniformRand() ;
475 for(i=0;i<200;i++)
476 {
477 if(random >= fIntegralProbabilityOfSR[i]) break ;
478 }
479 energyOfSR = 0.0001*i*i*fEnergyConst*gamma*gamma*perpB ;
480
481 // check against insufficient energy
482
483 if(energyOfSR <= 0.0)
484 {
485 return -1.0 ;
486 }
487 //G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
488 //G4ParticleMomentum
489 //particleDirection = aDynamicParticle->GetMomentumDirection();
490
491 // Gamma production cut in this material
492 //G4double
493 //gammaEnergyCut = (G4Gamma::GetCutsInEnergy())[aMaterial->GetIndex()];
494
495 // SR photon has energy more than the current material cut
496 // M-C of its direction
497
498 //G4double Teta = G4UniformRand()/gamma ; // Very roughly
499
500 //G4double Phi = twopi * G4UniformRand() ;
501 }
502 else
503 {
504 return -1.0 ;
505 }
506 }
507 return energyOfSR ;
508}
#define G4UniformRand()
Definition: Randomize.hh:53
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
const G4DynamicParticle * GetDynamicParticle() const

◆ GetProbSpectrumSRforEnergy()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforEnergy ( G4double  t)

Definition at line 584 of file G4SynchrotronRadiationInMat.cc.

585{
586 G4double result, hypCos=std::cosh(t);
587
588 result = std::cosh(5.*t/3.)*std::exp(t - fKsi*hypCos); // fKsi > 0. !
589 result /= hypCos;
590 return result;
591}

Referenced by GetEnergyProbSR().

◆ GetProbSpectrumSRforInt()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforInt ( G4double  t)

Definition at line 544 of file G4SynchrotronRadiationInMat.cc.

545{
546 G4double result, hypCos2, hypCos=std::cosh(t);
547
548 hypCos2 = hypCos*hypCos;
549 result = std::cosh(5.*t/3.)*std::exp(t-fKsi*hypCos); // fKsi > 0. !
550 result /= hypCos2;
551 return result;
552}

Referenced by GetIntProbSR().

◆ GetRandomEnergySR()

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

Definition at line 514 of file G4SynchrotronRadiationInMat.cc.

515{
516 G4int i, iMax;
517 G4double energySR, random, position;
518
519 iMax = 200;
520 random = G4UniformRand();
521
522 for( i = 0; i < iMax; i++ )
523 {
524 if( random >= fIntegralProbabilityOfSR[i] ) break;
525 }
526 if(i <= 0 ) position = G4UniformRand(); // 0.
527 else if( i>= iMax) position = G4double(iMax);
528 else position = i + G4UniformRand(); // -1
529 //
530 // it was in initial implementation:
531 // energyOfSR = 0.0001*i*i*fEnergyConst*gamma*gamma*perpB ;
532
533 energySR = 0.0001*position*position*fEnergyConst*gamma*gamma*perpB;
534
535 if( energySR < 0. ) energySR = 0.;
536
537 return energySR;
538}
#define position
Definition: xmlparse.cc:605

Referenced by PostStepDoIt().

◆ IsApplicable()

G4bool G4SynchrotronRadiationInMat::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 157 of file G4SynchrotronRadiationInMat.cc.

158{
159
160 return ( ( &particle == (const G4ParticleDefinition *)theElectron ) ||
161 ( &particle == (const G4ParticleDefinition *)thePositron ) );
162
163}

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 259 of file G4SynchrotronRadiationInMat.cc.

262{
263 aParticleChange.Initialize(trackData);
264
265 const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
266
267 G4double gamma = aDynamicParticle->GetTotalEnergy()/
268 (aDynamicParticle->GetMass() );
269
270 if(gamma <= 1.0e3 )
271 {
272 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
273 }
274 G4double particleCharge = aDynamicParticle->GetDefinition()->GetPDGCharge();
275
276 G4ThreeVector FieldValue;
277 const G4Field* pField = 0 ;
278
279 G4FieldManager* fieldMgr=0;
280 G4bool fieldExertsForce = false;
281
282 if( (particleCharge != 0.0) )
283 {
284 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
285 if ( fieldMgr != 0 )
286 {
287 // If the field manager has no field, there is no field !
288
289 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
290 }
291 }
292 if ( fieldExertsForce )
293 {
294 pField = fieldMgr->GetDetectorField() ;
295 G4ThreeVector globPosition = trackData.GetPosition() ;
296 G4double globPosVec[4], FieldValueVec[6] ;
297 globPosVec[0] = globPosition.x() ;
298 globPosVec[1] = globPosition.y() ;
299 globPosVec[2] = globPosition.z() ;
300 globPosVec[3] = trackData.GetGlobalTime();
301
302 pField->GetFieldValue( globPosVec, FieldValueVec ) ;
303 FieldValue = G4ThreeVector( FieldValueVec[0],
304 FieldValueVec[1],
305 FieldValueVec[2] );
306
307 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
308 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
309 G4double perpB = unitMcrossB.mag() ;
310 if(perpB > 0.0)
311 {
312 // M-C of synchrotron photon energy
313
314 G4double energyOfSR = GetRandomEnergySR(gamma,perpB);
315
316 if(fVerboseLevel > 0)
317 {
318 G4cout<<"SR photon energy = "<<energyOfSR/keV<<" keV"<<G4endl;
319 }
320 // check against insufficient energy
321
322 if( energyOfSR <= 0.0 )
323 {
324 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
325 }
326 G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
328 particleDirection = aDynamicParticle->GetMomentumDirection();
329
330 // M-C of its direction, simplified dipole busted approach
331
332 // G4double Teta = G4UniformRand()/gamma ; // Very roughly
333
334 G4double cosTheta, sinTheta, fcos, beta;
335
336 do
337 {
338 cosTheta = 1. - 2.*G4UniformRand();
339 fcos = (1 + cosTheta*cosTheta)*0.5;
340 }
341 while( fcos < G4UniformRand() );
342
343 beta = std::sqrt(1. - 1./(gamma*gamma));
344
345 cosTheta = (cosTheta + beta)/(1. + beta*cosTheta);
346
347 if( cosTheta > 1. ) cosTheta = 1.;
348 if( cosTheta < -1. ) cosTheta = -1.;
349
350 sinTheta = std::sqrt(1. - cosTheta*cosTheta );
351
352 G4double Phi = twopi * G4UniformRand() ;
353
354 G4double dirx = sinTheta*std::cos(Phi) ,
355 diry = sinTheta*std::sin(Phi) ,
356 dirz = cosTheta;
357
358 G4ThreeVector gammaDirection ( dirx, diry, dirz);
359 gammaDirection.rotateUz(particleDirection);
360
361 // polarization of new gamma
362
363 // G4double sx = std::cos(Teta)*std::cos(Phi);
364 // G4double sy = std::cos(Teta)*std::sin(Phi);
365 // G4double sz = -std::sin(Teta);
366
367 G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
368 gammaPolarization = gammaPolarization.unit();
369
370 // (sx, sy, sz);
371 // gammaPolarization.rotateUz(particleDirection);
372
373 // create G4DynamicParticle object for the SR photon
374
376 gammaDirection,
377 energyOfSR );
378 aGamma->SetPolarization( gammaPolarization.x(),
379 gammaPolarization.y(),
380 gammaPolarization.z() );
381
382
385
386 // Update the incident particle
387
388 G4double newKinEnergy = kineticEnergy - energyOfSR ;
389
390 if (newKinEnergy > 0.)
391 {
392 aParticleChange.ProposeMomentumDirection( particleDirection );
393 aParticleChange.ProposeEnergy( newKinEnergy );
395 }
396 else
397 {
400 G4double charge = aDynamicParticle->GetDefinition()->GetPDGCharge();
401 if (charge<0.)
402 {
404 }
405 else
406 {
408 }
409 }
410 }
411 else
412 {
413 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
414 }
415 }
416 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
417}
@ fStopAndKill
@ fStopButAlive
G4double fcos(G4double arg)
Hep3Vector unit() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4double GetRandomEnergySR(G4double, G4double)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

◆ SetEta()

void G4SynchrotronRadiationInMat::SetEta ( G4double  eta)
inline

Definition at line 115 of file G4SynchrotronRadiationInMat.hh.

115{ fEta = eta; };

◆ SetKsi()

void G4SynchrotronRadiationInMat::SetKsi ( G4double  ksi)
inline

Definition at line 114 of file G4SynchrotronRadiationInMat.hh.

114{ fKsi = ksi; };

◆ SetOrderAngleK()

void G4SynchrotronRadiationInMat::SetOrderAngleK ( G4double  ord)
inline

Definition at line 117 of file G4SynchrotronRadiationInMat.hh.

117{ fOrderAngleK = ord; }; // should be 1/3 or 2/3

◆ SetPsiGamma()

void G4SynchrotronRadiationInMat::SetPsiGamma ( G4double  psg)
inline

Definition at line 116 of file G4SynchrotronRadiationInMat.hh.

116{ fPsiGamma = psg; };

◆ SetRootNumber()

void G4SynchrotronRadiationInMat::SetRootNumber ( G4int  rn)
inline

Definition at line 112 of file G4SynchrotronRadiationInMat.hh.

112{ fRootNumber = rn; };

◆ SetVerboseLevel()

void G4SynchrotronRadiationInMat::SetVerboseLevel ( G4int  v)
inline

Definition at line 113 of file G4SynchrotronRadiationInMat.hh.

113{ fVerboseLevel = v; };

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