#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) |
G4VParticleChange * | PostStepDoIt (const G4Track &track, const G4Step &step) |
G4bool | IsApplicable (const G4ParticleDefinition &pd) |
void | SetMinGamma (G4double ming) |
G4double | GetMinGamma () |
G4double | GetRadius (const G4Track &track) |
Definition at line 53 of file G4QSynchRad.hh.
G4QSynchRad::G4QSynchRad | ( | const G4String & | processName = "CHIPS_SynchrotronRadiation" , |
|
G4ProcessType | type = fElectromagnetic | |||
) |
Definition at line 51 of file G4QSynchRad.cc.
References G4HadronicDeprecate.
00051 : 00052 G4VDiscreteProcess (Name, Type), minGamma(227.), Polarization(0.,0.,1.) { 00053 G4HadronicDeprecate("G4QSynchRad"); 00054 }
virtual G4QSynchRad::~G4QSynchRad | ( | ) | [inline, virtual] |
G4double G4QSynchRad::GetMeanFreePath | ( | const G4Track & | track, | |
G4double | step, | |||
G4ForceCondition * | fCond | |||
) | [virtual] |
Implements G4VDiscreteProcess.
Definition at line 57 of file G4QSynchRad.cc.
References DBL_MAX, G4cout, G4endl, G4Track::GetDynamicParticle(), G4DynamicParticle::GetMass(), GetRadius(), G4DynamicParticle::GetTotalEnergy(), and NotForced.
00058 { 00059 static const G4double coef = 0.4*std::sqrt(3.)/fine_structure_const; 00060 const G4DynamicParticle* particle = track.GetDynamicParticle(); 00061 *cond = NotForced ; 00062 G4double gamma = particle->GetTotalEnergy() / particle->GetMass(); 00063 #ifdef debug 00064 G4cout<<"G4QSynchRad::MeanFreePath: gamma = "<<gamma<<G4endl; 00065 #endif 00066 G4double MFP = DBL_MAX; 00067 if( gamma > minGamma ) // For smalle gamma neglect the process 00068 { 00069 G4double R = GetRadius(track); 00070 #ifdef debug 00071 G4cout<<"G4QSynchRad::MeanFreePath: Radius = "<<R/meter<<" [m]"<<G4endl; 00072 #endif 00073 if(R > 0.) MFP= coef*R/gamma; 00074 } 00075 #ifdef debug 00076 G4cout<<"G4QSynchRad::MeanFreePath = "<<MFP/centimeter<<" [cm]"<<G4endl; 00077 #endif 00078 return MFP; 00079 }
G4double G4QSynchRad::GetMinGamma | ( | ) | [inline] |
Definition at line 149 of file G4QSynchRad.cc.
References G4PropagatorInField::FindAndSetFieldManager(), G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4FieldManager::GetDetectorField(), G4Track::GetDynamicParticle(), G4Field::GetFieldValue(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGCharge(), G4Track::GetPosition(), G4TransportationManager::GetPropagatorInField(), G4DynamicParticle::GetTotalMomentum(), G4TransportationManager::GetTransportationManager(), G4Track::GetVolume(), and position.
Referenced by GetMeanFreePath(), and PostStepDoIt().
00150 { 00151 static const G4double unk = meter*tesla/0.3/gigaelectronvolt; 00152 const G4DynamicParticle* particle = track.GetDynamicParticle(); 00153 G4double z = particle->GetDefinition()->GetPDGCharge(); 00154 if(z == 0.) return 0.; // --> neutral particle 00155 if(z < 0.) z=-z; 00156 G4TransportationManager* transMan = G4TransportationManager::GetTransportationManager(); 00157 G4PropagatorInField* Field = transMan->GetPropagatorInField(); 00158 G4FieldManager* fMan = Field->FindAndSetFieldManager(track.GetVolume()); 00159 if(!fMan || !fMan->GetDetectorField()) return 0.; // --> no field at all 00160 const G4Field* pField = fMan->GetDetectorField(); 00161 G4ThreeVector position = track.GetPosition(); 00162 G4double PosArray[3]={position.x(), position.y(), position.z()}; 00163 G4double BArray[3]; 00164 pField->GetFieldValue(PosArray, BArray); 00165 G4ThreeVector B3D(BArray[0], BArray[1], BArray[2]); 00166 #ifdef debug 00167 G4cout<<"G4QSynchRad::GetRadius: Pos="<<position/meter<<", B(tesla)="<<B3D/tesla<<G4endl; 00168 #endif 00169 G4ThreeVector MomDir = particle->GetMomentumDirection(); 00170 G4ThreeVector Ort = B3D.cross(MomDir); 00171 G4double OrtB = Ort.mag(); // not negative (independent units) 00172 if(OrtB == 0.) return 0.; // --> along the field line 00173 Polarization = Ort/OrtB; // Polarization unit vector 00174 G4double mom = particle->GetTotalMomentum(); // Momentum of the particle 00175 #ifdef debug 00176 G4cout<<"G4QSynchRad::GetRadius: P(GeV)="<<mom/GeV<<", B(tesla)="<<OrtB/tesla<<G4endl; 00177 #endif 00178 // R [m]= mom [GeV]/(0.3 * z * OrtB [tesla]) 00179 return mom * unk / z / OrtB; 00180 }
G4bool G4QSynchRad::IsApplicable | ( | const G4ParticleDefinition & | pd | ) | [inline, virtual] |
Reimplemented from G4VProcess.
Definition at line 72 of file G4QSynchRad.hh.
References G4ParticleDefinition::GetPDGCharge().
00072 { return (pd.GetPDGCharge() != 0.); }
G4VParticleChange * G4QSynchRad::PostStepDoIt | ( | const G4Track & | track, | |
const G4Step & | step | |||
) | [virtual] |
Reimplemented from G4VDiscreteProcess.
Definition at line 81 of file G4QSynchRad.cc.
References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, fStopButAlive, G4cout, G4endl, G4UniformRand, G4Gamma::Gamma(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4DynamicParticle::GetMomentumDirection(), GetRadius(), G4DynamicParticle::GetTotalEnergy(), G4INCL::PhysicalConstants::hc, G4ParticleChange::Initialize(), G4VDiscreteProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::SetNumberOfSecondaries(), and G4DynamicParticle::SetPolarization().
00083 { 00084 static const G4double hc = 1.5 * c_light * hbar_Planck; // E_c=h*w_c=1.5*(hc)*(gamma^3)/R 00085 aParticleChange.Initialize(track); 00086 const G4DynamicParticle* particle=track.GetDynamicParticle(); 00087 G4double gamma = particle->GetTotalEnergy() / particle->GetMass(); 00088 if(gamma <= minGamma ) 00089 { 00090 #ifdef debug 00091 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt is called for small gamma="<<gamma<<G4endl; 00092 #endif 00093 return G4VDiscreteProcess::PostStepDoIt(track,step); 00094 } 00095 // Photon energy calculation (E < 8.1*Ec restriction) 00096 G4double R = GetRadius(track); 00097 if(R <= 0.) 00098 { 00099 #ifdef debug 00100 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ radius =" 00101 <<R/meter<<" [m]"<<G4endl; 00102 #endif 00103 return G4VDiscreteProcess::PostStepDoIt(track, step); 00104 } 00105 G4double EPhoton = hc * gamma * gamma * gamma / R; // E_c 00106 G4double dd=5.e-8; 00107 G4double rnd=G4UniformRand()*(1.+dd); 00108 if (rnd < 0.5 ) EPhoton *= .65 * rnd * rnd * rnd; 00109 else if(rnd > .997) EPhoton *= 15.-1.03*std::log((1.-rnd)/dd+1.); 00110 else 00111 { 00112 G4double r2=rnd*rnd; 00113 G4double dr=1.-rnd; 00114 EPhoton*=(2806.+28./rnd)/(1.+500./r2/r2+6500.*(std::sqrt(dr)+28.*dr*dr*dr)); 00115 } 00116 #ifdef debug 00117 G4cout<<"G4SynchRad::PostStepDoIt: PhotonEnergy = "<<EPhoton/keV<<" [keV]"<<G4endl; 00118 #endif 00119 if(EPhoton <= 0.) 00120 { 00121 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: zero or negativ photon energy=" 00122 <<EPhoton/keV<<" [keV]"<<G4endl; 00123 return G4VDiscreteProcess::PostStepDoIt(track, step); 00124 } 00125 G4double kinEn = particle->GetKineticEnergy(); 00126 G4double newEn = kinEn - EPhoton ; 00127 if (newEn > 0.) 00128 { 00129 aParticleChange.ProposeEnergy(newEn); 00130 aParticleChange.ProposeLocalEnergyDeposit (0.); 00131 } 00132 else // Very low probable event 00133 { 00134 G4cout<<"-Warning-G4QSynchRad::PostStepDoIt: PhotonEnergy > TotalKinEnergy"<<G4endl; 00135 EPhoton = kinEn; 00136 aParticleChange.ProposeEnergy(0.); 00137 aParticleChange.ProposeLocalEnergyDeposit(0.); 00138 aParticleChange.ProposeTrackStatus(fStopButAlive) ; 00139 } 00140 G4ThreeVector MomDir = particle->GetMomentumDirection(); 00141 G4DynamicParticle* Photon = new G4DynamicParticle(G4Gamma::Gamma(), MomDir, EPhoton); 00142 Photon->SetPolarization(Polarization.x(), Polarization.y(), Polarization.z()); 00143 aParticleChange.SetNumberOfSecondaries(1); 00144 aParticleChange.AddSecondary(Photon); 00145 return G4VDiscreteProcess::PostStepDoIt(track,step); 00146 }
void G4QSynchRad::SetMinGamma | ( | G4double | ming | ) | [inline] |