Geant4  10.02.p03
G4SynchrotronRadiationInMat Class Reference

#include <G4SynchrotronRadiationInMat.hh>

Inheritance diagram for G4SynchrotronRadiationInMat:
Collaboration 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)
 
G4VParticleChange * PostStepDoIt (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 G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (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
 
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 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
 
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)
 

Private Member Functions

G4SynchrotronRadiationInMatoperator= (const G4SynchrotronRadiationInMat &right)
 
 G4SynchrotronRadiationInMat (const G4SynchrotronRadiationInMat &)
 

Private Attributes

const G4double LowestKineticEnergy
 
const G4double HighestKineticEnergy
 
G4int TotBin
 
G4double CutInRange
 
const G4ParticleDefinitiontheGamma
 
const G4ParticleDefinitiontheElectron
 
const G4ParticleDefinitionthePositron
 
const G4doubleGammaCutInKineticEnergy
 
const G4doubleElectronCutInKineticEnergy
 
const G4doublePositronCutInKineticEnergy
 
const G4doubleParticleCutInKineticEnergy
 
G4double GammaCutInKineticEnergyNow
 
G4double ElectronCutInKineticEnergyNow
 
G4double PositronCutInKineticEnergyNow
 
G4double ParticleCutInKineticEnergyNow
 
G4double fAlpha
 
G4int fRootNumber
 
G4double fKsi
 
G4double fPsiGamma
 
G4double fEta
 
G4double fOrderAngleK
 
G4int fVerboseLevel
 
G4PropagatorInFieldfFieldPropagator
 

Static Private Attributes

static const G4double fLambdaConst
 
static const G4double fEnergyConst
 
static const G4double fIntegralProbabilityOfSR [200]
 

Additional Inherited Members

- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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() [1/2]

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),
135  fAlpha(0.0), fRootNumber(80),
137 {
139 
140  fFieldPropagator = transportMgr->GetPropagatorInField();
144  fPsiGamma = fEta = fOrderAngleK = 0.0;
145 }
const G4ParticleDefinition * theGamma
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4ParticleDefinition * theElectron
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4TransportationManager * GetTransportationManager()
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4ParticleDefinition * thePositron
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
G4PropagatorInField * GetPropagatorInField() const
Here is the call graph for this function:

◆ ~G4SynchrotronRadiationInMat()

G4SynchrotronRadiationInMat::~G4SynchrotronRadiationInMat ( )
virtual

Definition at line 152 of file G4SynchrotronRadiationInMat.cc.

153 {}

◆ G4SynchrotronRadiationInMat() [2/2]

G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat ( const G4SynchrotronRadiationInMat )
private

Member Function Documentation

◆ GetAngleK()

G4double G4SynchrotronRadiationInMat::GetAngleK ( G4double  eta)

Definition at line 637 of file G4SynchrotronRadiationInMat.cc.

638 {
639  fEta = eta; // should be > 0. !
640  G4int n;
641  G4double result, a;
642 
643  a = fAlpha; // always = 0.
644  n = fRootNumber; // around default = 80
645 
647 
648  result = integral.Laguerre(this,
650 
651  return result;
652 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
double G4double
Definition: G4Types.hh:76
G4double Laguerre(T &typeT, F f, G4double alpha, G4int n)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAngleNumberAtGammaKsi()

G4double G4SynchrotronRadiationInMat::GetAngleNumberAtGammaKsi ( G4double  gpsi)

Definition at line 658 of file G4SynchrotronRadiationInMat.cc.

659 {
660  G4double result, funK, funK2, gpsi2 = gpsi*gpsi;
661 
662  fPsiGamma = gpsi;
663  fEta = 0.5*fKsi*(1 + gpsi2)*std::sqrt(1 + gpsi2);
664 
665  fOrderAngleK = 1./3.;
666  funK = GetAngleK(fEta);
667  funK2 = funK*funK;
668 
669  result = gpsi2*funK2/(1 + gpsi2);
670 
671  fOrderAngleK = 2./3.;
672  funK = GetAngleK(fEta);
673  funK2 = funK*funK;
674 
675  result += funK2;
676  result *= (1 + gpsi2)*fKsi;
677 
678  return result;
679 }
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ 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 600 of file G4SynchrotronRadiationInMat.cc.

601 {
602  if (ksi <= 0.) return 1.0;
603  fKsi = ksi; // should be > 0. !
604  G4int n;
605  G4double result, a;
606 
607  a = fAlpha; // always = 0.
608  n = fRootNumber; // around default = 80
609 
611 
612  result = integral.Laguerre(this,
614 
615  result *= 9.*std::sqrt(3.)*ksi/8./pi;
616 
617  return result;
618 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
G4double Laguerre(T &typeT, F f, G4double alpha, G4int n)
Here is the call graph for this function:

◆ GetIntegrandForAngleK()

G4double G4SynchrotronRadiationInMat::GetIntegrandForAngleK ( G4double  t)

Definition at line 624 of file G4SynchrotronRadiationInMat.cc.

625 {
626  G4double result, hypCos=std::cosh(t);
627 
628  result = std::cosh(fOrderAngleK*t)*std::exp(t - fEta*hypCos); // fEta > 0. !
629  result /= hypCos;
630  return result;
631 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetIntProbSR()

G4double G4SynchrotronRadiationInMat::GetIntProbSR ( G4double  ksi)

Definition at line 561 of file G4SynchrotronRadiationInMat.cc.

562 {
563  if (ksi <= 0.) return 1.0;
564  fKsi = ksi; // should be > 0. !
565  G4int n;
566  G4double result, a;
567 
568  a = fAlpha; // always = 0.
569  n = fRootNumber; // around default = 80
570 
572 
573  result = integral.Laguerre(this,
575 
576  result *= 3./5./pi;
577 
578  return result;
579 }
int G4int
Definition: G4Types.hh:78
Char_t n[5]
static const double pi
Definition: G4SIunits.hh:74
double G4double
Definition: G4Types.hh:76
G4double Laguerre(T &typeT, F f, G4double alpha, G4int n)
Here is the call graph for this function:

◆ 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 
186  *condition = NotForced ;
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)
G4double GetMass() const
CLHEP::Hep3Vector G4ThreeVector
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
G4double GetTotalMomentum() const
G4double GetTotalEnergy() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
bool G4bool
Definition: G4Types.hh:79
double mag() const
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
double z() const
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
double G4double
Definition: G4Types.hh:76
static const G4double e3
#define DBL_MAX
Definition: templates.hh:83
G4double GetPDGCharge() const
const G4Field * GetDetectorField() const
Here is the call graph for this function:

◆ GetPhotonEnergy()

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

Definition at line 422 of file G4SynchrotronRadiationInMat.cc.

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

◆ GetProbSpectrumSRforEnergy()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforEnergy ( G4double  t)

Definition at line 585 of file G4SynchrotronRadiationInMat.cc.

586 {
587  G4double result, hypCos=std::cosh(t);
588 
589  result = std::cosh(5.*t/3.)*std::exp(t - fKsi*hypCos); // fKsi > 0. !
590  result /= hypCos;
591  return result;
592 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetProbSpectrumSRforInt()

G4double G4SynchrotronRadiationInMat::GetProbSpectrumSRforInt ( G4double  t)

Definition at line 545 of file G4SynchrotronRadiationInMat.cc.

546 {
547  G4double result, hypCos2, hypCos=std::cosh(t);
548 
549  hypCos2 = hypCos*hypCos;
550  result = std::cosh(5.*t/3.)*std::exp(t-fKsi*hypCos); // fKsi > 0. !
551  result /= hypCos2;
552  return result;
553 }
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetRandomEnergySR()

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

Definition at line 515 of file G4SynchrotronRadiationInMat.cc.

516 {
517  G4int i, iMax;
518  G4double energySR, random, position;
519 
520  iMax = 200;
521  random = G4UniformRand();
522 
523  for( i = 0; i < iMax; i++ )
524  {
525  if( random >= fIntegralProbabilityOfSR[i] ) break;
526  }
527  if(i <= 0 ) position = G4UniformRand(); // 0.
528  else if( i>= iMax) position = G4double(iMax);
529  else position = i + G4UniformRand(); // -1
530  //
531  // it was in initial implementation:
532  // energyOfSR = 0.0001*i*i*fEnergyConst*gamma*gamma*perpB ;
533 
534  energySR = 0.0001*position*position*fEnergyConst*gamma*gamma*perpB;
535 
536  if( energySR < 0. ) energySR = 0.;
537 
538  return energySR;
539 }
int G4int
Definition: G4Types.hh:78
static const G4double fIntegralProbabilityOfSR[200]
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ 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 }
const G4ParticleDefinition * theElectron
const G4ParticleDefinition * thePositron

◆ operator=()

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

◆ 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  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
342  while( fcos < G4UniformRand() );
343 
344  beta = std::sqrt(1. - 1./(gamma*gamma));
345 
346  cosTheta = (cosTheta + beta)/(1. + beta*cosTheta);
347 
348  if( cosTheta > 1. ) cosTheta = 1.;
349  if( cosTheta < -1. ) cosTheta = -1.;
350 
351  sinTheta = std::sqrt(1. - cosTheta*cosTheta );
352 
353  G4double Phi = twopi * G4UniformRand() ;
354 
355  G4double dirx = sinTheta*std::cos(Phi) ,
356  diry = sinTheta*std::sin(Phi) ,
357  dirz = cosTheta;
358 
359  G4ThreeVector gammaDirection ( dirx, diry, dirz);
360  gammaDirection.rotateUz(particleDirection);
361 
362  // polarization of new gamma
363 
364  // G4double sx = std::cos(Teta)*std::cos(Phi);
365  // G4double sy = std::cos(Teta)*std::sin(Phi);
366  // G4double sz = -std::sin(Teta);
367 
368  G4ThreeVector gammaPolarization = FieldValue.cross(gammaDirection);
369  gammaPolarization = gammaPolarization.unit();
370 
371  // (sx, sy, sz);
372  // gammaPolarization.rotateUz(particleDirection);
373 
374  // create G4DynamicParticle object for the SR photon
375 
377  gammaDirection,
378  energyOfSR );
379  aGamma->SetPolarization( gammaPolarization.x(),
380  gammaPolarization.y(),
381  gammaPolarization.z() );
382 
383 
384  aParticleChange.SetNumberOfSecondaries(1);
385  aParticleChange.AddSecondary(aGamma);
386 
387  // Update the incident particle
388 
389  G4double newKinEnergy = kineticEnergy - energyOfSR ;
390 
391  if (newKinEnergy > 0.)
392  {
393  aParticleChange.ProposeMomentumDirection( particleDirection );
394  aParticleChange.ProposeEnergy( newKinEnergy );
395  aParticleChange.ProposeLocalEnergyDeposit (0.);
396  }
397  else
398  {
399  aParticleChange.ProposeEnergy( 0. );
400  aParticleChange.ProposeLocalEnergyDeposit (0.);
401  G4double charge = aDynamicParticle->GetDefinition()->GetPDGCharge();
402  if (charge<0.)
403  {
404  aParticleChange.ProposeTrackStatus(fStopAndKill) ;
405  }
406  else
407  {
408  aParticleChange.ProposeTrackStatus(fStopButAlive) ;
409  }
410  }
411  }
412  else
413  {
414  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
415  }
416  }
417  return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
418 }
G4double GetMass() const
CLHEP::Hep3Vector G4ThreeVector
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
G4double fcos(G4double arg)
G4double GetTotalEnergy() const
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
bool G4bool
Definition: G4Types.hh:79
double mag() const
Hep3Vector unit() const
static const double twopi
Definition: G4SIunits.hh:75
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4FieldManager * FindAndSetFieldManager(G4VPhysicalVolume *pCurrentPhysVol)
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
double z() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
static const G4double e3
G4double GetRandomEnergySR(G4double, G4double)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetPDGCharge() const
const G4Field * GetDetectorField() const
Here is the call graph for this function:

◆ SetEta()

void G4SynchrotronRadiationInMat::SetEta ( G4double  eta)
inline

Definition at line 115 of file G4SynchrotronRadiationInMat.hh.

◆ SetKsi()

void G4SynchrotronRadiationInMat::SetKsi ( G4double  ksi)
inline

Definition at line 114 of file G4SynchrotronRadiationInMat.hh.

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

◆ SetRootNumber()

void G4SynchrotronRadiationInMat::SetRootNumber ( G4int  rn)
inline

◆ SetVerboseLevel()

void G4SynchrotronRadiationInMat::SetVerboseLevel ( G4int  v)
inline

Member Data Documentation

◆ CutInRange

G4double G4SynchrotronRadiationInMat::CutInRange
private

Definition at line 135 of file G4SynchrotronRadiationInMat.hh.

◆ ElectronCutInKineticEnergy

const G4double* G4SynchrotronRadiationInMat::ElectronCutInKineticEnergy
private

Definition at line 142 of file G4SynchrotronRadiationInMat.hh.

◆ ElectronCutInKineticEnergyNow

G4double G4SynchrotronRadiationInMat::ElectronCutInKineticEnergyNow
private

Definition at line 148 of file G4SynchrotronRadiationInMat.hh.

◆ fAlpha

G4double G4SynchrotronRadiationInMat::fAlpha
private

Definition at line 152 of file G4SynchrotronRadiationInMat.hh.

◆ fEnergyConst

const G4double G4SynchrotronRadiationInMat::fEnergyConst
staticprivate
Initial value:

Definition at line 123 of file G4SynchrotronRadiationInMat.hh.

◆ fEta

G4double G4SynchrotronRadiationInMat::fEta
private

Definition at line 156 of file G4SynchrotronRadiationInMat.hh.

◆ fFieldPropagator

G4PropagatorInField* G4SynchrotronRadiationInMat::fFieldPropagator
private

Definition at line 161 of file G4SynchrotronRadiationInMat.hh.

◆ fIntegralProbabilityOfSR

const G4double G4SynchrotronRadiationInMat::fIntegralProbabilityOfSR
staticprivate

Definition at line 125 of file G4SynchrotronRadiationInMat.hh.

◆ fKsi

G4double G4SynchrotronRadiationInMat::fKsi
private

Definition at line 154 of file G4SynchrotronRadiationInMat.hh.

◆ fLambdaConst

const G4double G4SynchrotronRadiationInMat::fLambdaConst
staticprivate
Initial value:

Definition at line 117 of file G4SynchrotronRadiationInMat.hh.

◆ fOrderAngleK

G4double G4SynchrotronRadiationInMat::fOrderAngleK
private

Definition at line 157 of file G4SynchrotronRadiationInMat.hh.

◆ fPsiGamma

G4double G4SynchrotronRadiationInMat::fPsiGamma
private

Definition at line 155 of file G4SynchrotronRadiationInMat.hh.

◆ fRootNumber

G4int G4SynchrotronRadiationInMat::fRootNumber
private

Definition at line 153 of file G4SynchrotronRadiationInMat.hh.

◆ fVerboseLevel

G4int G4SynchrotronRadiationInMat::fVerboseLevel
private

Definition at line 160 of file G4SynchrotronRadiationInMat.hh.

◆ GammaCutInKineticEnergy

const G4double* G4SynchrotronRadiationInMat::GammaCutInKineticEnergy
private

Definition at line 141 of file G4SynchrotronRadiationInMat.hh.

◆ GammaCutInKineticEnergyNow

G4double G4SynchrotronRadiationInMat::GammaCutInKineticEnergyNow
private

Definition at line 147 of file G4SynchrotronRadiationInMat.hh.

◆ HighestKineticEnergy

const G4double G4SynchrotronRadiationInMat::HighestKineticEnergy
private

Definition at line 131 of file G4SynchrotronRadiationInMat.hh.

◆ LowestKineticEnergy

const G4double G4SynchrotronRadiationInMat::LowestKineticEnergy
private

Definition at line 128 of file G4SynchrotronRadiationInMat.hh.

◆ ParticleCutInKineticEnergy

const G4double* G4SynchrotronRadiationInMat::ParticleCutInKineticEnergy
private

Definition at line 144 of file G4SynchrotronRadiationInMat.hh.

◆ ParticleCutInKineticEnergyNow

G4double G4SynchrotronRadiationInMat::ParticleCutInKineticEnergyNow
private

Definition at line 150 of file G4SynchrotronRadiationInMat.hh.

◆ PositronCutInKineticEnergy

const G4double* G4SynchrotronRadiationInMat::PositronCutInKineticEnergy
private

Definition at line 143 of file G4SynchrotronRadiationInMat.hh.

◆ PositronCutInKineticEnergyNow

G4double G4SynchrotronRadiationInMat::PositronCutInKineticEnergyNow
private

Definition at line 149 of file G4SynchrotronRadiationInMat.hh.

◆ theElectron

const G4ParticleDefinition* G4SynchrotronRadiationInMat::theElectron
private

Definition at line 138 of file G4SynchrotronRadiationInMat.hh.

◆ theGamma

const G4ParticleDefinition* G4SynchrotronRadiationInMat::theGamma
private

Definition at line 137 of file G4SynchrotronRadiationInMat.hh.

◆ thePositron

const G4ParticleDefinition* G4SynchrotronRadiationInMat::thePositron
private

Definition at line 139 of file G4SynchrotronRadiationInMat.hh.

◆ TotBin

G4int G4SynchrotronRadiationInMat::TotBin
private

Definition at line 133 of file G4SynchrotronRadiationInMat.hh.


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