G4hImpactIonisation.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 //
00027 // ------------------------------------------------------------
00028 // G4hImpactIonisation
00029 //
00030 // $Id$
00031 //
00032 // Author: Maria Grazia Pia (MariaGrazia.Pia@ge.infn.it)
00033 //
00034 // 08 Sep 2008 - MGP - Created (initially based on G4hLowEnergyIonisation) 
00035 //                     Added PIXE capabilities
00036 //                     Partial clean-up of the implementation (more needed)
00037 //                     Calculation of MicroscopicCrossSection delegated to specialised class 
00038 //
00039 // ------------------------------------------------------------
00040  
00041 // Class Description:
00042 // Impact Ionisation process of charged hadrons and ions
00043 // Initially based on G4hLowEnergyIonisation, to be subject to redesign
00044 // and further evolution of physics capabilities
00045 //
00046 // The physics model of G4hLowEnergyIonisation is described in 
00047 // CERN-OPEN-99-121 and CERN-OPEN-99-300. 
00048 //
00049 // Documentation available in:
00050 // M.G. Pia et al., PIXE Simulation With Geant4,
00051 // IEEE Trans. Nucl. Sci., vol. 56, no. 6, pp. 3614-3649, Dec. 2009.
00052 
00053 // ------------------------------------------------------------
00054 
00055 #ifndef G4HIMPACTIONISATION
00056 #define G4HIMPACTIONISATION 1
00057 
00058 #include <map>
00059 #include <CLHEP/Units/PhysicalConstants.h>
00060 
00061 #include "globals.hh"
00062 #include "G4hRDEnergyLoss.hh"
00063 #include "G4DataVector.hh"
00064 #include "G4AtomicDeexcitation.hh"
00065 #include "G4PixeCrossSectionHandler.hh"
00066 
00067 class G4VLowEnergyModel;
00068 class G4VParticleChange;
00069 class G4ParticleDefinition;
00070 class G4PhysicsTable;
00071 class G4MaterialCutsCouple;
00072 class G4Track;
00073 class G4Step;
00074 
00075 class G4hImpactIonisation : public G4hRDEnergyLoss
00076 {
00077 public: // With description
00078   
00079   G4hImpactIonisation(const G4String& processName = "hImpactIoni"); 
00080   // The ionisation process for hadrons/ions to be include in the
00081   // UserPhysicsList
00082 
00083   ~G4hImpactIonisation();
00084   // Destructor
00085   
00086   G4bool IsApplicable(const G4ParticleDefinition&);
00087   // True for all charged hadrons/ions
00088     
00089   void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ;
00090   // Build physics table during initialisation
00091 
00092   G4double GetMeanFreePath(const G4Track& track,
00093                            G4double previousStepSize,
00094                            enum G4ForceCondition* condition );
00095   // Return MeanFreePath until delta-electron production
00096   
00097   void PrintInfoDefinition() const;
00098   // Print out of the class parameters
00099 
00100   void SetHighEnergyForProtonParametrisation(G4double energy) {protonHighEnergy = energy;} ;
00101   // Definition of the boundary proton energy. For higher energies
00102   // Bethe-Bloch formula is used, for lower energies a parametrisation
00103   // of the energy losses is performed. Default is 2 MeV.
00104 
00105   void SetLowEnergyForProtonParametrisation(G4double energy) {protonLowEnergy = energy;} ;
00106   // Set of the boundary proton energy. For lower energies
00107   // the Free Electron Gas model is used for the energy losses.
00108   // Default is 1 keV.
00109 
00110   void SetHighEnergyForAntiProtonParametrisation(G4double energy) {antiprotonHighEnergy = energy;} ;
00111   // Set of the boundary antiproton energy. For higher energies
00112   // Bethe-Bloch formula is used, for lower energies parametrisation
00113   // of the energy losses is performed. Default is 2 MeV.
00114 
00115   void SetLowEnergyForAntiProtonParametrisation(G4double energy) {antiprotonLowEnergy = energy;} ;
00116   // Set of the boundary antiproton energy. For lower energies
00117   // the Free Electron Gas model is used for the energy losses. 
00118   // Default is 1 keV.
00119 
00120   G4double GetContinuousStepLimit(const G4Track& track,
00121                                   G4double previousStepSize,
00122                                   G4double currentMinimumStep,
00123                                   G4double& currentSafety); 
00124   // Calculation of the step limit due to ionisation losses
00125 
00126   void SetElectronicStoppingPowerModel(const G4ParticleDefinition* aParticle,
00127                                        const G4String& dedxTable);
00128   // This method defines the electron ionisation parametrisation method 
00129   // via the name of the table. Default is "ICRU_49p".
00130 
00131   void SetNuclearStoppingPowerModel(const G4String& dedxTable)
00132   {theNuclearTable = dedxTable; SetNuclearStoppingOn();};
00133   // This method defines the nuclear ionisation parametrisation method 
00134   // via the name of the table. Default is "ICRU_49".
00135 
00136   // ---- MGP ---- The following design of On/Off is nonsense; to be modified
00137   // in a following design iteration
00138 
00139   void SetNuclearStoppingOn() {nStopping = true;};
00140   // This method switch on calculation of the nuclear stopping power.
00141   
00142   void SetNuclearStoppingOff() {nStopping = false;};
00143   // This method switch off calculation of the nuclear stopping power.
00144 
00145   void SetBarkasOn() {theBarkas = true;};
00146   // This method switch on calculation of the Barkas and Bloch effects.
00147 
00148   void SetBarkasOff() {theBarkas = false;};
00149   // This method switch off calculation of the Barkas and Bloch effects.
00150 
00151   void SetPixe(const G4bool /* val */ ) {pixeIsActive = true;};
00152   // This method switches atomic relaxation on/off; currently always on
00153 
00154   G4VParticleChange* AlongStepDoIt(const G4Track& trackData ,
00155                                    const G4Step& stepData ) ;
00156   // Function to determine total energy deposition on the step
00157 
00158   G4VParticleChange* PostStepDoIt(const G4Track& track,
00159                                   const G4Step& Step  ) ;
00160   // Simulation of delta-ray production.
00161 
00162   G4double ComputeDEDX(const G4ParticleDefinition* aParticle,
00163                        const G4MaterialCutsCouple* couple,
00164                        G4double kineticEnergy);
00165   // This method returns electronic dE/dx for protons or antiproton
00166 
00167   void SetCutForSecondaryPhotons(G4double cut);
00168   // Set threshold energy for fluorescence
00169 
00170   void SetCutForAugerElectrons(G4double cut);
00171   // Set threshold energy for Auger electron production
00172 
00173   void ActivateAugerElectronProduction(G4bool val);
00174   // Set Auger electron production flag on/off
00175 
00176   // Accessors to configure PIXE
00177   void SetPixeCrossSectionK(const G4String& name) { modelK = name; }
00178   void SetPixeCrossSectionL(const G4String& name) { modelL = name; }
00179   void SetPixeCrossSectionM(const G4String& name) { modelM = name; }
00180   void SetPixeProjectileMinEnergy(G4double energy) { eMinPixe = energy; }
00181   void SetPixeProjectileMaxEnergy(G4double energy) { eMaxPixe = energy; }
00182 
00183 protected:
00184 
00185 private:
00186 
00187   void InitializeMe();
00188   void InitializeParametrisation();
00189   void BuildLossTable(const G4ParticleDefinition& aParticleType);
00190   // void BuildDataForFluorescence(const G4ParticleDefinition& aParticleType);
00191   void BuildLambdaTable(const G4ParticleDefinition& aParticleType);
00192   void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable)
00193   {protonTable = dedxTable ;};
00194   // This method defines the ionisation parametrisation method via its name
00195 
00196   void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable)
00197   {antiprotonTable = dedxTable;};
00198 
00199   G4double MicroscopicCrossSection(const G4ParticleDefinition& aParticleType,
00200                                    G4double kineticEnergy,
00201                                    G4double atomicNumber,
00202                                    G4double deltaCutInEnergy) const;
00203 
00204   G4double GetConstraints(const G4DynamicParticle* particle,
00205                           const G4MaterialCutsCouple* couple);
00206   // Function to determine StepLimit
00207 
00208   G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
00209                                   G4double kineticEnergy) const;
00210 
00211   G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
00212                                       G4double kineticEnergy) const;
00213 
00214   G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple,
00215                            G4double kineticEnergy,
00216                            G4double particleMass) const;
00217   // This method returns average energy loss due to delta-rays emission with
00218   // energy higher than the cut energy for given material.
00219 
00220   G4double BarkasTerm(const G4Material* material,
00221                       G4double kineticEnergy) const;
00222   // Function to compute the Barkas term for protons
00223 
00224   G4double BlochTerm(const G4Material* material,
00225                      G4double kineticEnergy,
00226                      G4double cSquare) const;
00227   // Function to compute the Bloch term for protons
00228 
00229   G4double ElectronicLossFluctuation(const G4DynamicParticle* particle,
00230                                      const G4MaterialCutsCouple* material,
00231                                      G4double meanLoss,
00232                                      G4double step) const;
00233   // Function to sample electronic losses
00234 
00235   // hide assignment operator
00236   G4hImpactIonisation & operator=(const G4hImpactIonisation &right);
00237   G4hImpactIonisation(const G4hImpactIonisation&);
00238 
00239 private:
00240   //  private data members ...............................
00241   G4VLowEnergyModel* betheBlochModel;
00242   G4VLowEnergyModel* protonModel;
00243   G4VLowEnergyModel* antiprotonModel;
00244   G4VLowEnergyModel* theIonEffChargeModel;
00245   G4VLowEnergyModel* theNuclearStoppingModel;
00246   G4VLowEnergyModel* theIonChuFluctuationModel;
00247   G4VLowEnergyModel* theIonYangFluctuationModel;
00248 
00249   // std::map<G4int,G4double,std::less<G4int> > totalCrossSectionMap;
00250 
00251   // name of parametrisation table of electron stopping power
00252   G4String protonTable;
00253   G4String antiprotonTable;
00254   G4String theNuclearTable;
00255 
00256   // interval of parametrisation of electron stopping power
00257   G4double protonLowEnergy;
00258   G4double protonHighEnergy;
00259   G4double antiprotonLowEnergy;
00260   G4double antiprotonHighEnergy;
00261 
00262   // flag of parametrisation of nucleus stopping power
00263   G4bool nStopping;
00264   G4bool theBarkas;
00265 
00266   G4DataVector cutForDelta;
00267   G4DataVector cutForGamma;
00268   G4double minGammaEnergy;
00269   G4double minElectronEnergy;
00270   G4PhysicsTable* theMeanFreePathTable;
00271 
00272   const G4double paramStepLimit; // parameter limits the step at low energy
00273 
00274   G4double fdEdx;        // computed in GetContraints
00275   G4double fRangeNow ;   //
00276   G4double charge;       //
00277   G4double chargeSquare; //
00278   G4double initialMass;  // mass to calculate Lambda tables
00279   G4double fBarkas;
00280 
00281   G4PixeCrossSectionHandler* pixeCrossSectionHandler;
00282   G4AtomicDeexcitation atomicDeexcitation;
00283   G4String modelK;
00284   G4String modelL;
00285   G4String modelM;
00286   G4double eMinPixe;
00287   G4double eMaxPixe;
00288  
00289   G4bool pixeIsActive;
00290 
00291 };
00292 
00293 
00294 inline G4double G4hImpactIonisation::GetContinuousStepLimit(const G4Track& track,
00295                                                             G4double,
00296                                                             G4double currentMinimumStep,
00297                                                             G4double&)
00298 {
00299   G4double step = GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
00300 
00301   // ---- MGP ---- The following line, taken as is from G4hLowEnergyIonisation,
00302   // is meaningless: currentMinimumStep is passed by value,
00303   // therefore any local modification to it has no effect
00304 
00305   if ((step > 0.) && (step < currentMinimumStep)) currentMinimumStep = step ;
00306 
00307   return step ;
00308 }
00309 
00310 
00311 inline G4bool G4hImpactIonisation::IsApplicable(const G4ParticleDefinition& particle)
00312 {
00313   // ---- MGP ---- Better criterion for applicability to be defined;
00314   // now hard-coded particle mass > 0.1 * proton_mass
00315 
00316   return (particle.GetPDGCharge() != 0.0 && particle.GetPDGMass() > CLHEP::proton_mass_c2*0.1);
00317 }
00318 
00319 #endif
00320 
00321 
00322 
00323 
00324 
00325 
00326 
00327 

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