source: trunk/source/processes/electromagnetic/polarisation/include/G4PolarizedAnnihilationCrossSection.hh @ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 15 years ago

update CVS release candidate geant4.9.3.01

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