Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ePolarizedBremsstrahlungModel.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//
29// GEANT4 Class file
30//
31//
32// File name: G4ePolarizedBremsstrahlungModel
33//
34// Author: Karim Laihem
35//
36// Creation date: 12.03.2005
37//
38// Modifications:
39// 19-08-06 addapted to accomodate geant481 structure
40//
41//
42// Class Description:
43//
44//
45// -------------------------------------------------------------------
46//
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
54
56 const G4ParticleDefinition* p, const G4String& nam)
57 : G4SeltzerBergerModel(p,nam),
58 crossSectionCalculator(nullptr)
59{
60}
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63
65{
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72 const G4DataVector& d)
73{
77}
78
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
82
83void G4ePolarizedBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
84 const G4MaterialCutsCouple* couple,
85 const G4DynamicParticle* dp,
86 G4double tmin,
87 G4double maxEnergy)
88{
89 G4SeltzerBergerModel::SampleSecondaries(vdp,couple,dp,tmin,maxEnergy);
90 G4int num = vdp->size();
91
92 if(num > 0) {
93 G4double lepEnergy0 = dp->GetKineticEnergy();
94 G4double gamEnergy1 = (*vdp)[0]->GetKineticEnergy();
95 G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
96 if (sintheta>1.) sintheta=1.;
97
98
99 G4StokesVector beamPol = dp->GetPolarization();
100
101 // determine interaction plane
102 G4ThreeVector nInteractionFrame =
105
106 // transform polarization into interaction frame
107 beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
108
109 // calulcate polarization transfer
110 crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
111 GetCurrentElement()->GetZ(),
112 GetCurrentElement()->GetfCoulomb());
113 crossSectionCalculator->Initialize(lepEnergy0, gamEnergy1, sintheta,
114 beamPol, G4StokesVector::ZERO);
115
116 // deterimine final state polarization
118 newBeamPol.RotateAz(nInteractionFrame,
121
122 if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized bremsstrahlung not supported!\n";
123 for (G4int i=0; i<num; i++) {
125 photonPol.SetPhoton();
126 photonPol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
127 (*vdp)[i]->SetPolarization(photonPol.p1(),
128 photonPol.p2(),
129 photonPol.p3());
130 }
131 }
132 return;
133}
134// The emitted gamma energy is sampled using a parametrized formula
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
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 &)
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
Definition: G4VEmModel.hh:495
void SetMaterial(G4double A, G4double Z, G4double coul)
virtual G4StokesVector GetPol3()
virtual G4StokesVector GetPol2()
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
G4ParticleChangeForLoss * fParticleChange
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ePolarizedBremsstrahlungModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="PolBrem")