Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedPhotoElectricModel.cc
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 file
29//
30// File name: G4PolarizedPhotoElectricModel
31//
32// Author: Andreas Schaelicke & Karim Laihem
33//
34// Class Description:
35// Implementation of Photo electric effect
36// including polarization transfer from circularly polarised gammas
37
39
40#include "G4DataVector.hh"
44
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
47 const G4ParticleDefinition*, const G4String& nam)
49 , fCrossSectionCalculator(nullptr)
50 , fVerboseLevel(0)
51{}
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55{
56 if(fCrossSectionCalculator)
57 delete fCrossSectionCalculator;
58}
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 const G4DataVector& dv)
63{
65 if(!fCrossSectionCalculator)
66 fCrossSectionCalculator = new G4PolarizedPhotoElectricXS();
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71 std::vector<G4DynamicParticle*>* vdp, const G4MaterialCutsCouple* couple,
72 const G4DynamicParticle* dp, G4double tmin, G4double maxEnergy)
73{
74 G4PEEffectFluoModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
75
76 if(fVerboseLevel >= 1)
77 {
78 G4cout << "G4PolarizedPhotoElectricModel::SampleSecondaries" << G4endl;
79 }
80
81 if(vdp && !vdp->empty())
82 {
83 G4double gamEnergy0 = dp->GetKineticEnergy();
84 G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
85 G4double sintheta =
86 dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
87 if(sintheta > 1.)
88 sintheta = 1.;
89
91 beamPol.SetPhoton();
92
93 // determine interaction plane
95 dp->GetMomentumDirection(), (*vdp)[0]->GetMomentumDirection());
97 .cross((*vdp)[0]->GetMomentumDirection())
98 .mag() < 1.e-10)
99 {
100 nInteractionFrame =
102 }
103
104 // transform polarization into interaction frame
105 beamPol.InvRotateAz(nInteractionFrame, dp->GetMomentumDirection());
106
107 // calulcate polarization transfer
108 fCrossSectionCalculator->SetMaterial(
109 GetCurrentElement()->GetN(), // number of nucleons
110 GetCurrentElement()->GetZ(), GetCurrentElement()->GetfCoulomb());
111 fCrossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
112 beamPol, G4StokesVector::ZERO);
113
114 // determine final state polarization
115 G4StokesVector lep1Pol = fCrossSectionCalculator->GetPol2();
116 lep1Pol.RotateAz(nInteractionFrame, (*vdp)[0]->GetMomentumDirection());
117 (*vdp)[0]->SetPolarization(lep1Pol.p1(), lep1Pol.p2(), lep1Pol.p3());
118
119 if(vdp->size() != 1)
120 {
122 ed << " WARNING " << vdp->size()
123 << " secondaries in polarized photo electric effect not supported!\n";
124 G4Exception("G4PolarizedPhotoElectricModel::SampleSecondaries", "pol024",
125 JustWarning, ed);
126 }
127 }
128}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
double mag() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
static G4ThreeVector GetRandomFrame(const G4ThreeVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *pd, const G4DataVector &dv) override
G4PolarizedPhotoElectricModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Polarized-PhotoElectric")
void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0) override
G4double p3() const
G4double p1() const
static const G4StokesVector ZERO
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4double p2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4Element * GetCurrentElement(const G4Material *mat=nullptr) const
Definition: G4VEmModel.cc:242
void SetMaterial(G4double A, G4double Z, G4double coul)