Geant4  10.03
G4eplusPolarizedAnnihilation.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: G4eplusPolarizedAnnihilation.hh 96114 2016-03-16 18:51:33Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class header file
31 //
32 //
33 // File name: G4eplusPolarizedAnnihilation
34 //
35 // Author: A. Schaelicke on base of Vladimir Ivanchenko / Michel Maire code
36 //
37 // Creation date: 02.07.2006
38 //
39 // Modifications:
40 // 26-07-06 modified cross section (P. Starovoitov)
41 // 21-08-06 interface updated (A. Schaelicke)
42 // 11-06-07, add PostStepGetPhysicalInteractionLength (A.Schalicke)
43 //
44 //
45 // Class Description:
46 //
47 // Polarized process of e+ annihilation into 2 gammas
48 //
49 
50 // -------------------------------------------------------------------
51 //
52 
53 #ifndef G4eplusPolarizedAnnihilation_h
54 #define G4eplusPolarizedAnnihilation_h 1
55 
56 #include "G4eplusAnnihilation.hh"
57 #include "G4Positron.hh"
58 #include "G4VEmModel.hh"
59 
60 
62 
64 {
65 
66 public:
67 
68  explicit G4eplusPolarizedAnnihilation(const G4String& name = "pol-annihil");
69 
71 
72  // Print out of the class parameters
73  virtual void PrintInfo() override;
74 
75  virtual G4double GetMeanFreePath(const G4Track& track,
76  G4double previousStepSize,
77  G4ForceCondition* condition) override;
78 
80  const G4Track& track,
81  G4double previousStepSize,
82  G4ForceCondition* condition) override;
83 
84  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
85 
86 private:
87 
88  void CleanTables();
89 
91 
93  const G4MaterialCutsCouple* couple,
95  G4double cut,
96  G4double &tasm);
97 
99 
102 
104 
105  // for polarization:
108 
109  G4PhysicsTable* theAsymmetryTable; // table for cross section assym.
110  G4PhysicsTable* theTransverseAsymmetryTable; // table for transverse cross section assym.
111 };
112 
113 #endif
G4double condition(const G4ErrorSymMatrix &m)
G4double ComputeSaturationFactor(const G4Track &aTrack)
CLHEP::Hep3Vector G4ThreeVector
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const char * name(G4int ptype)
G4eplusPolarizedAnnihilation & operator=(const G4eplusPolarizedAnnihilation &right)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
bool G4bool
Definition: G4Types.hh:79
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
const G4ParticleDefinition * particle
G4double energy(const ThreeVector &p, const G4double m)
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4PolarizedAnnihilationModel * emModel
double G4double
Definition: G4Types.hh:76
G4ForceCondition
void BuildAsymmetryTables(const G4ParticleDefinition &part)
G4eplusPolarizedAnnihilation(const G4String &name="pol-annihil")