Back to home page

EIC code displayed by LXR

 
 

    


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

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 // Calculation of the nucleus-nucleus total, inelastic, production, 
0027 // elastic and quasi-elastic  cross-sections
0028 // based on parametrisations of nucleon-nucleon
0029 // cross-sections  in 
0030 // the framework of simplified Glauber-Gribov approach
0031 //
0032 //
0033 // 24.11.08 V. Grichine - first implementation based on 
0034 //                        G4GlauberGribovCrossSection
0035 //
0036 // 04.09.18 V. Ivantchenko Major revision of interfaces and implementation
0037 // 27.05.19 V. Ivantchenko Removed obsolete methods and members 
0038 //
0039 
0040 #ifndef G4ComponentGGNuclNuclXsc_h
0041 #define G4ComponentGGNuclNuclXsc_h
0042 
0043 #include "globals.hh"
0044 #include "G4VComponentCrossSection.hh"
0045 #include "G4DynamicParticle.hh"
0046 
0047 class G4ParticleDefinition;
0048 class G4HadronNucleonXsc;
0049 class G4ComponentGGHadronNucleusXsc;
0050 class G4Material;
0051 
0052 class G4ComponentGGNuclNuclXsc : public G4VComponentCrossSection
0053 {
0054 public:
0055 
0056   G4ComponentGGNuclNuclXsc ();
0057   virtual ~G4ComponentGGNuclNuclXsc ();
0058 
0059   // virtual interface methods
0060   G4double GetTotalElementCrossSection(const G4ParticleDefinition* aParticle,
0061                        G4double kinEnergy, 
0062                        G4int Z, G4double A) final;
0063 
0064   G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0065                        G4double kinEnergy,
0066                        G4int Z, G4int A) final;
0067 
0068   G4double GetInelasticElementCrossSection(const G4ParticleDefinition* aParticle,
0069                        G4double kinEnergy, 
0070                        G4int Z, G4double A) final;
0071 
0072   G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0073                        G4double kinEnergy, 
0074                        G4int Z, G4int A) final;
0075 
0076   G4double GetElasticElementCrossSection(const G4ParticleDefinition* aParticle,
0077                      G4double kinEnergy, 
0078                      G4int Z, G4double A) final;
0079 
0080   G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition* aParticle,
0081                      G4double kinEnergy, 
0082                      G4int Z, G4int A) final;
0083  
0084   G4double ComputeQuasiElasticRatio(const G4ParticleDefinition* aParticle,
0085                     G4double kinEnergy, 
0086                     G4int Z, G4int A) final;
0087 
0088   void BuildPhysicsTable(const G4ParticleDefinition&) final;
0089 
0090   void DumpPhysicsTable(const G4ParticleDefinition&) final;
0091 
0092   void Description(std::ostream&) const final;
0093    
0094   // Extra methods
0095   //  inline G4double GetElementCrossSection(const G4DynamicParticle*, 
0096   //                         G4int Z, const G4Material*);
0097 
0098   inline G4double GetZandACrossSection(const G4DynamicParticle*, 
0099                        G4int Z, G4int A);
0100 
0101   inline G4double GetCoulombBarier(const G4DynamicParticle*, 
0102                        G4double Z, G4double A, 
0103                                    G4double pR, G4double tR);
0104 
0105   G4double ComputeCoulombBarier(const G4ParticleDefinition* aParticle,
0106                 G4double kinEnergy, G4int Z, G4int A,
0107                 G4double pR, G4double tR);
0108 
0109   G4double GetRatioSD(const G4DynamicParticle*, G4double At, G4double Zt);
0110   G4double GetRatioQE(const G4DynamicParticle*, G4double At, G4double Zt);
0111 
0112   inline G4double GetElasticGlauberGribov(const G4DynamicParticle*,G4int Z, G4int A);
0113   inline G4double GetInelasticGlauberGribov(const G4DynamicParticle*,G4int Z, G4int A);
0114 
0115   inline G4double GetTotalGlauberGribovXsc() const       { return fTotalXsc;     }; 
0116   inline G4double GetElasticGlauberGribovXsc() const     { return fElasticXsc;   }; 
0117   inline G4double GetInelasticGlauberGribovXsc() const   { return fInelasticXsc; }; 
0118   inline G4double GetProductionGlauberGribovXsc() const  { return fProductionXsc; }; 
0119   inline G4double GetDiffractionGlauberGribovXsc() const { return fDiffractionXsc; }; 
0120 
0121 private:
0122 
0123   // Glauber-Gribov cross section
0124   void ComputeCrossSections(const G4ParticleDefinition* aParticle,
0125                 G4double kinEnergy, G4int Z, G4int A);
0126 
0127   G4double fTotalXsc, fElasticXsc, fInelasticXsc;
0128   G4double fProductionXsc, fDiffractionXsc;
0129   // Cache
0130   G4double fEnergy;
0131  
0132   const G4ParticleDefinition* theProton;
0133   const G4ParticleDefinition* theNeutron;
0134   const G4ParticleDefinition* theLambda;
0135 
0136   G4ComponentGGHadronNucleusXsc* fHadrNucl; 
0137   G4HadronNucleonXsc* fHNXsc;
0138 
0139   // Cache
0140   const G4ParticleDefinition* fParticle;
0141   G4int fZ, fA;    
0142 };
0143 
0144 inline G4double
0145 G4ComponentGGNuclNuclXsc::GetElasticGlauberGribov(const G4DynamicParticle* dp,
0146                                                   G4int Z, G4int A)
0147 {
0148   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0149   return fElasticXsc;
0150 }
0151 
0152 inline G4double
0153 G4ComponentGGNuclNuclXsc::GetInelasticGlauberGribov(const G4DynamicParticle* dp,
0154                                                     G4int Z, G4int A)
0155 {
0156   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0157   return fInelasticXsc;
0158 }
0159 
0160 /*
0161 inline G4double
0162 G4ComponentGGNuclNuclXsc::GetElementCrossSection(const G4DynamicParticle* dp,
0163                          G4int Z, const G4Material*)
0164 {
0165   G4int A = G4lrint(fNist->GetAtomicMassAmu(Z));
0166   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0167   return fInelasticXsc;
0168 }
0169 */
0170 inline G4double
0171 G4ComponentGGNuclNuclXsc::GetZandACrossSection(const G4DynamicParticle* dp,
0172                            G4int Z, G4int A)
0173 {
0174   ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z, A);
0175   return fInelasticXsc;
0176 }
0177 
0178 inline G4double 
0179 G4ComponentGGNuclNuclXsc::GetCoulombBarier(const G4DynamicParticle* dp, 
0180                        G4double Z, G4double A, 
0181                        G4double pR, G4double tR)
0182 {
0183   return ComputeCoulombBarier(dp->GetDefinition(), dp->GetKineticEnergy(),
0184                               G4lrint(Z), G4lrint(A), pR, tR);
0185 }
0186 
0187 #endif