Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PolarizedBremsstrahlungModel.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: G4PolarizedBremsstrahlungModel
31//
32// Author: Karim Laihem
33
35
39
41 const G4ParticleDefinition* p, const G4String& nam)
42 : G4SeltzerBergerModel(p, nam)
43 , fCrossSectionCalculator(nullptr)
44{}
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48{
49 if(fCrossSectionCalculator)
50 delete fCrossSectionCalculator;
51}
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 const G4DataVector& d)
56{
58 if(!fCrossSectionCalculator)
59 fCrossSectionCalculator = new G4PolarizedBremsstrahlungXS();
60}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63// The emitted gamma energy is sampled using a parametrized formula
65 std::vector<G4DynamicParticle*>* vdp, const G4MaterialCutsCouple* couple,
66 const G4DynamicParticle* dp, G4double tmin, G4double maxEnergy)
67{
68 G4SeltzerBergerModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
69 std::size_t num = vdp->size();
70
71 if(num > 0)
72 {
73 G4double lepEnergy0 = dp->GetKineticEnergy();
74 G4double gamEnergy1 = (*vdp)[0]->GetKineticEnergy();
75 G4double sintheta =
76 dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
77 if(sintheta > 1.)
78 sintheta = 1.;
79
81
82 // determine interaction plane
86
87 // transform polarization into interaction frame
88 beamPol.InvRotateAz(nInteractionFrame, dp->GetMomentumDirection());
89
90 // calculate polarization transfer
91 fCrossSectionCalculator->SetMaterial(
92 GetCurrentElement()->GetN(), // number of nucleons
93 GetCurrentElement()->GetZ(), GetCurrentElement()->GetfCoulomb());
94 fCrossSectionCalculator->Initialize(lepEnergy0, gamEnergy1, sintheta,
95 beamPol, G4StokesVector::ZERO);
96
97 // determine final state polarization
98 G4StokesVector newBeamPol = fCrossSectionCalculator->GetPol2();
99 newBeamPol.RotateAz(nInteractionFrame,
102
103 if(num != 1)
104 {
106 ed << num << " secondaries in polarized bremsstrahlung not supported!\n";
107 G4Exception("G4PolarizedBremsstrahlungModel::SampleSecondaries", "pol001",
108 JustWarning, ed);
109 }
110 for(std::size_t i = 0; i < num; ++i)
111 {
112 G4StokesVector photonPol = fCrossSectionCalculator->GetPol3();
113 photonPol.SetPhoton();
114 photonPol.RotateAz(nInteractionFrame, (*vdp)[i]->GetMomentumDirection());
115 (*vdp)[i]->SetPolarization(photonPol.p1(), photonPol.p2(),
116 photonPol.p3());
117 }
118 }
119 return;
120}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
Hep3Vector cross(const Hep3Vector &) const
double mag() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
const G4ThreeVector & GetProposedMomentumDirection() const
void ProposePolarization(const G4ThreeVector &dir)
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4PolarizedBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="PolBrem")
G4ParticleChangeForLoss * fParticleChange
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double cutEnergy, G4double maxEnergy) 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
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)=0
virtual G4StokesVector GetPol2()
virtual G4StokesVector GetPol3()
void SetMaterial(G4double A, G4double Z, G4double coul)