/home/cern/BDSIM_new/include/BDSLaserCompton.hh

00001 /* BDSIM code.    Version 1.0
00002    Author: Grahame A. Blair, Royal Holloway, Univ. of London.
00003    Last modified 24.7.2002
00004    Copyright (c) 2002 by G.A.Blair.  ALL RIGHTS RESERVED. 
00005 */
00006 #ifndef BDSLaserCompton_h
00007 #define BDSLaserCompton_h 1
00008 
00009 #include "G4ios.hh" 
00010 #include "globals.hh"
00011 #include "Randomize.hh" 
00012 #if G4VERSION > 8
00013 #include "G4VDiscreteProcess.hh"
00014 #else
00015 #include "G4VeEnergyLoss.hh"
00016 #endif
00017 #include "G4VPhysicsConstructor.hh"
00018 #include "G4Track.hh"
00019 #include "G4Step.hh"
00020 #include "G4Gamma.hh"
00021 #include "G4Electron.hh"
00022 #include "G4Positron.hh"
00023 #include "G4OrderedTable.hh" 
00024 #include "G4PhysicsTable.hh"
00025 #include "G4PhysicsLogVector.hh"
00026 #include "BDSComptonEngine.hh"
00027 #include "BDSMaterials.hh"
00028 
00029 
00030 extern BDSMaterials* theMaterials;
00031 // flag initiated in BDSEventAction
00032 extern G4bool FireLaserCompton;
00033 
00034 #if G4VERSION > 8
00035 class BDSLaserCompton : public G4VDiscreteProcess
00036 #else
00037 class BDSLaserCompton : public G4VeEnergyLoss
00038 #endif
00039 { 
00040   public:
00041    
00042   BDSLaserCompton(const G4String& processName = "eLaser");
00043   
00044   ~BDSLaserCompton();
00045 
00046 #if G4VERSION > 8  
00047 //  virtual void PrintInfo();
00048 #endif
00049 
00050   G4bool IsApplicable(const G4ParticleDefinition&);
00051   
00052   G4double GetMeanFreePath(const G4Track& track,
00053                            G4double previousStepSize,
00054                            G4ForceCondition* condition );
00055   
00056   G4VParticleChange *PostStepDoIt(const G4Track& track,         
00057                                   const G4Step&  step);                 
00058   
00059   inline void SetLaserDirection(G4ThreeVector aDirection);
00060   inline G4ThreeVector GetLaserDirection();
00061   
00062   inline void SetLaserWavelength(G4double aWavelength);
00063   inline G4double GetLaserWavelength();
00064   
00065 protected:
00066   
00067   G4double ComputeMeanFreePath( const G4ParticleDefinition* ParticleType,
00068                                 G4double KineticEnergy, 
00069                                 const G4Material* aMaterial);
00070 
00071   virtual G4double SecondaryEnergyThreshold(size_t index);
00072 
00073 protected:
00074 #if G4VERSION > 8
00075   G4bool isInitialised;  
00076   const G4ParticleDefinition* particle;
00077 
00078 //  virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*, const G4ParticleDefinition*);
00079 #endif
00080   
00081 private:
00082   
00083   BDSLaserCompton & operator=(const BDSLaserCompton &right);
00084   
00085   BDSLaserCompton(const BDSLaserCompton&);
00086 
00087   const std::vector<G4double>* secondaryEnergyCuts;
00088   
00089 private:
00090   G4double itsLaserWavelength;
00091   G4ThreeVector itsLaserDirection;
00092   G4double itsLaserEnergy;
00093   BDSComptonEngine* itsComptonEngine;
00094   G4Material* itsLastMaterial;
00095 
00096 };
00097 inline G4bool BDSLaserCompton::IsApplicable(
00098                                             const G4ParticleDefinition& particle)
00099 {
00100   return(  (&particle == G4Electron::Electron())
00101            ||(&particle == G4Positron::Positron()) );
00102 }
00103 
00104 inline G4double BDSLaserCompton::GetMeanFreePath(const G4Track& track,
00105                                                  G4double PreviousStepSize,
00106                                                  G4ForceCondition* ForceCondition)
00107 {
00108   //  G4cout<<" FireLaserCompton="<<FireLaserCompton<<G4endl;
00109   if( track.GetMaterial()==theMaterials->GetMaterial("LaserVac") &&
00110       //      itsLastMaterial!=theMaterials->GetMaterial("LaserVac") )
00111 // flag initiated in BDSEventActionxem
00112       FireLaserCompton)
00113     {
00114       //return 0;}
00115       *ForceCondition=Forced;}
00116   // itsLastMaterial=track.GetMaterial();
00117   return DBL_MAX;
00118 }
00119 
00120 
00121 inline void BDSLaserCompton::SetLaserDirection(G4ThreeVector aDirection)
00122 {itsLaserDirection=aDirection;}
00123 
00124 inline G4ThreeVector BDSLaserCompton::GetLaserDirection()
00125 {return itsLaserDirection;}
00126 
00127 inline void BDSLaserCompton::SetLaserWavelength(G4double aWavelength)
00128 {itsLaserWavelength=aWavelength;}
00129 
00130 inline G4double BDSLaserCompton::GetLaserWavelength()
00131 {return itsLaserWavelength;}
00132 
00133 inline G4double BDSLaserCompton::SecondaryEnergyThreshold(size_t index)
00134 {
00135   return (*secondaryEnergyCuts)[index];
00136 }
00137 
00138 #endif

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