Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermorePolarizedGammaConversionModel.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// Authors: G.Depaola & F.Longo
29//
30
31#ifndef G4LivermorePolarizedGammaConversionModel_h
32#define G4LivermorePolarizedGammaConversionModel_h 1
33
34#include "G4VEmModel.hh"
35#include "G4Electron.hh"
36#include "G4Positron.hh"
42#include "G4ForceCondition.hh"
43
44
46{
47
48public:
49
51 const G4String& nam = "LivermorePolarizedGammaConversion");
52
54
55 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
56
59 G4double kinEnergy,
60 G4double Z,
61 G4double A=0,
62 G4double cut=0,
63 G4double emax=DBL_MAX);
64
65 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
67 const G4DynamicParticle*,
68 G4double tmin,
69 G4double maxEnergy);
70
71protected:
72
74
76 G4double previousStepSize,
78private:
79
80 G4double lowEnergyLimit;
81 G4double highEnergyLimit;
82 G4bool isInitialised;
83 G4int verboseLevel;
84
85 G4VEMDataSet* meanFreePathTable;
86 G4VCrossSectionHandler* crossSectionHandler;
87
88
89 // specific methods for polarization
90
91 G4ThreeVector GetRandomPolarization(G4ThreeVector& direction0); // Random Polarization
92 G4ThreeVector GetPerpendicularPolarization(const G4ThreeVector& direction0, const G4ThreeVector& polarization0) const;
93 G4ThreeVector SetPerpendicularVector(G4ThreeVector& a); // temporary
94 void SystemOfRefChange(G4ThreeVector& direction0, G4ThreeVector& direction1,
95 G4ThreeVector& polarization0);
96
97 // Gamma Conversion methods
98
99
100 G4double SetPhi(G4double);
101 G4double SetPsi(G4double, G4double); //
102
103 G4double ScreenFunction1(G4double screenVariable);
104 G4double ScreenFunction2(G4double screenVariable);
105 void SetTheta(G4double*, G4double*, G4double);
106
109
111
112 G4double Flor(G4double*, G4double);
113 G4double Glor(G4double*, G4double);
114
115 G4double Fdlor(G4double*, G4double);
116 G4double Fintlor(G4double*, G4double);
118
119 G4double Ftan(G4double*, G4double);
120
121 //G4double Gtan(G4double*, G4double);
122
123 G4double Fdtan(G4double*, G4double);
124 G4double Finttan(G4double*, G4double);
126
127 G4double smallEnergy;
128 G4double Psi, Phi;
129
130
133
134};
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
#define DBL_MAX
Definition: templates.hh:83