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
G4PolarizedGammaConversionModel.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: G4PolarizedGammaConversionModel
31//
32// Author: Karim Laihem
33//
34// Class Description:
35// Implementation of gamma conversion to e+e- in the field of a nucleus
36// including polarization transfer
37
39
40#include "G4DataVector.hh"
41#include "G4DynamicParticle.hh"
45#include "G4StokesVector.hh"
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 const G4ParticleDefinition* pd, const G4String& nam)
50 : G4BetheHeitlerModel(pd, nam)
51 , fCrossSectionCalculator(nullptr)
52{}
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56{
57 if(fCrossSectionCalculator)
58 delete fCrossSectionCalculator;
59}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 const G4DataVector& dv)
64{
66 if(!fCrossSectionCalculator)
67 fCrossSectionCalculator = new G4PolarizedGammaConversionXS();
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 std::vector<G4DynamicParticle*>* vdp, const G4MaterialCutsCouple* couple,
73 const G4DynamicParticle* dp, G4double tmin, G4double maxEnergy)
74{
75 G4BetheHeitlerModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
76
77 if(vdp && vdp->size() > 0)
78 {
79 G4double gamEnergy0 = dp->GetKineticEnergy();
80 G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
81 G4double sintheta =
82 dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
83 if(sintheta > 1.)
84 sintheta = 1.;
85
87 beamPol.SetPhoton();
88
89 // determine interaction plane
91 dp->GetMomentumDirection(), (*vdp)[0]->GetMomentumDirection());
92
93 // transform polarization into interaction frame
94 beamPol.InvRotateAz(nInteractionFrame, dp->GetMomentumDirection());
95
96 // calculate polarization transfer
97 fCrossSectionCalculator->SetMaterial(
98 GetCurrentElement()->GetN(), // number of nucleons
99 GetCurrentElement()->GetZ(), GetCurrentElement()->GetfCoulomb());
100 fCrossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
101 beamPol, G4StokesVector::ZERO);
102
103 // determine final state polarization
104 G4StokesVector lep1Pol = fCrossSectionCalculator->GetPol2();
105 lep1Pol.RotateAz(nInteractionFrame, (*vdp)[0]->GetMomentumDirection());
106 (*vdp)[0]->SetPolarization(lep1Pol.p1(), lep1Pol.p2(), lep1Pol.p3());
107
108 size_t num = vdp->size();
109 if(num != 2)
110 {
112 ed << " WARNING " << num
113 << " secondaries in polarized pairproduction not supported!\n";
114 G4Exception("G4PolarizedGammaConversionModel::SampleSecondaries",
115 "pol018", JustWarning, ed);
116 }
117 for(size_t i = 1; i < num; ++i)
118 {
119 G4StokesVector lep2Pol = fCrossSectionCalculator->GetPol3();
120 lep2Pol.RotateAz(nInteractionFrame, (*vdp)[i]->GetMomentumDirection());
121 (*vdp)[i]->SetPolarization(lep2Pol.p1(), lep2Pol.p2(), lep2Pol.p3());
122 }
123 }
124}
@ 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
Hep3Vector cross(const Hep3Vector &) const
double mag() const
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4PolarizedGammaConversionModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="polConv")
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void Initialize(G4double eps, G4double X, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, 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)