Back to home page

EIC code displayed by LXR

 
 

    


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

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 #ifndef G4DNAChampionElasticModel_h
0029 #define G4DNAChampionElasticModel_h 1
0030 
0031 #include <map>
0032 #include "G4DNACrossSectionDataSet.hh"
0033 #include "G4VEmModel.hh"
0034 #include "G4Electron.hh"
0035 #include "G4ParticleChangeForGamma.hh"
0036 #include "G4LogLogInterpolation.hh"
0037 #include "G4ProductionCutsTable.hh"
0038 #include "G4NistManager.hh"
0039 
0040 class G4DNAChampionElasticModel : public G4VEmModel
0041 {
0042 
0043 public:
0044 
0045   G4DNAChampionElasticModel(const G4ParticleDefinition* p = nullptr,
0046                             const G4String& nam = "DNAChampionElasticModel");
0047 
0048   ~G4DNAChampionElasticModel() override;
0049 
0050   G4DNAChampionElasticModel & operator=(const G4DNAChampionElasticModel &right) = delete;
0051   G4DNAChampionElasticModel(const G4DNAChampionElasticModel&) = delete;
0052 
0053   void Initialise(const G4ParticleDefinition*,
0054                           const G4DataVector&) override;
0055 
0056   G4double CrossSectionPerVolume(const G4Material* material,
0057                                          const G4ParticleDefinition* p,
0058                                          G4double ekin,
0059                                          G4double emin,
0060                                          G4double emax) override;
0061 
0062   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0063                                  const G4MaterialCutsCouple*,
0064                                  const G4DynamicParticle*,
0065                                  G4double tmin,
0066                                  G4double maxEnergy) override;
0067 
0068   void SetKillBelowThreshold(G4double threshold);
0069 
0070   inline G4double GetKillBelowThreshold()
0071   {
0072     G4ExceptionDescription errMsg;
0073     errMsg << "The method G4DNAChampionElasticModel::"
0074               "GetKillBelowThreshold is deprecated";
0075     
0076     G4Exception("G4DNAChampionElasticModel::GetKillBelowThreshold",
0077                 "deprecated",
0078                 JustWarning,
0079                 errMsg);
0080     return 0.;
0081   }
0082 
0083 private:
0084   // Cross section
0085   using VecMap = std::map<G4double, std::vector<G4double>>;
0086   VecMap eVecm;
0087   using TriDimensionMap = std::map<G4double, std::map<G4double, G4double>>;
0088   TriDimensionMap eDiffCrossSectionData;
0089   std::vector<G4double> eTdummyVec;
0090 
0091   // Water density table
0092   const std::vector<G4double>* fpMolWaterDensity;
0093 
0094   // Cross section
0095   G4DNACrossSectionDataSet* fpData;
0096 
0097 protected:
0098   G4ParticleChangeForGamma* fParticleChangeForGamma;
0099 
0100 private:
0101   G4int verboseLevel;
0102   G4bool isInitialised{false};
0103   
0104   // Final state
0105 
0106   //G4double DifferentialCrossSection(G4ParticleDefinition* aParticle,
0107   //                                  G4double k, G4double theta);
0108 
0109   G4double Theta(//G4ParticleDefinition * aParticleDefinition,
0110                  G4double k,
0111                  G4double integrDiff);
0112 
0113   G4double LinLinInterpolate(G4double e1,
0114                              G4double e2,
0115                              G4double e,
0116                              G4double xs1,
0117                              G4double xs2);
0118 
0119   G4double LinLogInterpolate(G4double e1,
0120                              G4double e2,
0121                              G4double e,
0122                              G4double xs1,
0123                              G4double xs2);
0124 
0125   G4double LogLogInterpolate(G4double e1,
0126                              G4double e2,
0127                              G4double e,
0128                              G4double xs1,
0129                              G4double xs2);
0130 
0131   G4double QuadInterpolator(G4double e11,
0132                             G4double e12,
0133                             G4double e21,
0134                             G4double e22,
0135                             G4double x11,
0136                             G4double x12,
0137                             G4double x21,
0138                             G4double x22,
0139                             G4double t1,
0140                             G4double t2,
0141                             G4double t,
0142                             G4double e);
0143 
0144   G4double RandomizeCosTheta(G4double k);
0145 
0146 };
0147 
0148 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0149 
0150 #endif