Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
processes
electromagnetic
lowenergy
include
G4LivermorePolarizedRayleighModel.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: Sebastien Incerti
29
// 30 October 2008
30
// on base of G4LowEnergyPolarizedRayleigh developed by R. Capra
31
//
32
33
#ifndef G4LivermorePolarizedRayleighModel_h
34
#define G4LivermorePolarizedRayleighModel_h 1
35
36
#include "
G4VEmModel.hh
"
37
#include "
G4ParticleChangeForGamma.hh
"
38
#include "
G4CrossSectionHandler.hh
"
39
#include "
G4LogLogInterpolation.hh
"
40
#include "
G4CompositeEMDataSet.hh
"
41
#include "
G4Gamma.hh
"
42
43
class
G4LivermorePolarizedRayleighModel
:
public
G4VEmModel
44
{
45
46
public
:
47
48
G4LivermorePolarizedRayleighModel
(
const
G4ParticleDefinition
*
p
= 0,
49
const
G4String
& nam =
"LivermorePolarizedRayleigh"
);
50
51
virtual
~G4LivermorePolarizedRayleighModel
();
52
53
virtual
void
Initialise
(
const
G4ParticleDefinition
*,
const
G4DataVector
&);
54
55
virtual
G4double
ComputeCrossSectionPerAtom
(
56
const
G4ParticleDefinition
*,
57
G4double
kinEnergy,
58
G4double
Z
,
59
G4double
A=0,
60
G4double
cut=0,
61
G4double
emax=
DBL_MAX
);
62
63
virtual
void
SampleSecondaries
(std::vector<G4DynamicParticle*>*,
64
const
G4MaterialCutsCouple
*,
65
const
G4DynamicParticle
*,
66
G4double
tmin,
67
G4double
maxEnergy);
68
69
protected
:
70
71
G4ParticleChangeForGamma
*
fParticleChange
;
72
73
private
:
74
75
G4double
lowEnergyLimit;
76
G4double
highEnergyLimit;
77
78
G4int
verboseLevel;
79
G4bool
isInitialised;
80
81
G4VCrossSectionHandler
* crossSectionHandler;
82
G4VEMDataSet
* formFactorData;
83
84
// Generates \f$cos \left ( \theta\right )\f$ of the scattered photon
85
// incomingPhotonEnergy The energy of the incoming photon
86
// zAtom Atomic number
87
// \f$cos \left ( \theta\right )\f$
88
G4double
GenerateCosTheta(
G4double
incomingPhotonEnergy,
G4int
zAtom)
const
;
89
90
// Generates \f$\phi\f$ of the scattered photon
91
// cosTheta \f$cos \left ( \theta\right )\f$ of the scattered photon
92
// \f$\phi\f$
93
G4double
GeneratePhi(
G4double
cosTheta)
const
;
94
95
// Generates the polarization direction \f$\beta\f$ in the plane x, y relative to the x direction
96
// \f$\beta\f$
97
G4double
GeneratePolarizationAngle(
void
)
const
;
98
99
G4ThreeVector
GetPhotonPolarization(
const
G4DynamicParticle
&
photon
);
100
101
G4LivermorePolarizedRayleighModel
& operator=(
const
G4LivermorePolarizedRayleighModel
&
right
);
102
G4LivermorePolarizedRayleighModel
(
const
G4LivermorePolarizedRayleighModel
&);
103
104
};
105
106
//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
108
#endif
Generated on Sat May 25 2013 14:33:32 for Geant4 by
1.8.4