source: trunk/source/processes/electromagnetic/polarisation/include/G4PolarizedAnnihilationModel.hh @ 1350

Last change on this file since 1350 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 5.4 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// $Id: G4PolarizedAnnihilationModel.hh,v 1.3 2007/07/10 09:38:17 schaelic Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4PolarizedAnnihilationModel
35//
36// Author:        Andreas Schaelicke and Pavel Starovoitov
37//
38// Creation date: 01.05.2005
39//
40// Modifications:
41// 18-07-06 use newly calculated cross sections (P. Starovoitov)
42// 21-08-06 update interface to geant4.8.1 (A. Schaelicke)
43// 10-07-07 copied Initialise() method from G4eeToTwoGammaModel to provide a 
44//          ParticleChangeForGamma object (A. Schaelicke)
45//
46//
47// Class Description:
48//
49// Implementation of polarized gamma Annihilation scattering on free electron
50//
51
52// -------------------------------------------------------------------
53//
54
55#ifndef G4PolarizedAnnihilationModel_h
56#define G4PolarizedAnnihilationModel_h 1
57
58#include "G4eeToTwoGammaModel.hh"
59#include "G4StokesVector.hh"
60
61class G4ParticleChangeForGamma;
62class G4PolarizedAnnihilationCrossSection;
63
64class G4PolarizedAnnihilationModel : public G4eeToTwoGammaModel
65{
66
67public:
68
69  G4PolarizedAnnihilationModel(const G4ParticleDefinition* p = 0, 
70                        const G4String& nam = "Polarized-Annihilation");
71
72  virtual ~G4PolarizedAnnihilationModel();
73
74  virtual void Initialise(const G4ParticleDefinition*, 
75                          const G4DataVector&);
76  virtual G4double ComputeCrossSectionPerElectron(
77                                const G4ParticleDefinition*,
78                                      G4double kinEnergy, 
79                                      G4double cut,
80                                      G4double emax);
81  void ComputeAsymmetriesPerElectron(G4double gammaEnergy,
82                                     G4double & valueX,
83                                     G4double & valueA,
84                                     G4double & valueT);
85
86  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
87                                 const G4MaterialCutsCouple*,
88                                 const G4DynamicParticle*,
89                                 G4double tmin,
90                                 G4double maxEnergy);
91
92  // polarized routines
93  inline void SetTargetPolarization(const G4ThreeVector & pTarget);
94  inline void SetBeamPolarization(const G4ThreeVector & pBeam);
95  inline const G4ThreeVector & GetTargetPolarization() const;
96  inline const G4ThreeVector & GetBeamPolarization() const;
97  inline const G4ThreeVector & GetFinalGamma1Polarization() const;
98  inline const G4ThreeVector & GetFinalGamma2Polarization() const;
99private:
100
101  // hide assignment operator
102  G4PolarizedAnnihilationModel & operator=(const  G4PolarizedAnnihilationModel &right);
103  G4PolarizedAnnihilationModel(const  G4PolarizedAnnihilationModel&);
104
105  G4PolarizedAnnihilationCrossSection * crossSectionCalculator;
106  // incomming
107  G4StokesVector theBeamPolarization;    // positron
108  G4StokesVector theTargetPolarization;  // electron
109  // outgoing
110  G4StokesVector finalGamma1Polarization;
111  G4StokesVector finalGamma2Polarization;
112
113  G4int verboseLevel;
114
115  G4ParticleChangeForGamma* gParticleChange;
116  G4bool gIsInitialised;
117};
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
121inline void G4PolarizedAnnihilationModel::SetTargetPolarization(const G4ThreeVector & pTarget)
122{
123  theTargetPolarization = pTarget;
124}
125inline void G4PolarizedAnnihilationModel::SetBeamPolarization(const G4ThreeVector & pBeam)
126{
127  theBeamPolarization = pBeam;
128}
129inline const G4ThreeVector & G4PolarizedAnnihilationModel::GetTargetPolarization() const
130{
131  return theTargetPolarization;
132}
133inline const G4ThreeVector & G4PolarizedAnnihilationModel::GetBeamPolarization() const
134{
135  return theBeamPolarization;
136}
137inline const G4ThreeVector &  G4PolarizedAnnihilationModel::GetFinalGamma1Polarization() const
138{
139  return finalGamma1Polarization;
140}
141inline const G4ThreeVector & G4PolarizedAnnihilationModel::GetFinalGamma2Polarization() const
142{
143  return finalGamma2Polarization;
144}
145
146#endif
Note: See TracBrowser for help on using the repository browser.