include/BDSSynchrotronRadiation.hh

00001 
00002 #ifndef BDSSynchrotronRadiation_h
00003 #define BDSSynchrotronRadiation_h 1
00004 
00005 #include "BDSGlobalConstants.hh"
00006 
00007 #include "G4ios.hh" 
00008 #include "globals.hh"
00009 #include "Randomize.hh" 
00010 #include "G4VDiscreteProcess.hh"
00011 #include "G4Track.hh"
00012 #include "G4Step.hh"
00013 #include "G4Gamma.hh"
00014 #include "G4Electron.hh"
00015 #include "G4Positron.hh"
00016 #include "G4OrderedTable.hh" 
00017 #include "G4PhysicsTable.hh"
00018 #include "G4PhysicsLogVector.hh"
00019 #include "BDSComptonEngine.hh"
00020 #include "BDSMaterials.hh"
00021 #include "Randomize.hh"
00022 
00023 #include "G4ChordFinder.hh"
00024 #include "G4FieldManager.hh"
00025 #include "G4MagIntegratorDriver.hh"
00026 #include "G4MagIntegratorStepper.hh"
00027 #include "G4FieldTrack.hh"
00028 
00029 #include "G4MagIntegratorDriver.hh"
00030 
00031 #include "G4Navigator.hh"
00032 #include "G4AffineTransform.hh"
00033 
00034 extern G4double BDSLocalRadiusOfCurvature;
00035  
00036 class BDSSynchrotronRadiation : public G4VDiscreteProcess 
00037 { 
00038 public:
00039  
00040   BDSSynchrotronRadiation(const G4String& processName = "BDSSynchRad");
00041  
00042   ~BDSSynchrotronRadiation();
00043 
00044   G4bool IsApplicable(const G4ParticleDefinition&);
00045      
00046   G4double GetMeanFreePath(const G4Track& track,
00047                            G4double previousStepSize,
00048                            G4ForceCondition* condition );
00049  
00050   G4VParticleChange *PostStepDoIt(const G4Track& track,         
00051                                   const G4Step&  step);                 
00052 
00053   G4double SynGenC(G4double xmin);
00054   G4double SynRadC(G4double x);
00055 
00056 protected:
00057 
00058 private:
00059 
00060   // assignment and copy constructor not implemented nor used
00061   BDSSynchrotronRadiation & operator=(const BDSSynchrotronRadiation &right);
00062   BDSSynchrotronRadiation(const BDSSynchrotronRadiation&);
00063 
00064   G4double nExpConst;
00065   G4double CritEngFac;
00066   G4int MeanFreePathCounter;
00067 
00068 private:
00069 };
00070 
00071 inline G4bool 
00072 BDSSynchrotronRadiation::IsApplicable(const G4ParticleDefinition& particle)
00073 {
00074   return(  (&particle == G4Electron::Electron())
00075            ||(&particle == G4Positron::Positron()) );
00076 }
00077 
00078 inline G4double 
00079 BDSSynchrotronRadiation::GetMeanFreePath(const G4Track& track,
00080                                          G4double /*PreviousStepSize*/,
00081                                          G4ForceCondition* ForceCondition)
00082 {  
00083   *ForceCondition = NotForced ;
00084 
00085   //   static G4FieldManager* lastFieldManager;
00086 
00087   G4double MeanFreePath;
00088   G4FieldManager* TheFieldManager=
00089     track.GetVolume()->GetLogicalVolume()->GetFieldManager();
00090 
00091   if(track.GetTotalEnergy()<BDSGlobalConstants::Instance()->GetThresholdCutCharged())
00092     return DBL_MAX;
00093   /*
00094   G4double SynchOnZPos = (7.184+4.0) * m;
00095   if(track.GetPosition().z() + BDSGlobalConstants::Instance()->GetWorldSizeZ() < SynchOnZPos)
00096     return DBL_MAX;
00097   */
00098   if(TheFieldManager)
00099     {
00100       const G4Field* pField = TheFieldManager->GetDetectorField() ;
00101       G4ThreeVector  globPosition = track.GetPosition() ;
00102       G4double  globPosVec[3], FieldValueVec[3]={0.} ;
00103       globPosVec[0] = globPosition.x() ;
00104       globPosVec[1] = globPosition.y() ;
00105       globPosVec[2] = globPosition.z() ;
00106       
00107       if(pField)
00108         pField->GetFieldValue( globPosVec, FieldValueVec ) ; 
00109         
00110       G4double Blocal= FieldValueVec[1];
00111       if ( FieldValueVec[0]!=0.)
00112         Blocal=sqrt(Blocal*Blocal+FieldValueVec[0]*FieldValueVec[0]);
00113 
00114       
00115  
00116       if(track.GetMaterial()==BDSMaterials::Instance()->GetMaterial("Vacuum")&&Blocal !=0 )
00117         {
00118           G4ThreeVector InitMag=track.GetMomentum();
00119 
00120           // transform to the local coordinate frame
00121 
00122           G4AffineTransform tf(track.GetTouchable()->GetHistory()->GetTopTransform().Inverse());
00123           const G4RotationMatrix Rot=tf.NetRotation();
00124           const G4ThreeVector Trans=-tf.NetTranslation();
00125           InitMag=Rot*InitMag; 
00126 
00127 
00128           G4double Rlocal=(InitMag.z()/GeV)/(0.299792458 * Blocal/tesla) *m;
00129           
00130           MeanFreePath=
00131             fabs(Rlocal)/(track.GetTotalEnergy()*nExpConst);
00132 
00133           MeanFreePath /= BDSGlobalConstants::Instance()->GetSynchMeanFreeFactor();
00134 
00135           if(MeanFreePathCounter==BDSGlobalConstants::Instance()->GetSynchMeanFreeFactor())
00136             MeanFreePathCounter=0;
00137 
00138 #ifdef DEBUG
00139           G4cout<<"*****************SR*************************"<<G4endl;
00140           G4cout<<"Track momentum: "<<InitMag<<G4endl;
00141           G4cout<<"Blocal="<<Blocal/tesla<<"  Rlocal="<<Rlocal/m<<G4endl;
00142           G4cout<<track.GetVolume()->GetName()<<" mfp="<<MeanFreePath/m<<G4endl;
00143           G4cout<<"********************************************"<<G4endl;
00144 #endif
00145 
00146           return MeanFreePath;
00147         }
00148       else
00149         return DBL_MAX;
00150     }
00151   else
00152     return DBL_MAX;
00153   
00154 }
00155   
00156 
00157 
00158 #endif

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7