/home/cern/BDSIM_new/include/BDSContinuousSR.hh

00001 
00002 #ifndef BDSContinuousSR_h
00003 #define BDSContinuousSR_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 BDSMaterials* theMaterials;
00035 extern G4double BDSLocalRadiusOfCurvature;
00036  
00037 class BDSContinuousSR : public G4VDiscreteProcess 
00038 { 
00039 public:
00040  
00041   BDSContinuousSR(const G4String& processName = "contSR");
00042  
00043   ~BDSContinuousSR();
00044 
00045   G4bool IsApplicable(const G4ParticleDefinition&);
00046      
00047   G4double GetMeanFreePath(const G4Track& track,
00048                            G4double previousStepSize,
00049                            G4ForceCondition* condition );
00050  
00051   G4VParticleChange *PostStepDoIt(const G4Track& track,         
00052                                   const G4Step&  step);                 
00053 
00054   G4double SynGenC(G4double xmin);
00055   G4double SynRadC(G4double x);
00056 
00057 protected:
00058 
00059 private:
00060 
00061   BDSContinuousSR & operator=(const BDSContinuousSR &right);
00062      
00063   BDSContinuousSR(const BDSContinuousSR&);
00064 
00065   G4double nExpConst;
00066   G4double CritEngFac;
00067 
00068 private:
00069 };
00070 
00071 inline G4bool 
00072 BDSContinuousSR::IsApplicable(const G4ParticleDefinition& particle)
00073 {
00074   return(  (&particle == G4Electron::Electron())
00075            ||(&particle == G4Positron::Positron()) );
00076 }
00077 
00078 inline G4double 
00079 BDSContinuousSR::GetMeanFreePath(const G4Track& track,
00080                                         G4double PreviousStepSize,
00081                                         G4ForceCondition* ForceCondition)
00082 {  
00083   *ForceCondition = NotForced ;
00084 
00085   G4double rfact = 1.; // mean free path reduction factor 
00086                         // to reduce fluctuations 
00087 
00088   //   static G4FieldManager* lastFieldManager;
00089 
00090   G4double MeanFreePath;
00091   G4FieldManager* TheFieldManager=
00092     track.GetVolume()->GetLogicalVolume()->GetFieldManager();
00093 
00094   if(track.GetTotalEnergy()<BDSGlobals->GetThresholdCutCharged())
00095     return DBL_MAX;
00096   /*
00097   G4double SynchOnZPos = (7.184+4.0) * m;
00098   if(track.GetPosition().z() + BDSGlobals->GetWorldSizeZ() < SynchOnZPos)
00099     return DBL_MAX;
00100   */
00101   if(TheFieldManager)
00102     {
00103       const G4Field* pField = TheFieldManager->GetDetectorField() ;
00104       G4ThreeVector  globPosition = track.GetPosition() ;
00105       G4double  globPosVec[3], FieldValueVec[3] ;
00106       globPosVec[0] = globPosition.x() ;
00107       globPosVec[1] = globPosition.y() ;
00108       globPosVec[2] = globPosition.z() ;
00109       
00110       if(pField)
00111         pField->GetFieldValue( globPosVec, FieldValueVec ) ; 
00112         
00113       G4double Blocal= FieldValueVec[1];
00114       if ( FieldValueVec[0]!=0.)
00115         Blocal=sqrt(Blocal*Blocal+FieldValueVec[0]*FieldValueVec[0]);
00116 
00117       
00118  
00119       if(track.GetMaterial()==theMaterials->GetMaterial("Vacuum")&&Blocal !=0 )
00120         {
00121           G4ThreeVector InitMag=track.GetMomentum();
00122 
00123           // transform to the local coordinate frame
00124 
00125           G4AffineTransform tf(track.GetTouchable()->GetHistory()->GetTopTransform().Inverse());
00126           const G4RotationMatrix Rot=tf.NetRotation();
00127           const G4ThreeVector Trans=-tf.NetTranslation();
00128           InitMag=Rot*InitMag; 
00129 
00130 
00131           G4double Rlocal=(InitMag.z()/GeV)/(0.299792458 * Blocal/tesla) *m;
00132           
00133           MeanFreePath=
00134             fabs(Rlocal)/(track.GetTotalEnergy()*nExpConst);
00135 
00136           
00137 
00138          //  G4cout<<"*****************SR*************************"<<G4endl;
00139 //        G4cout<<"Track momentum: "<<InitMag<<G4endl;;
00140 //        G4cout<<"Blocal="<<Blocal/tesla<<"  Rlocal="<<Rlocal/m<<G4endl;
00141 //        G4cout<<track.GetVolume()->GetName()<<" mfp="<<MeanFreePath/m<<G4endl;
00142 //        G4cout<<"********************************************"<<G4endl;
00143 
00144           return rfact * MeanFreePath;
00145         }
00146       else
00147         return DBL_MAX;
00148     }
00149   else
00150     return DBL_MAX;
00151   
00152 }
00153   
00154 
00155 
00156 #endif

Generated on Wed Mar 5 17:25:21 2008 for BDSIM by  doxygen 1.5.3