source: trunk/source/processes/electromagnetic/polarisation/src/G4ePolarizedBremsstrahlungModel.cc

Last change on this file was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 5.3 KB
Line 
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: G4ePolarizedBremsstrahlungModel.cc,v 1.4 2007/11/01 17:32:34 schaelic Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4ePolarizedBremsstrahlungModel
35//
36// Author:        Karim Laihem
37//
38// Creation date: 12.03.2005
39//
40// Modifications:
41//    19-08-06 addapted to accomodate geant481 structure
42//
43//
44// Class Description:
45//
46//
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
52#include "G4ePolarizedBremsstrahlungModel.hh"
53#include "G4PolarizedBremsstrahlungCrossSection.hh"
54#include "G4ParticleChangeForLoss.hh"
55#include "G4PolarizationHelper.hh"
56
57G4ePolarizedBremsstrahlungModel::G4ePolarizedBremsstrahlungModel(const G4ParticleDefinition* p,
58                                               const G4String& nam)
59  : G4eBremsstrahlungModel(p,nam),
60    crossSectionCalculator(0)
61{
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66G4ePolarizedBremsstrahlungModel::~G4ePolarizedBremsstrahlungModel()
67{
68  if (crossSectionCalculator) delete crossSectionCalculator;
69}
70 
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72
73void G4ePolarizedBremsstrahlungModel::Initialise(const G4ParticleDefinition* p, 
74                                                         const G4DataVector& d)
75{
76  G4eBremsstrahlungModel::Initialise(p,d);
77  if (!crossSectionCalculator)
78    crossSectionCalculator = new G4PolarizedBremsstrahlungCrossSection();
79}
80
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
84
85void G4ePolarizedBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
86                                                        const G4MaterialCutsCouple* couple,
87                                                        const G4DynamicParticle* dp,
88                                                        G4double tmin,
89                                                        G4double maxEnergy)
90{
91  G4eBremsstrahlungModel::SampleSecondaries(vdp,couple,dp,tmin,maxEnergy);
92  G4int num = vdp->size();
93
94  if(num > 0) {
95    G4double lepEnergy0 = dp->GetKineticEnergy();
96    G4double gamEnergy1 = (*vdp)[0]->GetKineticEnergy();
97    G4double sintheta   = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
98    if (sintheta>1.) sintheta=1.;
99
100
101    G4StokesVector beamPol = dp->GetPolarization();
102
103    // determine interaction plane
104    G4ThreeVector  nInteractionFrame = 
105      G4PolarizationHelper::GetFrame(dp->GetMomentumDirection(), 
106                 fParticleChange->GetProposedMomentumDirection());
107
108    // transform polarization into interaction frame
109     beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
110
111    // calulcate polarization transfer
112    crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(),  // number of nucleons
113                                        GetCurrentElement()->GetZ(),
114                                        GetCurrentElement()->GetfCoulomb());
115    crossSectionCalculator->Initialize(lepEnergy0, gamEnergy1, sintheta,
116                                       beamPol, G4StokesVector::ZERO);
117
118    // deterimine final state polarization
119    G4StokesVector newBeamPol = crossSectionCalculator->GetPol2();
120    newBeamPol.RotateAz(nInteractionFrame,
121        fParticleChange->GetProposedMomentumDirection());
122    fParticleChange->ProposePolarization(newBeamPol);
123
124    if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized bremsstrahlung not supported!\n"; 
125    for (G4int i=0; i<num; i++) {
126      G4StokesVector photonPol = crossSectionCalculator->GetPol3();
127      photonPol.SetPhoton();
128      photonPol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
129      (*vdp)[i]->SetPolarization(photonPol.p1(),
130                                 photonPol.p2(),
131                                 photonPol.p3());
132    }
133  }
134  return;
135}
136// The emitted gamma energy is sampled using a parametrized formula
Note: See TracBrowser for help on using the repository browser.