Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedAnnihilationCrossSection.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// -------------------------------------------------------------------
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: PolarizedAnnihilationCrossSection
33//
34// Author: Andreas Schaelicke and Pavel Starovoitov
35//
36// Creation date: 22.03.2006
37//
38// Modifications:
39// 15.10.07 introduced a more general framework for cross sections
40//
41// Class Description:
42// * calculates the differential cross section (ME squared,
43// without phase space) incoming positron (along positive z direction)
44// annihilations with an electron at rest
45// * phi denotes the angle between the scattering plane and
46// X axis of incoming partice reference frame (PRF)
47//
48#ifndef G4PolarizedAnnihilationCrossSection_h
49#define G4PolarizedAnnihilationCrossSection_h 1
50
51#include "G4StokesVector.hh"
53
55{
56public:
59public:
60 virtual void Initialize(G4double eps, G4double gamma, G4double phi,
61 const G4StokesVector & p0,const G4StokesVector & p1,
62 G4int flag=0) override;
63
65 virtual G4double XSection(const G4StokesVector & pol2,
66 const G4StokesVector & pol3) override;
67 virtual G4double TotalXSection(G4double xmin, G4double xmax,
68 G4double y,
69 const G4StokesVector & pol0,
70 const G4StokesVector & pol1) override;
71
72 // return expected mean polarisation
73 G4StokesVector GetPol2() override;
74 G4StokesVector GetPol3() override;
75
76 virtual G4double GetXmin(G4double y) override; // minimal energy fraction in TotalXSection
77 virtual G4double GetXmax(G4double y) override; // maximal energy fraction in TotalXSection
78
80 // test routine
81 void getCoeff();
82private:
83 void TotalXS();
84 void DefineCoefficients(const G4StokesVector & pol0,
85 const G4StokesVector & pol1);
86
87
88 G4double polxx, polyy, polzz, polxz, polzx, polxy, polyx, polyz, polzy;
89
90 G4double re2, diffXSFactor, totalXSFactor;
91 // - unpolarised + part depending on the polarization of the initial pair
92 G4double phi0;
93 // - part depending on the polarization of the final positron
94 G4ThreeVector phi2;
95 // - part depending on the polarization of the final electron
96 G4ThreeVector phi3;
97 G4double dice;
98 G4double polXS, unpXS;
99 G4double ISPxx, ISPyy, ISPzz, ISPnd;
100};
101#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1) override
virtual void Initialize(G4double eps, G4double gamma, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
virtual G4double XSection(const G4StokesVector &pol2, const G4StokesVector &pol3) override