Geant4 10.7.0
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 form factors
38// 27 Aug 2020 G. Paternò Further improvement of MI implementation
39//
40// -------------------------------------------------------------------
41// Class description:
42// Low Energy Electromagnetic Physics, Rayleigh Scattering
43// with the model from Penelope, version 2008
44// extended for Molecular Interference Effects
45// -------------------------------------------------------------------
46//
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49#ifndef G4PenelopeRayleighModelMI_HH
50#define G4PenelopeRayleighModelMI_HH 1
51
52#include "globals.hh"
53#include "G4VEmModel.hh"
54#include "G4DataVector.hh"
56
57#include "G4ExtendedMaterial.hh"
58#include "G4MIData.hh"
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
65class G4Material;
68
70{
71
72public:
74 const G4String& processName = "PenRayleighMI");
75
77
78 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
79 virtual void InitialiseLocal(const G4ParticleDefinition*,
80 G4VEmModel *masterModel) override;
81
83 G4double kinEnergy,
84 G4double Z,
85 G4double A = 0,
86 G4double cut = 0,
87 G4double emax = DBL_MAX) override;
88
89 //Overriding of parent's (G4VEmModel) method
92 G4double kineticEnergy,
93 G4double cutEnergy = 0.,
94 G4double maxEnergy = DBL_MAX) override;
95
96 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
98 const G4DynamicParticle*,
99 G4double tmin,
100 G4double maxEnergy) override;
101
102 void SetVerbosityLevel(G4int lev) {verboseLevel = lev;};
103 G4int GetVerbosityLevel() {return verboseLevel;};
104
105 //Testing purposes
106 void DumpFormFactorTable(const G4Material*);
107
108 //Settings
109 void SetMIActive(G4bool val){fIsMIActive = val;};
110 G4bool IsMIActive(){return fIsMIActive;};
111
112private:
115 void SetParticle(const G4ParticleDefinition*);
116
117 //Helper methods
118 void ReadDataFile(G4int);
119 void ClearTables();
120 void BuildFormFactorTable(const G4Material*);
121 void GetPMaxTable(const G4Material*);
122 G4double GetFSquared(const G4Material*,const G4double);
123 void InitializeSamplingAlgorithm(const G4Material*);
124 void ReadMolInterferenceData(const G4String&,const G4String& filename="NULL");
125 G4MIData* GetMIData(const G4Material*);
126 void CalculateThetaAndAngFun();
127 G4double CalculateQSquared(G4double angle, G4double energy);
128 G4double IntegrateFun(G4double y[], G4int n, G4double dTheta);
129 void LoadKnownMIFFMaterials();
130
131
132 /// Data members
133 G4ParticleChangeForGamma* fParticleChange;
134 const G4ParticleDefinition* fParticle;
135
136 //Intrinsic energy limits of the model: cannot be extended by the parent process
137 G4double fIntrinsicLowEnergyLimit;
138 G4double fIntrinsicHighEnergyLimit;
139
140 G4int verboseLevel;
141 G4bool isInitialised;
142
143 //Internal tables and manager methods
144 std::map<G4int,G4PhysicsFreeVector*> *logAtomicCrossSection;
145 std::map<G4int,G4PhysicsFreeVector*> *atomicFormFactor;
146 std::map<G4String,G4PhysicsFreeVector*> *MolInterferenceData; //G. Paternò
147
148 G4DataVector logQSquareGrid; //log(Q^2) grid for interpolation
149 std::map<const G4Material*,G4PhysicsFreeVector*> *logFormFactorTable; //log(Q^2) vs. log(F^2)
150
151 G4DataVector logEnergyGridPMax; //energy grid for PMax (and originally for the x-section)
152 std::map<const G4Material*,G4PhysicsFreeVector*> *pMaxTable; //E vs. Pmax
153 std::map<const G4Material*,G4PenelopeSamplingData*> *samplingTable;
154
155 //Used only for G4EmCalculator and Unit Tests
156 G4bool fLocalTable;
157 static const G4int Ntheta = 31415;
158 G4double fDTheta = {0.0001};
159 G4bool fIsMIActive;
160 G4PhysicsFreeVector* angularFunction;
161 std::map<G4String,G4String> *fKnownMaterials;
162
163};
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
166
167#endif
168
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void DumpFormFactorTable(const G4Material *)
#define DBL_MAX
Definition: templates.hh:62