G4VLowEnergyDiscretePhotonProcess.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // --------------------------------------------------------------
00029 //
00030 // File name:     G4VLowEnergyDiscretePhotonProcess.hh
00031 //
00032 // Author:        Capra Riccardo
00033 //
00034 // Creation date: May 2005
00035 //
00036 // History:
00037 // -----------
00038 // 02 May 2005  R. Capra         1st implementation
00039 // 20 May 2005  MGP              Modified protected member functions to return a const pointer
00040 //
00041 //----------------------------------------------------------------
00042 
00043 
00044 
00045 #ifndef   G4VLowEnergyDiscretePhotonProcess_hh
00046 #define  G4VLowEnergyDiscretePhotonProcess_hh
00047  
00048 // Base class
00049 #include "G4VLowEnergyTestableDiscreteProcess.hh"
00050 
00051 // Forward declaration
00052 class G4String;
00053 class G4ParticleDefinition;
00054 class G4DynamicParticle;
00055 class G4Track;
00056 class G4VCrossSectionHandler;
00057 class G4VEMDataSet;
00058 class G4VDataSetAlgorithm;
00059  
00060 // G4VLowEnergyDiscretePhotonProcess
00061 // A common class for Rayleigh and Compton processes
00062 class G4VLowEnergyDiscretePhotonProcess : public G4VLowEnergyTestableDiscreteProcess
00063 {
00064 public:
00065   //   Class constructor
00066   //
00067   //   Creates crossSectionHandler and scatterFunctionData
00068   //   scatterFunctionData are loaded from the scatterFile file, and are interpolated with scatterInterpolation algorithm
00069   //   processName The name of the process
00070   //   aCrossSectionFileName The name of the cross-section data file
00071   //   aScatterFileName The name of the scatter function data
00072   //   aScatterInterpolation The interpolation algorithm
00073   //   aLowEnergyLimit The lower energy limit of validity against data
00074   //   aHighEnergyLimit The higher energy limit of validity against data
00075   G4VLowEnergyDiscretePhotonProcess(const G4String& processName, 
00076                                     const G4String& aCrossSectionFileName, 
00077                                     const G4String& aScatterFileName, 
00078                                     G4VDataSetAlgorithm* aScatterInterpolation, 
00079                                     G4double aLowEnergyLimit, 
00080                                     G4double aHighEnergyLimit);
00081 
00082   //   Class destructor
00083   //
00084   //        Deletes crossSectionHandler, meanFreePathTable and scatterFunctionData
00085   virtual ~G4VLowEnergyDiscretePhotonProcess(void);
00086 
00087 
00088 
00089   //   Checks if the process is applicable to the particle
00090   //
00091   //   For processes inheriting from this class the only applicable particle is the photon
00092   //   particleDefinition Is the particle to be checked for the applicability of the process
00093   //   true only if the particle is a photon
00094   virtual G4bool IsApplicable(const G4ParticleDefinition& particleDefinition);
00095 
00096   //   Updates the crossSectionHandler and builds the meanFreePathTable from it
00097   //
00098   //        crossSectionHandler data is loaded from crossSectionFile file
00099   //   photon particle is always a photon
00100   virtual void BuildPhysicsTable(const G4ParticleDefinition& photon);
00101    
00102    
00103    
00104   //   Returns the low energy limit of the process
00105   //   The low energy limit of the process
00106   inline G4double GetLowEnergyLimit(void) const {return lowEnergyLimit;}
00107    
00108   //   Returns the high energy limit of the process
00109   //   The high energy limit of the process
00110   inline G4double GetHighEnergyLimit(void) const {return highEnergyLimit;}
00111 
00112 
00113 
00114 protected:
00115   //   Evaluates the process mean free path
00116   //
00117   //    Mean free path evaluation is based on meanFreePathTable generated by the crossSectionHandler with data taken from the crossSectionFile.
00118   //   The method uses lowEnergyLimit and highEnergyLimit
00119   //   aTrack the particle momentum for which the mean free path must be evaluates (a photon)
00120   //   previousStepSize the size of the prevous step (not used by this implementation)
00121   //   condition the consition to be updated (not used by this implementation)
00122   //   The mean free path for the process
00123   virtual G4double GetMeanFreePath(const G4Track &aTrack, 
00124                                    G4double previousStepSize, 
00125                                    G4ForceCondition* condition);
00126    
00127   //   Returns the cross-section handler
00128   //   The cross-section handler
00129   inline const G4VCrossSectionHandler* GetCrossSectionHandler(void) const {return crossSectionHandler;}
00130    
00131   //   Returns the mean free path table
00132   //   The mean free path table
00133   inline const G4VEMDataSet* GetMeanFreePathTable(void) const {return meanFreePathTable;}
00134    
00135   //   Returns the scatter function data
00136   //   The scatter function data
00137   inline const G4VEMDataSet* GetScatterFunctionData(void) const {return scatterFunctionData;}
00138 
00139 
00140 
00141   //   Verifies if the polarization vector is orthonormal to the photon direction
00142   //
00143   //        This method is used by polarized processes. It returns always a vector orthonormal to the photon direction.
00144   //        When G4DynamicParticle::GetPolarization() is well defined the returned vector matches this vector.
00145   //        When G4DynamicParticle::GetPolarization() is sensibly different from a null-vector, the returned vector is
00146   //        orthonormal to the photon direction and on the plane made with G4DynamicParticle::GetPolarization()
00147   //        When G4DynamicParticle::GetPolarization() is almost a null-vector, the returned vector is a random
00148   //        vector orthonormal to the photon direction
00149   //   photon The incoming photon
00150   //  s Always a vector ortogonal to the photon direction.
00151   G4ThreeVector GetPhotonPolarization(const G4DynamicParticle& photon);
00152 
00153 
00154 
00155 private:
00156 
00157   //   Hides copy constructor
00158   G4VLowEnergyDiscretePhotonProcess(const G4VLowEnergyDiscretePhotonProcess& );
00159 
00160   //   Hides assignment operator
00161   G4VLowEnergyDiscretePhotonProcess& operator=(const G4VLowEnergyDiscretePhotonProcess& );
00162 
00163   //   Low energy limit
00164   G4double lowEnergyLimit;
00165 
00166   //   High energy limit
00167   G4double highEnergyLimit;
00168 
00169   //   Cross-section file name
00170   G4String crossSectionFileName;
00171 
00172   //   Cross-section handler
00173   G4VCrossSectionHandler* crossSectionHandler;
00174 
00175   //   The mean free path table based on cross-section data
00176   G4VEMDataSet* meanFreePathTable;
00177 
00178   //   Scatter function table
00179   G4VEMDataSet* scatterFunctionData;
00180 };
00181 
00182 #endif /* G4VLowEnergyDiscretePhotonProcess_hh */

Generated on Mon May 27 17:50:17 2013 for Geant4 by  doxygen 1.4.7