Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeRayleighModelMI.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4PenelopeRayleighModelMI.hh 75573 2013-11-04 11:48:15Z gcosmo $
27//
28// Author: Luciano Pandola and Gianfranco Paternò
29//
30// -------------------------------------------------------------------
31// History:
32// 03 Dec 2009 L. Pandola 1st implementation
33// 25 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
34// 27 Sep 2013 L. Pandola Migration to MT paradigm
35// 20 Aug 2017 G. Paternò Molecular Interference implementation
36// 24 Mar 2019 G. Paternò Improved Molecular Interference implementation
37// 20 Jun 2020 G. Paternò Read qext separately and leave original atomic
38// form factors
39// 27 Aug 2020 G. Paternò Further improvement of MI implementation
40// 04 Mar 2021 L. Pandola Replace maps with arrays
41//
42// -------------------------------------------------------------------
43// Class description:
44// Low Energy Electromagnetic Physics, Rayleigh Scattering
45// with the model from Penelope, version 2008
46// extended for Molecular Interference Effects
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
51#ifndef G4PenelopeRayleighModelMI_HH
52#define G4PenelopeRayleighModelMI_HH 1
53
54#include "globals.hh"
55#include "G4VEmModel.hh"
56#include "G4DataVector.hh"
58
59#include "G4ExtendedMaterial.hh"
60#include "G4MIData.hh"
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
67class G4Material;
70
72{
73
74public:
75 explicit G4PenelopeRayleighModelMI(const G4ParticleDefinition* p = nullptr,
76 const G4String& processName = "PenRayleighMI");
77
79
80 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
82 G4VEmModel *masterModel) override;
83
85 G4double kinEnergy,
86 G4double Z,
87 G4double A = 0,
88 G4double cut = 0,
89 G4double emax = DBL_MAX) override;
90
91 //Overriding of parent's (G4VEmModel) method
94 G4double kineticEnergy,
95 G4double cutEnergy = 0.,
96 G4double maxEnergy = DBL_MAX) override;
97
98 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
100 const G4DynamicParticle*,
101 G4double tmin,
102 G4double maxEnergy) override;
103
104 void SetVerbosityLevel(G4int lev) {fVerboseLevel = lev;};
105 G4int GetVerbosityLevel() {return fVerboseLevel;};
106
107 //Testing purposes
108 void DumpFormFactorTable(const G4Material*);
109
110 //Settings
111 void SetMIActive(G4bool val){fIsMIActive = val;};
112 G4bool IsMIActive(){return fIsMIActive;};
113
116
117private:
118 void SetParticle(const G4ParticleDefinition*);
119
120 //Helper methods
121 void ReadDataFile(G4int);
122 void ClearTables();
123 void BuildFormFactorTable(const G4Material*);
124 void GetPMaxTable(const G4Material*);
125 G4double GetFSquared(const G4Material*,const G4double);
126 void InitializeSamplingAlgorithm(const G4Material*);
127 void ReadMolInterferenceData(const G4String&,const G4String& filename="NULL");
128 G4MIData* GetMIData(const G4Material*);
129 void CalculateThetaAndAngFun();
130 G4double CalculateQSquared(G4double angle, G4double energy);
131 G4double IntegrateFun(G4double y[], G4int n, G4double dTheta);
132 void LoadKnownMIFFMaterials();
133
134 /// Data members
135 G4ParticleChangeForGamma* fParticleChange;
136 const G4ParticleDefinition* fParticle;
137
138 G4DataVector fLogQSquareGrid; //log(Q^2) grid for interpolation
139 std::map<const G4Material*,G4PhysicsFreeVector*> *fLogFormFactorTable; //log(Q^2) vs. log(F^2)
140
141 G4DataVector fLogEnergyGridPMax; //energy grid for PMax (and originally for the x-section)
142 std::map<const G4Material*,G4PhysicsFreeVector*> *fPMaxTable; //E vs. Pmax
143 std::map<const G4Material*,G4PenelopeSamplingData*> *fSamplingTable;
144 //Internal tables and manager methods
145 std::map<G4String,G4PhysicsFreeVector*> *fMolInterferenceData;
146 G4PhysicsFreeVector* fAngularFunction;
147 std::map<G4String,G4String> *fKnownMaterials;
148 static const G4int fMaxZ =99;
149 static G4PhysicsFreeVector* fLogAtomicCrossSection[fMaxZ+1];
150 static G4PhysicsFreeVector* fAtomicFormFactor[fMaxZ+1];
151
152 //Intrinsic energy limits of the model: cannot be extended by the parent process
153 G4double fIntrinsicLowEnergyLimit;
154 G4double fIntrinsicHighEnergyLimit;
155 G4double fDTheta = {0.0001};
156
157 static const G4int fNtheta = 31415;
158 G4int fVerboseLevel;
159 G4bool fIsInitialised;
160 //Used only for G4EmCalculator and Unit Tests
161 G4bool fLocalTable;
162 G4bool fIsMIActive;
163};
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167#endif
168
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4PenelopeRayleighModelMI(const G4PenelopeRayleighModelMI &)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
G4PenelopeRayleighModelMI & operator=(const G4PenelopeRayleighModelMI &right)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void DumpFormFactorTable(const G4Material *)
#define DBL_MAX
Definition: templates.hh:62