Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedIonisationModel.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: G4PolarizedIonisationModel
31//
32// Author: A.Schaelicke on base of Vladimir Ivanchenko code
33//
34// Class Description:
35// Physics implementation of polarized Bhabha/Moller scattering
36//
37// -------------------------------------------------------------------
38
39#ifndef G4PolarizedIonisationModel_h
40#define G4PolarizedIonisationModel_h 1
41
42#include "globals.hh"
44#include "G4StokesVector.hh"
45#include "G4ThreeVector.hh"
46
50class G4VPolarizedXS;
51
53{
54 public:
56 const G4ParticleDefinition* p = nullptr,
57 const G4String& nam = "PolarizedMollerBhabha");
58
59 virtual ~G4PolarizedIonisationModel() override;
60
62 G4double kinEnergy,
63 G4double cut,
64 G4double emax) override;
65
66 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
68 const G4DynamicParticle*, G4double tmin,
69 G4double maxEnergy) override;
70
73 const G4PolarizedIonisationModel& right) = delete;
74
76 {
77 fTargetPolarization = G4StokesVector(pTarget);
78 }
80 {
81 fBeamPolarization = G4StokesVector(pBeam);
82 }
83 const G4StokesVector& GetTargetPolarization() { return fTargetPolarization; }
84 const G4StokesVector& GetBeamPolarization() { return fBeamPolarization; }
86 {
87 return fElectronPolarization;
88 }
90 {
91 return fPositronPolarization;
92 }
93
94 private:
95 G4VPolarizedXS* fCrossSectionCalculator;
96
97 G4StokesVector fBeamPolarization;
98 G4StokesVector fTargetPolarization;
99
100 G4StokesVector fPositronPolarization;
101 G4StokesVector fElectronPolarization;
102};
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105
106#endif
double G4double
Definition: G4Types.hh:83
virtual ~G4PolarizedIonisationModel() override
G4PolarizedIonisationModel(G4PolarizedIonisationModel &)=delete
const G4StokesVector & GetFinalPositronPolarization()
const G4StokesVector & GetTargetPolarization()
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SetTargetPolarization(const G4ThreeVector &pTarget)
G4PolarizedIonisationModel & operator=(const G4PolarizedIonisationModel &right)=delete
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax) override
const G4StokesVector & GetBeamPolarization()
const G4StokesVector & GetFinalElectronPolarization()
void SetBeamPolarization(const G4ThreeVector &pBeam)