Geant4  10.03
G4PolarizedCompton.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 // GEANT4 Class header file
29 //
30 //
31 // File name: G4PolarizedCompton
32 //
33 // Author: Andreas Schaelicke
34 // based on code by Michel Maire / Vladimir IVANTCHENKO
35 // Class description
36 //
37 // modified version respecting media and beam polarization
38 // using the stokes formalism
39 //
40 // Creation date: 01.05.2005
41 //
42 // Modifications:
43 //
44 // 01-01-05, include polarization description (A.Stahl)
45 // 01-01-05, create asymmetry table and determine interactionlength (A.Stahl)
46 // 01-05-05, update handling of media polarization (A.Schalicke)
47 // 01-05-05, update polarized differential cross section (A.Schalicke)
48 // 26-07-06, cross section recalculated (P.Starovoitov)
49 // 09-08-06, make it work under current geant4 release (A.Schalicke)
50 // 11-06-07, add PostStepGetPhysicalInteractionLength (A.Schalicke)
51 //
52 // -----------------------------------------------------------------------------
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56 
57 #ifndef PolarizedComptonScattering_h
58 #define PolarizedComptonScattering_h 1
59 
60 #include "globals.hh"
61 #include "G4VEmProcess.hh"
62 #include "G4Gamma.hh"
63 
64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65 
68 class G4DynamicParticle;
69 
71 
73 
74 {
75 public: // with description
76 
77  explicit G4PolarizedCompton(const G4String& processName ="pol-compt",
79 
80  virtual ~G4PolarizedCompton();
81 
82  // true for Gamma only.
83  virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
84 
85  // Print few lines of informations about the process: validity range,
86  virtual void PrintInfo() override;
87 
88  void SetModel(const G4String& name);
89 
90 protected:
91 
92  virtual void InitialiseProcess(const G4ParticleDefinition*) override;
93 
94  // added for polarization treatment of polarized media:
95  virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
96 
97  virtual G4double GetMeanFreePath(const G4Track& aTrack,
98  G4double previousStepSize,
99  G4ForceCondition* condition) override;
100 
102  const G4Track& track,
103  G4double previousStepSize,
104  G4ForceCondition* condition) override;
105 
106 private:
107 
108  void CleanTable();
109 
110  void BuildAsymmetryTable(const G4ParticleDefinition& part);
111 
113  const G4MaterialCutsCouple* couple,
115  G4double cut,
116  G4double & tAsymmetry);
117 
119 
121  G4PolarizedCompton(const G4PolarizedCompton& ) = delete;
122 
125 
128 
129  // added for polarization treatment:
131  static G4PhysicsTable* theAsymmetryTable; // table for crosssection assymmetry
133 };
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
137 #endif
138 
G4double condition(const G4ErrorSymMatrix &m)
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
CLHEP::Hep3Vector G4ThreeVector
void BuildAsymmetryTable(const G4ParticleDefinition &part)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
virtual void PrintInfo() override
G4PolarizedComptonModel * emModel
const char * name(G4int ptype)
int G4int
Definition: G4Types.hh:78
G4PolarizedCompton & operator=(const G4PolarizedCompton &right)=delete
virtual void InitialiseProcess(const G4ParticleDefinition *) override
G4double ComputeSaturationFactor(const G4Track &aTrack)
bool G4bool
Definition: G4Types.hh:79
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4ParticleDefinition * particle
G4PolarizedCompton(const G4String &processName="pol-compt", G4ProcessType type=fElectromagnetic)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
G4double energy(const ThreeVector &p, const G4double m)
G4ThreeVector targetPolarization
double G4double
Definition: G4Types.hh:76
static G4PhysicsTable * theAsymmetryTable
G4ForceCondition
void SetModel(const G4String &name)
G4ProcessType