Geant4 10.7.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) 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 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 ()
 
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 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

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 67 of file G4SynchrotronRadiationInMat.hh.

Constructor & Destructor Documentation

◆ G4SynchrotronRadiationInMat()

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

Definition at line 122 of file G4SynchrotronRadiationInMat.cc.

123 :G4VDiscreteProcess (processName, type),
124 LowestKineticEnergy (10.*keV),
125 theGamma (G4Gamma::Gamma() ),
126 theElectron ( G4Electron::Electron() ),
127 thePositron ( G4Positron::Positron() ),
128 fAlpha(0.0), fRootNumber(80),
129 fVerboseLevel( verboseLevel )
130{
132
133 fFieldPropagator = transportMgr->GetPropagatorInField();
135 CutInRange = GammaCutInKineticEnergyNow = ElectronCutInKineticEnergyNow =
136 PositronCutInKineticEnergyNow = ParticleCutInKineticEnergyNow = fKsi =
137 fPsiGamma = fEta = fOrderAngleK = 0.0;
138}
@ fSynchrotronRadiation
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4SynchrotronRadiationInMat()

G4SynchrotronRadiationInMat::~G4SynchrotronRadiationInMat ( )
virtual

Definition at line 145 of file G4SynchrotronRadiationInMat.cc.

146{}

Member Function Documentation

◆ GetAngleK()

G4double G4SynchrotronRadiationInMat::GetAngleK ( G4double  eta)

Definition at line 638 of file G4SynchrotronRadiationInMat.cc.

639{
640 fEta = eta; // should be > 0. !
641 G4int n;
642 G4double result, a;
643
644 a = fAlpha; // always = 0.
645 n = fRootNumber; // around default = 80
646
648
649 result = integral.Laguerre(this,
651
652 return result;
653}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85

Referenced by GetAngleNumberAtGammaKsi().

◆ GetAngleNumberAtGammaKsi()

G4double G4SynchrotronRadiationInMat::GetAngleNumberAtGammaKsi ( G4double  gpsi)

Definition at line 659 of file G4SynchrotronRadiationInMat.cc.

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

◆ GetEnergyConst()

G4double G4SynchrotronRadiationInMat::GetEnergyConst ( )
static

Definition at line 161 of file G4SynchrotronRadiationInMat.cc.

162{
163 return fEnergyConst;
164}

◆ GetEnergyProbSR()

G4double G4SynchrotronRadiationInMat::GetEnergyProbSR ( G4double  ksi)

Definition at line 601 of file G4SynchrotronRadiationInMat.cc.

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

◆ GetIntegrandForAngleK()

G4double G4SynchrotronRadiationInMat::GetIntegrandForAngleK ( G4double  t)

Definition at line 625 of file G4SynchrotronRadiationInMat.cc.

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

Referenced by GetAngleK().

◆ GetIntProbSR()

G4double G4SynchrotronRadiationInMat::GetIntProbSR ( G4double  ksi)

Definition at line 562 of file G4SynchrotronRadiationInMat.cc.

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

◆ GetLambdaConst()

G4double G4SynchrotronRadiationInMat::GetLambdaConst ( )
static

Definition at line 156 of file G4SynchrotronRadiationInMat.cc.

157{
158 return fLambdaConst;
159}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 175 of file G4SynchrotronRadiationInMat.cc.

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

◆ GetPhotonEnergy()

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

Definition at line 423 of file G4SynchrotronRadiationInMat.cc.

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

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

Referenced by GetEnergyProbSR().

◆ GetProbSpectrumSRforInt()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforInt ( G4double  t)

Definition at line 546 of file G4SynchrotronRadiationInMat.cc.

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

Referenced by GetIntProbSR().

◆ GetRandomEnergySR()

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

Definition at line 516 of file G4SynchrotronRadiationInMat.cc.

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

Referenced by PostStepDoIt().

◆ IsApplicable()

G4bool G4SynchrotronRadiationInMat::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 150 of file G4SynchrotronRadiationInMat.cc.

151{
152 return ( ( &particle == (const G4ParticleDefinition *)theElectron ) ||
153 ( &particle == (const G4ParticleDefinition *)thePositron ));
154}

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 260 of file G4SynchrotronRadiationInMat.cc.

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

◆ 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: