Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4PenelopeRayleighModel.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$
27//
28// Author: Luciano Pandola
29//
30// History:
31// -----------
32// 03 Dec 2009 L. Pandola 1st implementation.
33// 25 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
34//
35// -------------------------------------------------------------------
36//
37// Class description:
38// Low Energy Electromagnetic Physics, Rayleigh Scattering
39// with the model from Penelope, version 2008
40// -------------------------------------------------------------------
41
42#ifndef G4PENELOPERAYLEIGHMODEL_HH
43#define G4PENELOPERAYLEIGHMODEL_HH 1
44
45#include "globals.hh"
46#include "G4VEmModel.hh"
47#include "G4DataVector.hh"
49
53class G4Material;
56
58{
59
60public:
62 const G4String& processName ="PenRayleigh");
63
65
66 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
67
69 G4double kinEnergy,
70 G4double Z,
71 G4double A=0,
72 G4double cut=0,
73 G4double emax=DBL_MAX);
74
75 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
77 const G4DynamicParticle*,
78 G4double tmin,
79 G4double maxEnergy);
80
81 void SetVerbosityLevel(G4int lev){verboseLevel = lev;};
82 G4int GetVerbosityLevel(){return verboseLevel;};
83
84 //Testing purposes
86
87protected:
89
90private:
91 G4PenelopeRayleighModel& operator=(const G4PenelopeRayleighModel &right);
93
94 //Intrinsic energy limits of the model: cannot be extended by
95 //the parent process
96 G4double fIntrinsicLowEnergyLimit;
97 G4double fIntrinsicHighEnergyLimit;
98
99 G4int verboseLevel;
100 G4bool isInitialised;
101
102 //Internal tables and manager methods
103 std::map<G4int,G4PhysicsFreeVector*> *logAtomicCrossSection;
104 std::map<G4int,G4PhysicsFreeVector*> *atomicFormFactor;
105
106
107 G4DataVector logQSquareGrid; //log(Q^2) grid for interpolation
108 std::map<const G4Material*,G4PhysicsFreeVector*> *logFormFactorTable; //log(Q^2) vs. log(F^2)
109
110 G4DataVector logEnergyGridPMax; //energy grid for PMac (and originally for the x-section)
111 std::map<const G4Material*,G4PhysicsFreeVector*> *pMaxTable; //E vs. Pmax
112
113 std::map<const G4Material*,G4PenelopeSamplingData*> *samplingTable;
114
115 //Helper methods
116 void ReadDataFile(G4int);
117 void ClearTables();
118 void BuildFormFactorTable(const G4Material*);
119 void GetPMaxTable(const G4Material*);
120
121 G4double GetFSquared(const G4Material*,const G4double);
122 void InitializeSamplingAlgorithm(const G4Material*);
123
124};
125
126#endif
127
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
void DumpFormFactorTable(const G4Material *)
G4ParticleChangeForGamma * fParticleChange
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
#define DBL_MAX
Definition: templates.hh:83