Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:24

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // ------------------------------------------------------------
0028 // G4hImpactIonisation
0029 //
0030 //
0031 // Author: Maria Grazia Pia (MariaGrazia.Pia@ge.infn.it)
0032 //
0033 // 08 Sep 2008 - MGP - Created (initially based on G4hLowEnergyIonisation) 
0034 //                     Added PIXE capabilities
0035 //                     Partial clean-up of the implementation (more needed)
0036 //                     Calculation of MicroscopicCrossSection delegated to specialised class 
0037 //
0038 // ------------------------------------------------------------
0039  
0040 // Class Description:
0041 // Impact Ionisation process of charged hadrons and ions
0042 // Initially based on G4hLowEnergyIonisation, to be subject to redesign
0043 // and further evolution of physics capabilities
0044 //
0045 // The physics model of G4hLowEnergyIonisation is described in 
0046 // CERN-OPEN-99-121 and CERN-OPEN-99-300. 
0047 //
0048 // Documentation available in:
0049 // M.G. Pia et al., PIXE Simulation With Geant4,
0050 // IEEE Trans. Nucl. Sci., vol. 56, no. 6, pp. 3614-3649, Dec. 2009.
0051 
0052 // ------------------------------------------------------------
0053 
0054 #ifndef G4HIMPACTIONISATION
0055 #define G4HIMPACTIONISATION 1
0056 
0057 #include <map>
0058 #include <CLHEP/Units/PhysicalConstants.h>
0059 
0060 #include "globals.hh"
0061 #include "G4hRDEnergyLoss.hh"
0062 #include "G4DataVector.hh"
0063 #include "G4AtomicDeexcitation.hh"
0064 #include "G4PixeCrossSectionHandler.hh"
0065 
0066 class G4VLowEnergyModel;
0067 class G4VParticleChange;
0068 class G4ParticleDefinition;
0069 class G4PhysicsTable;
0070 class G4MaterialCutsCouple;
0071 class G4Track;
0072 class G4Step;
0073 
0074 class G4hImpactIonisation : public G4hRDEnergyLoss
0075 {
0076 public: // With description
0077   
0078   G4hImpactIonisation(const G4String& processName = "hImpactIoni"); 
0079   // The ionisation process for hadrons/ions to be include in the
0080   // UserPhysicsList
0081 
0082   ~G4hImpactIonisation();
0083   // Destructor
0084   
0085   G4bool IsApplicable(const G4ParticleDefinition&);
0086   // True for all charged hadrons/ions
0087     
0088   void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ;
0089   // Build physics table during initialisation
0090 
0091   G4double GetMeanFreePath(const G4Track& track,
0092                G4double previousStepSize,
0093                enum G4ForceCondition* condition );
0094   // Return MeanFreePath until delta-electron production
0095   
0096   void PrintInfoDefinition() const;
0097   // Print out of the class parameters
0098 
0099   void SetHighEnergyForProtonParametrisation(G4double energy) {protonHighEnergy = energy;} ;
0100   // Definition of the boundary proton energy. For higher energies
0101   // Bethe-Bloch formula is used, for lower energies a parametrisation
0102   // of the energy losses is performed. Default is 2 MeV.
0103 
0104   void SetLowEnergyForProtonParametrisation(G4double energy) {protonLowEnergy = energy;} ;
0105   // Set of the boundary proton energy. For lower energies
0106   // the Free Electron Gas model is used for the energy losses.
0107   // Default is 1 keV.
0108 
0109   void SetHighEnergyForAntiProtonParametrisation(G4double energy) {antiprotonHighEnergy = energy;} ;
0110   // Set of the boundary antiproton energy. For higher energies
0111   // Bethe-Bloch formula is used, for lower energies parametrisation
0112   // of the energy losses is performed. Default is 2 MeV.
0113 
0114   void SetLowEnergyForAntiProtonParametrisation(G4double energy) {antiprotonLowEnergy = energy;} ;
0115   // Set of the boundary antiproton energy. For lower energies
0116   // the Free Electron Gas model is used for the energy losses. 
0117   // Default is 1 keV.
0118 
0119   G4double GetContinuousStepLimit(const G4Track& track,
0120                   G4double previousStepSize,
0121                   G4double currentMinimumStep,
0122                   G4double& currentSafety); 
0123   // Calculation of the step limit due to ionisation losses
0124 
0125   void SetElectronicStoppingPowerModel(const G4ParticleDefinition* aParticle,
0126                                        const G4String& dedxTable);
0127   // This method defines the electron ionisation parametrisation method 
0128   // via the name of the table. Default is "ICRU_49p".
0129 
0130   void SetNuclearStoppingPowerModel(const G4String& dedxTable)
0131   {theNuclearTable = dedxTable; SetNuclearStoppingOn();};
0132   // This method defines the nuclear ionisation parametrisation method 
0133   // via the name of the table. Default is "ICRU_49".
0134 
0135   // ---- MGP ---- The following design of On/Off is nonsense; to be modified
0136   // in a following design iteration
0137 
0138   void SetNuclearStoppingOn() {nStopping = true;};
0139   // This method switch on calculation of the nuclear stopping power.
0140   
0141   void SetNuclearStoppingOff() {nStopping = false;};
0142   // This method switch off calculation of the nuclear stopping power.
0143 
0144   void SetBarkasOn() {theBarkas = true;};
0145   // This method switch on calculation of the Barkas and Bloch effects.
0146 
0147   void SetBarkasOff() {theBarkas = false;};
0148   // This method switch off calculation of the Barkas and Bloch effects.
0149 
0150   void SetPixe(const G4bool /* val */ ) {pixeIsActive = true;};
0151   // This method switches atomic relaxation on/off; currently always on
0152 
0153   G4VParticleChange* AlongStepDoIt(const G4Track& trackData ,
0154                                    const G4Step& stepData ) ;
0155   // Function to determine total energy deposition on the step
0156 
0157   G4VParticleChange* PostStepDoIt(const G4Track& track,
0158                   const G4Step& Step  ) ;
0159   // Simulation of delta-ray production.
0160 
0161   G4double ComputeDEDX(const G4ParticleDefinition* aParticle,
0162                        const G4MaterialCutsCouple* couple,
0163                G4double kineticEnergy);
0164   // This method returns electronic dE/dx for protons or antiproton
0165 
0166   void SetCutForSecondaryPhotons(G4double cut);
0167   // Set threshold energy for fluorescence
0168 
0169   void SetCutForAugerElectrons(G4double cut);
0170   // Set threshold energy for Auger electron production
0171 
0172   void ActivateAugerElectronProduction(G4bool val);
0173   // Set Auger electron production flag on/off
0174 
0175   // Accessors to configure PIXE
0176   void SetPixeCrossSectionK(const G4String& name) { modelK = name; }
0177   void SetPixeCrossSectionL(const G4String& name) { modelL = name; }
0178   void SetPixeCrossSectionM(const G4String& name) { modelM = name; }
0179   void SetPixeProjectileMinEnergy(G4double energy) { eMinPixe = energy; }
0180   void SetPixeProjectileMaxEnergy(G4double energy) { eMaxPixe = energy; }
0181 
0182 protected:
0183 
0184 private:
0185 
0186   void InitializeMe();
0187   void InitializeParametrisation();
0188   void BuildLossTable(const G4ParticleDefinition& aParticleType);
0189   // void BuildDataForFluorescence(const G4ParticleDefinition& aParticleType);
0190   void BuildLambdaTable(const G4ParticleDefinition& aParticleType);
0191   void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable)
0192   {protonTable = dedxTable ;};
0193   // This method defines the ionisation parametrisation method via its name
0194 
0195   void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable)
0196   {antiprotonTable = dedxTable;};
0197 
0198   G4double MicroscopicCrossSection(const G4ParticleDefinition& aParticleType,
0199                    G4double kineticEnergy,
0200                    G4double atomicNumber,
0201                    G4double deltaCutInEnergy) const;
0202 
0203   G4double GetConstraints(const G4DynamicParticle* particle,
0204                           const G4MaterialCutsCouple* couple);
0205   // Function to determine StepLimit
0206 
0207   G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
0208                   G4double kineticEnergy) const;
0209 
0210   G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
0211                       G4double kineticEnergy) const;
0212 
0213   G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple,
0214                G4double kineticEnergy,
0215                G4double particleMass) const;
0216   // This method returns average energy loss due to delta-rays emission with
0217   // energy higher than the cut energy for given material.
0218 
0219   G4double BarkasTerm(const G4Material* material,
0220               G4double kineticEnergy) const;
0221   // Function to compute the Barkas term for protons
0222 
0223   G4double BlochTerm(const G4Material* material,
0224              G4double kineticEnergy,
0225              G4double cSquare) const;
0226   // Function to compute the Bloch term for protons
0227 
0228   G4double ElectronicLossFluctuation(const G4DynamicParticle* particle,
0229                                      const G4MaterialCutsCouple* material,
0230                      G4double meanLoss,
0231                      G4double step) const;
0232   // Function to sample electronic losses
0233 
0234   // hide assignment operator
0235   G4hImpactIonisation & operator=(const G4hImpactIonisation &right);
0236   G4hImpactIonisation(const G4hImpactIonisation&);
0237 
0238 private:
0239   //  private data members ...............................
0240   G4VLowEnergyModel* betheBlochModel;
0241   G4VLowEnergyModel* protonModel;
0242   G4VLowEnergyModel* antiprotonModel;
0243   G4VLowEnergyModel* theIonEffChargeModel;
0244   G4VLowEnergyModel* theNuclearStoppingModel;
0245   G4VLowEnergyModel* theIonChuFluctuationModel;
0246   G4VLowEnergyModel* theIonYangFluctuationModel;
0247 
0248   // std::map<G4int,G4double,std::less<G4int> > totalCrossSectionMap;
0249 
0250   // name of parametrisation table of electron stopping power
0251   G4String protonTable;
0252   G4String antiprotonTable;
0253   G4String theNuclearTable;
0254 
0255   // interval of parametrisation of electron stopping power
0256   G4double protonLowEnergy;
0257   G4double protonHighEnergy;
0258   G4double antiprotonLowEnergy;
0259   G4double antiprotonHighEnergy;
0260 
0261   // flag of parametrisation of nucleus stopping power
0262   G4bool nStopping;
0263   G4bool theBarkas;
0264 
0265   G4DataVector cutForDelta;
0266   G4DataVector cutForGamma;
0267   G4double minGammaEnergy;
0268   G4double minElectronEnergy;
0269   G4PhysicsTable* theMeanFreePathTable;
0270 
0271   const G4double paramStepLimit; // parameter limits the step at low energy
0272 
0273   G4double fdEdx;        // computed in GetContraints
0274   G4double fRangeNow ;   //
0275   G4double charge;       //
0276   G4double chargeSquare; //
0277   G4double initialMass;  // mass to calculate Lambda tables
0278   G4double fBarkas;
0279 
0280   G4PixeCrossSectionHandler* pixeCrossSectionHandler;
0281   G4AtomicDeexcitation atomicDeexcitation;
0282   G4String modelK;
0283   G4String modelL;
0284   G4String modelM;
0285   G4double eMinPixe;
0286   G4double eMaxPixe;
0287  
0288   G4bool pixeIsActive;
0289 
0290 };
0291 
0292 
0293 inline G4double G4hImpactIonisation::GetContinuousStepLimit(const G4Track& track,
0294                                 G4double,
0295                                 G4double currentMinimumStep,
0296                                 G4double&)
0297 {
0298   G4double step = GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
0299 
0300   // ---- MGP ---- The following line, taken as is from G4hLowEnergyIonisation,
0301   // is meaningless: currentMinimumStep is passed by value,
0302   // therefore any local modification to it has no effect
0303 
0304   if ((step > 0.) && (step < currentMinimumStep)) currentMinimumStep = step ;
0305 
0306   return step ;
0307 }
0308 
0309 
0310 inline G4bool G4hImpactIonisation::IsApplicable(const G4ParticleDefinition& particle)
0311 {
0312   // ---- MGP ---- Better criterion for applicability to be defined;
0313   // now hard-coded particle mass > 0.1 * proton_mass
0314 
0315   return (particle.GetPDGCharge() != 0.0 && particle.GetPDGMass() > CLHEP::proton_mass_c2*0.1);
0316 }
0317 
0318 #endif
0319 
0320 
0321 
0322 
0323 
0324 
0325 
0326