Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4ePolarizedBremsstrahlungModel
34//
35// Author: Karim Laihem
36//
37// Creation date: 12.03.2005
38//
39// Modifications:
40// 19-08-06 addapted to accomodate geant481 structure
41//
42//
43// Class Description:
44//
45//
46// -------------------------------------------------------------------
47//
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
55
57 const G4String& nam)
59 crossSectionCalculator(0)
60{
61}
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66{
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71
73 const G4DataVector& d)
74{
78}
79
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
83
84void G4ePolarizedBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
85 const G4MaterialCutsCouple* couple,
86 const G4DynamicParticle* dp,
87 G4double tmin,
88 G4double maxEnergy)
89{
90 G4eBremsstrahlungModel::SampleSecondaries(vdp,couple,dp,tmin,maxEnergy);
91 G4int num = vdp->size();
92
93 if(num > 0) {
94 G4double lepEnergy0 = dp->GetKineticEnergy();
95 G4double gamEnergy1 = (*vdp)[0]->GetKineticEnergy();
96 G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
97 if (sintheta>1.) sintheta=1.;
98
99
100 G4StokesVector beamPol = dp->GetPolarization();
101
102 // determine interaction plane
103 G4ThreeVector nInteractionFrame =
106
107 // transform polarization into interaction frame
108 beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
109
110 // calulcate polarization transfer
111 crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
112 GetCurrentElement()->GetZ(),
113 GetCurrentElement()->GetfCoulomb());
114 crossSectionCalculator->Initialize(lepEnergy0, gamEnergy1, sintheta,
115 beamPol, G4StokesVector::ZERO);
116
117 // deterimine final state polarization
119 newBeamPol.RotateAz(nInteractionFrame,
122
123 if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized bremsstrahlung not supported!\n";
124 for (G4int i=0; i<num; i++) {
126 photonPol.SetPhoton();
127 photonPol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
128 (*vdp)[i]->SetPolarization(photonPol.p1(),
129 photonPol.p2(),
130 photonPol.p3());
131 }
132 }
133 return;
134}
135// The emitted gamma energy is sampled using a parametrized formula
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT 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 &)
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:391
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)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4ParticleChangeForLoss * fParticleChange
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ePolarizedBremsstrahlungModel(const G4ParticleDefinition *p=0, const G4String &nam="PolBrem")
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)