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