Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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// File name: G4PolarizedCompton
31//
32// Author: Andreas Schaelicke
33// based on code by Michel Maire / Vladimir IVANTCHENKO
34//
35// Class description
36// polarized version of Compton scattering
37//
38// -----------------------------------------------------------------------------
39
40#ifndef PolarizedComptonScattering_h
41#define PolarizedComptonScattering_h 1
42
43#include "globals.hh"
44#include "G4Gamma.hh"
45#include "G4VEmProcess.hh"
46
51
53
54{
55 public:
56 explicit G4PolarizedCompton(const G4String& processName = "pol-compt",
58
59 virtual ~G4PolarizedCompton() override;
60
61 // true for Gamma only.
62 virtual G4bool IsApplicable(const G4ParticleDefinition&) override;
63
64 virtual void ProcessDescription(std::ostream&) const override;
65 virtual void DumpInfo() const override { ProcessDescription(G4cout); };
66
67 void SetModel(const G4String& name);
68
71
72 protected:
73 virtual void InitialiseProcess(const G4ParticleDefinition*) override;
74
75 virtual void BuildPhysicsTable(const G4ParticleDefinition&) override;
76
77 virtual G4double GetMeanFreePath(const G4Track& aTrack,
78 G4double previousStepSize,
79 G4ForceCondition* condition) override;
80
82 const G4Track& track, G4double previousStepSize,
83 G4ForceCondition* condition) override;
84
85 private:
86 static G4PhysicsTable* theAsymmetryTable; // table for crosssection asymmetry
87 void CleanTable();
88
89 void BuildAsymmetryTable(const G4ParticleDefinition& part);
90
91 G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple* couple,
92 const G4ParticleDefinition& particle, G4double cut,
93 G4double& tAsymmetry);
94
95 G4double ComputeSaturationFactor(const G4Track& aTrack);
96
98
99 G4int fType;
100
101 G4bool fBuildAsymmetryTable;
102 G4bool fUseAsymmetryTable;
103
104 G4bool fIsInitialised;
105};
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
109#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
void SetModel(const G4String &name)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
virtual void ProcessDescription(std::ostream &) const override
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
virtual ~G4PolarizedCompton() override
virtual void DumpInfo() const override
G4PolarizedCompton & operator=(const G4PolarizedCompton &right)=delete
virtual G4bool IsApplicable(const G4ParticleDefinition &) override
G4PolarizedCompton(const G4PolarizedCompton &)=delete
virtual void InitialiseProcess(const G4ParticleDefinition *) override