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

#include <G4QSynchRad.hh>

+ Inheritance diagram for G4QSynchRad:

Public Member Functions

 G4QSynchRad (const G4String &processName="CHIPS_SynchrotronRadiation", G4ProcessType type=fElectromagnetic)
 
virtual ~G4QSynchRad ()
 
G4double GetMeanFreePath (const G4Track &track, G4double step, G4ForceCondition *fCond)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
G4bool IsApplicable (const G4ParticleDefinition &pd)
 
void SetMinGamma (G4double ming)
 
G4double GetMinGamma ()
 
G4double GetRadius (const G4Track &track)
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double 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 53 of file G4QSynchRad.hh.

Constructor & Destructor Documentation

◆ G4QSynchRad()

G4QSynchRad::G4QSynchRad ( const G4String processName = "CHIPS_SynchrotronRadiation",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 51 of file G4QSynchRad.cc.

51 :
52 G4VDiscreteProcess (Name, Type), minGamma(227.), Polarization(0.,0.,1.) {
53 G4HadronicDeprecate("G4QSynchRad");
54}
#define G4HadronicDeprecate(name)

◆ ~G4QSynchRad()

virtual G4QSynchRad::~G4QSynchRad ( )
inlinevirtual

Definition at line 60 of file G4QSynchRad.hh.

60{;}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4QSynchRad::GetMeanFreePath ( const G4Track track,
G4double  step,
G4ForceCondition fCond 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 57 of file G4QSynchRad.cc.

58{
59 static const G4double coef = 0.4*std::sqrt(3.)/fine_structure_const;
60 const G4DynamicParticle* particle = track.GetDynamicParticle();
61 *cond = NotForced ;
62 G4double gamma = particle->GetTotalEnergy() / particle->GetMass();
63#ifdef debug
64 G4cout<<"G4QSynchRad::MeanFreePath: gamma = "<<gamma<<G4endl;
65#endif
66 G4double MFP = DBL_MAX;
67 if( gamma > minGamma ) // For smalle gamma neglect the process
68 {
69 G4double R = GetRadius(track);
70#ifdef debug
71 G4cout<<"G4QSynchRad::MeanFreePath: Radius = "<<R/meter<<" [m]"<<G4endl;
72#endif
73 if(R > 0.) MFP= coef*R/gamma;
74 }
75#ifdef debug
76 G4cout<<"G4QSynchRad::MeanFreePath = "<<MFP/centimeter<<" [cm]"<<G4endl;
77#endif
78 return MFP;
79}
@ NotForced
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetMass() const
G4double GetTotalEnergy() const
G4double GetRadius(const G4Track &track)
Definition: G4QSynchRad.cc:149
const G4DynamicParticle * GetDynamicParticle() const
#define DBL_MAX
Definition: templates.hh:83

◆ GetMinGamma()

G4double G4QSynchRad::GetMinGamma ( )
inline

Definition at line 75 of file G4QSynchRad.hh.

75{return minGamma;}

◆ GetRadius()

G4double G4QSynchRad::GetRadius ( const G4Track track)

Definition at line 149 of file G4QSynchRad.cc.

150{
151 static const G4double unk = meter*tesla/0.3/gigaelectronvolt;
152 const G4DynamicParticle* particle = track.GetDynamicParticle();
153 G4double z = particle->GetDefinition()->GetPDGCharge();
154 if(z == 0.) return 0.; // --> neutral particle
155 if(z < 0.) z=-z;
157 G4PropagatorInField* Field = transMan->GetPropagatorInField();
158 G4FieldManager* fMan = Field->FindAndSetFieldManager(track.GetVolume());
159 if(!fMan || !fMan->GetDetectorField()) return 0.; // --> no field at all
160 const G4Field* pField = fMan->GetDetectorField();
162 G4double PosArray[3]={position.x(), position.y(), position.z()};
163 G4double BArray[3];
164 pField->GetFieldValue(PosArray, BArray);
165 G4ThreeVector B3D(BArray[0], BArray[1], BArray[2]);
166#ifdef debug
167 G4cout<<"G4QSynchRad::GetRadius: Pos="<<position/meter<<", B(tesla)="<<B3D/tesla<<G4endl;
168#endif
169 G4ThreeVector MomDir = particle->GetMomentumDirection();
170 G4ThreeVector Ort = B3D.cross(MomDir);
171 G4double OrtB = Ort.mag(); // not negative (independent units)
172 if(OrtB == 0.) return 0.; // --> along the field line
173 Polarization = Ort/OrtB; // Polarization unit vector
174 G4double mom = particle->GetTotalMomentum(); // Momentum of the particle
175#ifdef debug
176 G4cout<<"G4QSynchRad::GetRadius: P(GeV)="<<mom/GeV<<", B(tesla)="<<OrtB/tesla<<G4endl;
177#endif
178 // R [m]= mom [GeV]/(0.3 * z * OrtB [tesla])
179 return mom * unk / z / OrtB;
180}
Hep3Vector cross(const Hep3Vector &) const
double mag() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() 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)
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
static G4TransportationManager * GetTransportationManager()
G4PropagatorInField * GetPropagatorInField() const

Referenced by GetMeanFreePath(), and PostStepDoIt().

◆ IsApplicable()

G4bool G4QSynchRad::IsApplicable ( const G4ParticleDefinition pd)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 72 of file G4QSynchRad.hh.

72{ return (pd.GetPDGCharge() != 0.); }

◆ PostStepDoIt()

G4VParticleChange * G4QSynchRad::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 81 of file G4QSynchRad.cc.

83{
84 static const G4double hc = 1.5 * c_light * hbar_Planck; // E_c=h*w_c=1.5*(hc)*(gamma^3)/R
86 const G4DynamicParticle* particle=track.GetDynamicParticle();
87 G4double gamma = particle->GetTotalEnergy() / particle->GetMass();
88 if(gamma <= minGamma )
89 {
90#ifdef debug
91 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt is called for small gamma="<<gamma<<G4endl;
92#endif
93 return G4VDiscreteProcess::PostStepDoIt(track,step);
94 }
95 // Photon energy calculation (E < 8.1*Ec restriction)
96 G4double R = GetRadius(track);
97 if(R <= 0.)
98 {
99#ifdef debug
100 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ radius ="
101 <<R/meter<<" [m]"<<G4endl;
102#endif
103 return G4VDiscreteProcess::PostStepDoIt(track, step);
104 }
105 G4double EPhoton = hc * gamma * gamma * gamma / R; // E_c
106 G4double dd=5.e-8;
107 G4double rnd=G4UniformRand()*(1.+dd);
108 if (rnd < 0.5 ) EPhoton *= .65 * rnd * rnd * rnd;
109 else if(rnd > .997) EPhoton *= 15.-1.03*std::log((1.-rnd)/dd+1.);
110 else
111 {
112 G4double r2=rnd*rnd;
113 G4double dr=1.-rnd;
114 EPhoton*=(2806.+28./rnd)/(1.+500./r2/r2+6500.*(std::sqrt(dr)+28.*dr*dr*dr));
115 }
116#ifdef debug
117 G4cout<<"G4SynchRad::PostStepDoIt: PhotonEnergy = "<<EPhoton/keV<<" [keV]"<<G4endl;
118#endif
119 if(EPhoton <= 0.)
120 {
121 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ photon energy="
122 <<EPhoton/keV<<" [keV]"<<G4endl;
123 return G4VDiscreteProcess::PostStepDoIt(track, step);
124 }
125 G4double kinEn = particle->GetKineticEnergy();
126 G4double newEn = kinEn - EPhoton ;
127 if (newEn > 0.)
128 {
131 }
132 else // Very low probable event
133 {
134 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: PhotonEnergy > TotalKinEnergy"<<G4endl;
135 EPhoton = kinEn;
139 }
140 G4ThreeVector MomDir = particle->GetMomentumDirection();
141 G4DynamicParticle* Photon = new G4DynamicParticle(G4Gamma::Gamma(), MomDir, EPhoton);
142 Photon->SetPolarization(Polarization.x(), Polarization.y(), Polarization.z());
145 return G4VDiscreteProcess::PostStepDoIt(track,step);
146}
@ fStopButAlive
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double y() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
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
const G4double hc
[MeV*fm]

◆ SetMinGamma()

void G4QSynchRad::SetMinGamma ( G4double  ming)
inline

Definition at line 74 of file G4QSynchRad.hh.

74{minGamma = ming;}

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