Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/Geant4/G4PenelopeRayleighModelMI.hh was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

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 // $Id: G4PenelopeRayleighModelMI.hh 75573 2013-11-04 11:48:15Z gcosmo $
0027 //
0028 // Author: Luciano Pandola and Gianfranco Paternò
0029 //
0030 // -------------------------------------------------------------------
0031 // History:
0032 // 03 Dec 2009   L. Pandola   1st implementation 
0033 // 25 May 2011   L. Pandola   Renamed (make v2008 as default Penelope)
0034 // 27 Sep 2013   L. Pandola   Migration to MT paradigm
0035 // 20 Aug 2017   G. Paternò   Molecular Interference implementation
0036 // 24 Mar 2019   G. Paternò   Improved Molecular Interference implementation
0037 // 20 Jun 2020   G. Paternò   Read qext separately and leave original atomic 
0038 //                            form factors
0039 // 27 Aug 2020   G. Paternò   Further improvement of MI implementation
0040 // 04 Mar 2021   L. Pandola   Replace maps with arrays
0041 //
0042 // -------------------------------------------------------------------
0043 // Class description:
0044 // Low Energy Electromagnetic Physics, Rayleigh Scattering
0045 // with the model from Penelope, version 2008
0046 // extended for Molecular Interference Effects
0047 // -------------------------------------------------------------------
0048 //
0049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0050 
0051 #ifndef G4PenelopeRayleighModelMI_HH
0052 #define G4PenelopeRayleighModelMI_HH 1
0053 
0054 #include "globals.hh"
0055 #include "G4VEmModel.hh"
0056 #include "G4DataVector.hh"
0057 #include "G4ParticleChangeForGamma.hh"
0058 
0059 #include "G4ExtendedMaterial.hh"
0060 #include "G4MIData.hh"
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0063 
0064 class G4ParticleDefinition;
0065 class G4DynamicParticle;
0066 class G4MaterialCutsCouple;
0067 class G4Material;
0068 class G4PhysicsFreeVector;
0069 class G4PenelopeSamplingData;
0070 
0071 class G4PenelopeRayleighModelMI : public G4VEmModel 
0072 {
0073 
0074 public:  
0075   explicit G4PenelopeRayleighModelMI(const G4ParticleDefinition* p = nullptr,
0076                      const G4String& processName = "PenRayleighMI");
0077 
0078   virtual ~G4PenelopeRayleighModelMI();
0079   
0080   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0081   void InitialiseLocal(const G4ParticleDefinition*,
0082                    G4VEmModel *masterModel) override;
0083 
0084   G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
0085                       G4double kinEnergy,
0086                       G4double Z,
0087                       G4double A = 0,
0088                       G4double cut = 0,
0089                       G4double emax = DBL_MAX) override;
0090   
0091   //Overriding of parent's (G4VEmModel) method                                           
0092   G4double CrossSectionPerVolume(const G4Material*,
0093                  const G4ParticleDefinition*,
0094                  G4double kineticEnergy,
0095                  G4double cutEnergy = 0.,
0096                  G4double maxEnergy = DBL_MAX) override;    
0097   
0098   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0099              const G4MaterialCutsCouple*,
0100              const G4DynamicParticle*,
0101              G4double tmin,
0102              G4double maxEnergy) override;
0103   
0104   void SetVerbosityLevel(G4int lev) {fVerboseLevel = lev;};
0105   G4int GetVerbosityLevel() {return fVerboseLevel;};
0106   
0107   //Testing purposes
0108   void DumpFormFactorTable(const G4Material*);
0109   
0110   //Settings
0111   void SetMIActive(G4bool val){fIsMIActive = val;};
0112   G4bool IsMIActive(){return fIsMIActive;};
0113 
0114   G4PenelopeRayleighModelMI& operator=(const G4PenelopeRayleighModelMI &right) = delete;
0115   G4PenelopeRayleighModelMI(const G4PenelopeRayleighModelMI&) = delete;
0116 
0117 private:
0118    void SetParticle(const G4ParticleDefinition*);
0119 
0120   //Helper methods
0121   void ReadDataFile(G4int);   
0122   void ClearTables();
0123   void BuildFormFactorTable(const G4Material*);
0124   void GetPMaxTable(const G4Material*); 
0125   G4double GetFSquared(const G4Material*,const G4double);
0126   void InitializeSamplingAlgorithm(const G4Material*);
0127   void ReadMolInterferenceData(const G4String&,const G4String& filename="NULL");
0128   G4MIData* GetMIData(const G4Material*);
0129   void CalculateThetaAndAngFun();  
0130   G4double CalculateQSquared(G4double angle, G4double energy);
0131   G4double IntegrateFun(G4double y[], G4int n, G4double dTheta);
0132   void LoadKnownMIFFMaterials();
0133 
0134   /// Data members
0135   G4ParticleChangeForGamma* fParticleChange;
0136   const G4ParticleDefinition* fParticle;
0137   
0138   G4DataVector fLogQSquareGrid; //log(Q^2) grid for interpolation
0139   std::map<const G4Material*,G4PhysicsFreeVector*> *fLogFormFactorTable; //log(Q^2) vs. log(F^2)
0140   
0141   G4DataVector fLogEnergyGridPMax; //energy grid for PMax (and originally for the x-section)
0142   std::map<const G4Material*,G4PhysicsFreeVector*> *fPMaxTable; //E vs. Pmax
0143   std::map<const G4Material*,G4PenelopeSamplingData*> *fSamplingTable;
0144   //Internal tables and manager methods
0145   std::map<G4String,G4PhysicsFreeVector*> *fMolInterferenceData; 
0146   G4PhysicsFreeVector* fAngularFunction;
0147   std::map<G4String,G4String> *fKnownMaterials;
0148   static const G4int fMaxZ =99;
0149   static G4PhysicsFreeVector* fLogAtomicCrossSection[fMaxZ+1];
0150   static G4PhysicsFreeVector* fAtomicFormFactor[fMaxZ+1];
0151 
0152   //Intrinsic energy limits of the model: cannot be extended by the parent process
0153   G4double fIntrinsicLowEnergyLimit;
0154   G4double fIntrinsicHighEnergyLimit;
0155   G4double fDTheta = {0.0001};
0156 
0157   static const G4int fNtheta = 31415;
0158   G4int fVerboseLevel;
0159   G4bool fIsInitialised;  
0160   //Used only for G4EmCalculator and Unit Tests
0161   G4bool fLocalTable;  
0162   G4bool fIsMIActive;           
0163 };
0164 
0165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0166 
0167 #endif
0168