G4DNAScreenedRutherfordElasticModel.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 #ifndef G4DNAScreenedRutherfordElasticModel_h
00030 #define G4DNAScreenedRutherfordElasticModel_h 1
00031 
00032 #include <CLHEP/Units/SystemOfUnits.h>
00033 
00034 #include "G4VEmModel.hh"
00035 #include "G4ParticleChangeForGamma.hh"
00036 #include "G4ProductionCutsTable.hh"
00037 #include "G4NistManager.hh"
00038 
00039 class G4DNAScreenedRutherfordElasticModel : public G4VEmModel
00040 {
00041 
00042 public:
00043 
00044   G4DNAScreenedRutherfordElasticModel(const G4ParticleDefinition* p = 0, 
00045                           const G4String& nam = "DNAScreenedRutherfordElasticModel");
00046 
00047   virtual ~G4DNAScreenedRutherfordElasticModel();
00048 
00049   virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00050 
00051   virtual G4double CrossSectionPerVolume(const G4Material* material,
00052                                            const G4ParticleDefinition* p,
00053                                            G4double ekin,
00054                                            G4double emin,
00055                                            G4double emax);
00056 
00057   virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00058                                  const G4MaterialCutsCouple*,
00059                                  const G4DynamicParticle*,
00060                                  G4double tmin,
00061                                  G4double maxEnergy);
00062 
00063   inline void SetKillBelowThreshold (G4double threshold);                
00064   G4double GetKillBelowThreshold () { return killBelowEnergy; }          
00065 
00066 protected:
00067 
00068   G4ParticleChangeForGamma* fParticleChangeForGamma;
00069 
00070 private:
00071   // Water density table
00072   const std::vector<G4double>* fpWaterDensity;
00073 
00074   G4double killBelowEnergy;  
00075   G4double lowEnergyLimit;  
00076   G4double intermediateEnergyLimit;
00077   G4double highEnergyLimit; 
00078   G4bool isInitialised;
00079   G4int verboseLevel;
00080   
00081   // Cross section
00082   
00083   G4double RutherfordCrossSection(G4double energy, G4double z);
00084   
00085   G4double ScreeningFactor(G4double energy, G4double z);
00086   
00087   // Final state according to Brenner & Zaider
00088 
00089   G4double BrennerZaiderRandomizeCosTheta(G4double k);
00090   G4double CalculatePolynomial(G4double k, std::vector<G4double>& vec);
00091   std::vector<G4double> betaCoeff;
00092   std::vector<G4double> deltaCoeff;
00093   std::vector<G4double> gamma035_10Coeff;
00094   std::vector<G4double> gamma10_100Coeff;
00095   std::vector<G4double> gamma100_200Coeff;
00096    
00097   // Final state according to Screened Rutherford
00098 
00099   G4double ScreenedRutherfordRandomizeCosTheta(G4double k, G4double z);
00100 
00101   //
00102    
00103   G4DNAScreenedRutherfordElasticModel & operator=(const  G4DNAScreenedRutherfordElasticModel &right);
00104   G4DNAScreenedRutherfordElasticModel(const  G4DNAScreenedRutherfordElasticModel&);
00105 
00106 };
00107 
00108 inline void G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold (G4double threshold) 
00109 { 
00110     killBelowEnergy = threshold; 
00111     if (threshold < 9*CLHEP::eV)
00112      G4Exception ("*** WARNING : the G4DNAScreenedRutherfordElasticModel class is not validated below 9 eV !","",JustWarning,"") ;   
00113 }                
00114 
00115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00116 
00117 #endif

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