src/BDSContinuousSR.cc

00001 
00008 //
00009 //  Synchrotron radiation energy loss process
00010 //
00011 
00012 #include <list>
00013 
00014 
00015 
00016 #include "BDSGlobalConstants.hh" 
00017 #include "BDSContinuousSR.hh"
00018 #include "G4ios.hh"
00019 #include "G4UnitsTable.hh"
00020 
00021 #include "BDSAcceleratorComponent.hh"
00022 #include "BDSBeamline.hh"
00023 
00024 
00025 BDSContinuousSR::BDSContinuousSR(const G4String& processName)
00026   : G4VDiscreteProcess(processName)
00027      // initialization
00028 {
00029 
00030   G4cout<<"initializing contSR"<<G4endl;
00031 
00032   nExpConst=5*fine_structure_const/(2*sqrt(3.0))/electron_mass_c2;
00033   CritEngFac=3./2.*hbarc/pow(electron_mass_c2,3);
00034 
00035 } 
00036  
00037 
00038 G4VParticleChange* 
00039 BDSContinuousSR::PostStepDoIt(const G4Track& trackData,
00040                                      const G4Step& stepData)
00041 {
00042   aParticleChange.Initialize(trackData);
00043 
00044   G4double eEnergy=trackData.GetTotalEnergy();
00045 
00046   G4double R=BDSLocalRadiusOfCurvature;
00047 
00048   G4double NewKinEnergy = trackData.GetKineticEnergy();
00049  
00050   G4double GamEnergy=0;
00051 
00052   aParticleChange.SetNumberOfSecondaries(0);
00053 
00054   if(fabs(R)>0) {
00055     G4double l = trackData.GetStep()->GetStepLength();
00056     const G4DynamicParticle* aParticle = trackData.GetDynamicParticle();
00057     G4double mass = aParticle->GetMass();
00058     G4double gamma = 1.e-3 * eEnergy / mass; // in 1.e3 units
00059     
00060     G4double r0 = 2.817940325e-3*m; // classical electron radius in 1..e-12 units
00061 
00062     // G4cout<<"mass="<<mass<<G4endl;
00063 //     G4cout<<"energy="<<eEnergy<<G4endl;
00064 //     G4cout<<"gamma="<<gamma<<G4endl;
00065 //     G4cout<<"R="<<R<<G4endl;  
00066    
00067     // energy loss (in MeV)
00068     GamEnergy = l * (gamma * gamma * gamma * gamma ) * 2. * r0 * mass / ( 3. * R * R) ;
00069     
00070     //G4cout<<"constSr : poststepdoit, l= "<<l<<" energy loss =" << GamEnergy  <<" MeV"<<G4endl;
00071   }
00072 
00073   NewKinEnergy -= GamEnergy;
00074 
00075   if(GamEnergy>0)
00076     {
00077       
00078       if (NewKinEnergy > 0.)
00079         {
00080           //
00081           // Update the incident particle 
00082           aParticleChange.ProposeEnergy(NewKinEnergy);
00083         } 
00084       else
00085         { 
00086           aParticleChange.ProposeEnergy( 0. );
00087           aParticleChange.ProposeLocalEnergyDeposit (0.);
00088           G4double charge= trackData.GetDynamicParticle()->GetCharge();
00089           if (charge<0.) aParticleChange.ProposeTrackStatus(fStopAndKill);
00090           else       aParticleChange.ProposeTrackStatus(fStopButAlive);
00091         }
00092       
00093     }
00094   
00095   return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00096 }
00097 
00098 
00099 BDSContinuousSR::~BDSContinuousSR(){
00100 }
00101 

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7