G4BetheBlochModel.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 // GEANT4 Class header file
00031 //
00032 //
00033 // File name:     G4BetheBlochModel
00034 //
00035 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
00036 //
00037 // Creation date: 03.01.2002
00038 //
00039 // Modifications:
00040 //
00041 // 23-12-02 V.Ivanchenko change interface in order to moveto cut per region
00042 // 24-01-03 Make models region aware (V.Ivanchenko)
00043 // 13-02-03 Add name (V.Ivanchenko)
00044 // 12-11-03 Fix for GenericIons (V.Ivanchenko)
00045 // 24-03-05 Add G4EmCorrections (V.Ivanchenko)
00046 // 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
00047 // 11-04-04 Move MaxSecondaryEnergy to models (V.Ivanchenko)
00048 // 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
00049 // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio, 
00050 //          CorrectionsAlongStep needed for ions(V.Ivanchenko)
00051 
00052 //
00053 // Class Description:
00054 //
00055 // Implementation of Bethe-Bloch model of energy loss and
00056 // delta-electron production by heavy charged particles
00057 
00058 // -------------------------------------------------------------------
00059 //
00060 
00061 #ifndef G4BetheBlochModel_h
00062 #define G4BetheBlochModel_h 1
00063 
00064 #include <CLHEP/Units/SystemOfUnits.h>
00065 
00066 #include "G4VEmModel.hh"
00067 #include "G4NistManager.hh"
00068 
00069 class G4EmCorrections;
00070 class G4ParticleChangeForLoss;
00071 
00072 class G4BetheBlochModel : public G4VEmModel
00073 {
00074 
00075 public:
00076 
00077   G4BetheBlochModel(const G4ParticleDefinition* p = 0,
00078                     const G4String& nam = "BetheBloch");
00079 
00080   virtual ~G4BetheBlochModel();
00081 
00082   virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00083 
00084   virtual G4double ComputeCrossSectionPerElectron(
00085                                  const G4ParticleDefinition*,
00086                                  G4double kineticEnergy,
00087                                  G4double cutEnergy,
00088                                  G4double maxEnergy);
00089                                  
00090   virtual G4double ComputeCrossSectionPerAtom(
00091                                  const G4ParticleDefinition*,
00092                                  G4double kineticEnergy,
00093                                  G4double Z, G4double A,
00094                                  G4double cutEnergy,
00095                                  G4double maxEnergy);
00096                                                                  
00097   virtual G4double CrossSectionPerVolume(const G4Material*,
00098                                  const G4ParticleDefinition*,
00099                                  G4double kineticEnergy,
00100                                  G4double cutEnergy,
00101                                  G4double maxEnergy);
00102                                  
00103   virtual G4double ComputeDEDXPerVolume(const G4Material*,
00104                                         const G4ParticleDefinition*,
00105                                         G4double kineticEnergy,
00106                                         G4double cutEnergy);
00107 
00108   virtual G4double GetChargeSquareRatio(const G4ParticleDefinition* p,
00109                                         const G4Material* mat,
00110                                         G4double kineticEnergy);
00111 
00112   virtual G4double GetParticleCharge(const G4ParticleDefinition* p,
00113                                      const G4Material* mat,
00114                                      G4double kineticEnergy);
00115 
00116   virtual void CorrectionsAlongStep(const G4MaterialCutsCouple* couple,
00117                                     const G4DynamicParticle* dp,
00118                                     G4double& eloss,
00119                                     G4double&,
00120                                     G4double length);
00121 
00122   virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00123                                  const G4MaterialCutsCouple*,
00124                                  const G4DynamicParticle*,
00125                                  G4double tmin,
00126                                  G4double maxEnergy);
00127 
00128 protected:
00129 
00130   virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
00131                                       G4double kinEnergy);
00132 
00133   inline G4double GetChargeSquareRatio() const;
00134 
00135   inline void SetChargeSquareRatio(G4double val);
00136 
00137 private:
00138 
00139   void SetupParameters();
00140 
00141   inline void SetParticle(const G4ParticleDefinition* p);
00142 
00143   inline void SetGenericIon(const G4ParticleDefinition* p);
00144 
00145   // hide assignment operator
00146   G4BetheBlochModel & operator=(const  G4BetheBlochModel &right);
00147   G4BetheBlochModel(const  G4BetheBlochModel&);
00148 
00149   const G4ParticleDefinition* particle;
00150   G4ParticleDefinition*       theElectron;
00151   G4EmCorrections*            corr;
00152   G4ParticleChangeForLoss*    fParticleChange;
00153   G4NistManager*              nist;
00154 
00155   G4double mass;
00156   G4double tlimit;
00157   G4double spin;
00158   G4double magMoment2;
00159   G4double chargeSquare;
00160   G4double ratio;
00161   G4double formfact;
00162   G4double twoln10;
00163   G4double bg2lim;
00164   G4double taulim;
00165   G4double corrFactor;
00166   G4bool   isIon;
00167   G4bool   isInitialised;
00168 };
00169 
00170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00171 
00172 inline void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p)
00173 {
00174   if(particle != p) {
00175     particle = p;
00176     if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > 1.5*CLHEP::eplus) 
00177       { isIon = true; }
00178     SetupParameters();
00179   }
00180 }
00181 
00182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00183 
00184 inline void G4BetheBlochModel::SetGenericIon(const G4ParticleDefinition* p)
00185 {
00186   if(p && particle != p) { 
00187     if(p->GetParticleName() == "GenericIon") { isIon = true; }
00188   }
00189 }
00190 
00191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00192 
00193 inline G4double G4BetheBlochModel::GetChargeSquareRatio() const
00194 {
00195   return chargeSquare;
00196 }
00197 
00198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00199 
00200 inline void G4BetheBlochModel::SetChargeSquareRatio(G4double val)
00201 {
00202   chargeSquare = val;
00203 }
00204 
00205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00206 
00207 #endif

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