source: trunk/source/processes/electromagnetic/polarisation/src/G4PolarizedAnnihilationModel.cc @ 846

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 15.9 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: G4PolarizedAnnihilationModel.cc,v 1.6 2007/07/10 09:38:17 schaelic Exp $
27// GEANT4 tag $Name:  $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4PolarizedAnnihilationModel
35//
36// Author:        Andreas Schaelicke
37//
38// Creation date: 01.05.2005
39//
40// Modifications:
41// 18-07-06 use newly calculated cross sections (P. Starovoitov)
42// 21-08-06 update interface (A. Schaelicke)
43// 17-11-06 add protection agaist e+ zero energy PostStep (V.Ivanchenko)
44// 10-07-07 copied Initialise() method from G4eeToTwoGammaModel to provide a 
45//          local ParticleChangeForGamma object and reduce overhead
46//          in SampleSecondaries()  (A. Schaelicke)
47//
48//
49// Class Description:
50//
51// Implementation of polarized gamma Annihilation scattering on free electron
52//
53
54// -------------------------------------------------------------------
55#include "G4PolarizedAnnihilationModel.hh"
56#include "G4PolarizationManager.hh"
57#include "G4PolarizationHelper.hh"
58#include "G4StokesVector.hh"
59#include "G4PolarizedAnnihilationCrossSection.hh"
60#include "G4ParticleChangeForGamma.hh"
61#include "G4TrackStatus.hh"
62#include "G4Gamma.hh"
63
64G4PolarizedAnnihilationModel::G4PolarizedAnnihilationModel(const G4ParticleDefinition* p, 
65                                                           const G4String& nam)
66  : G4eeToTwoGammaModel(p,nam),crossSectionCalculator(0),gParticleChange(0),
67    gIsInitialised(false)
68{
69  crossSectionCalculator=new G4PolarizedAnnihilationCrossSection();
70}
71
72G4PolarizedAnnihilationModel::~G4PolarizedAnnihilationModel()
73{
74  if (crossSectionCalculator) delete crossSectionCalculator;
75}
76
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79
80void G4PolarizedAnnihilationModel::Initialise(const G4ParticleDefinition*,
81                                     const G4DataVector&)
82{
83  //  G4eeToTwoGammaModel::Initialise(part,dv);
84  if(gIsInitialised) return;
85
86  if(pParticleChange)
87    gParticleChange =
88      reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
89  else
90    gParticleChange = new G4ParticleChangeForGamma();
91
92  gIsInitialised = true;
93}
94
95G4double G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron(
96                                const G4ParticleDefinition* pd,
97                                      G4double kinEnergy, 
98                                      G4double cut,
99                                      G4double emax)
100{
101  G4double xs = G4eeToTwoGammaModel::ComputeCrossSectionPerElectron(pd,kinEnergy,
102                                                                cut,emax);
103
104  G4double polzz = theBeamPolarization.z()*theTargetPolarization.z();
105  G4double poltt = theBeamPolarization.x()*theTargetPolarization.x() 
106                 + theBeamPolarization.y()*theTargetPolarization.y();
107  if (polzz!=0 || poltt!=0) {
108    G4double xval,lasym,tasym;
109    ComputeAsymmetriesPerElectron(kinEnergy,xval,lasym,tasym);
110    xs*=(1.+polzz*lasym+poltt*tasym);
111  }
112
113  return xs;
114}
115
116void G4PolarizedAnnihilationModel::ComputeAsymmetriesPerElectron(G4double ene,
117                                               G4double & valueX,
118                                               G4double & valueA,
119                                               G4double & valueT)
120{
121  // *** calculate asymmetries
122  G4double gam = 1. + ene/electron_mass_c2;
123  G4double xs0=crossSectionCalculator->TotalXSection(0.,1.,gam,
124                                               G4StokesVector::ZERO,
125                                               G4StokesVector::ZERO);
126  G4double xsA=crossSectionCalculator->TotalXSection(0.,1.,gam,
127                                               G4StokesVector::P3,
128                                               G4StokesVector::P3);
129  G4double xsT1=crossSectionCalculator->TotalXSection(0.,1.,gam,
130                                               G4StokesVector::P1,
131                                               G4StokesVector::P1);
132  G4double xsT2=crossSectionCalculator->TotalXSection(0.,1.,gam,
133                                               G4StokesVector::P2,
134                                               G4StokesVector::P2);
135  G4double xsT=0.5*(xsT1+xsT2);
136 
137  valueX=xs0;
138  valueA=xsA/xs0-1.;
139  valueT=xsT/xs0-1.;
140  //  G4cout<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
141  if ( (valueA < -1) || (1 < valueA)) {
142    G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
143    G4cout<< " something wrong in total cross section calculation (valueA)\n";
144    G4cout<<"*********** LONG "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
145  }
146  if ( (valueT < -1) || (1 < valueT)) {
147    G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
148    G4cout<< " something wrong in total cross section calculation (valueT)\n";
149    G4cout<<"****** TRAN "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
150  }
151}
152
153
154void G4PolarizedAnnihilationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
155                                                     const G4MaterialCutsCouple* /*couple*/,
156                                                     const G4DynamicParticle* dp,
157                                                     G4double /*tmin*/,
158                                                     G4double /*maxEnergy*/) 
159{
160//   G4ParticleChangeForGamma*  gParticleChange
161//     = dynamic_cast<G4ParticleChangeForGamma*>(pParticleChange);
162  const G4Track * aTrack = gParticleChange->GetCurrentTrack();
163
164  // kill primary
165  gParticleChange->SetProposedKineticEnergy(0.);
166  gParticleChange->ProposeTrackStatus(fStopAndKill);
167
168  // V.Ivanchenko add protection against zero kin energy
169  G4double PositKinEnergy = dp->GetKineticEnergy();
170
171  if(PositKinEnergy < DBL_MIN) {
172
173    G4double cosTeta = 2.*G4UniformRand()-1.;
174    G4double sinTeta = std::sqrt((1.0 - cosTeta)*(1.0 + cosTeta));
175    G4double phi     = twopi * G4UniformRand();
176    G4ThreeVector dir(sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
177    fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(), dir, electron_mass_c2));
178    fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(),-dir, electron_mass_c2));
179    return;
180  }
181
182  // *** obtain and save target and beam polarization ***
183  G4PolarizationManager * polarizationManager = G4PolarizationManager::GetInstance();
184
185  // obtain polarization of the beam
186  theBeamPolarization = aTrack->GetPolarization();
187
188  // obtain polarization of the media
189  G4VPhysicalVolume*  aPVolume  = aTrack->GetVolume();
190  G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
191  const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
192  theTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
193
194  // transfer target electron polarization in frame of positron
195  if (targetIsPolarized)
196      theTargetPolarization.rotateUz(dp->GetMomentumDirection());
197 
198  G4ParticleMomentum PositDirection = dp->GetMomentumDirection();
199
200  // polar asymmetry:
201  G4double polarization = theBeamPolarization.p3()*theTargetPolarization.p3();
202
203  G4double gamam1 = PositKinEnergy/electron_mass_c2;
204  G4double gama   = gamam1+1. , gamap1 = gamam1+2.;
205  G4double sqgrate = std::sqrt(gamam1/gamap1)/2. , sqg2m1 = std::sqrt(gamam1*gamap1);
206
207  // limits of the energy sampling
208  G4double epsilmin = 0.5 - sqgrate , epsilmax = 0.5 + sqgrate;
209  G4double epsilqot = epsilmax/epsilmin;
210 
211  //
212  // sample the energy rate of the created gammas
213  // note: for polarized partices, the actual dicing strategy
214  //       will depend on the energy, and the degree of polarization !!
215  //
216  G4double epsil;
217  G4double gmax=1. + std::fabs(polarization); // crude estimate
218
219  G4bool check_range=true;
220
221  crossSectionCalculator->Initialize(epsilmin, gama, 0.,  theBeamPolarization, theTargetPolarization);
222  if (crossSectionCalculator->DiceEpsilon()<0) {
223    G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
224          <<"epsilmin DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
225    check_range=false;
226  }
227
228  crossSectionCalculator->Initialize(epsilmax, gama, 0.,  theBeamPolarization, theTargetPolarization);
229  if (crossSectionCalculator->DiceEpsilon()<0) {
230    G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
231          <<"epsilmax DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
232    check_range=false;
233  }
234
235  G4int ncount=0;
236  G4double trejectmax=0.;
237  G4double treject;
238
239
240  do {
241    //
242    epsil = epsilmin*std::pow(epsilqot,G4UniformRand());
243
244    crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,1);
245
246    treject = crossSectionCalculator->DiceEpsilon(); 
247    treject*=epsil;
248
249    if (treject>gmax  || treject<0.) 
250      G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
251            <<" eps ("<<epsil<<") rejection does not work properly: "<<treject<<G4endl;
252    ++ncount;
253    if (treject>trejectmax) trejectmax=treject;
254    if (ncount>1000) {
255      G4cout<<"WARNING  in PolarizedAnnihilationPS::PostStepDoIt\n"
256            <<"eps dicing very inefficient ="<<trejectmax/gmax
257            <<", "<<treject/gmax<<".  For secondary energy = "<<epsil<<"   "<<ncount<<G4endl;
258      break;
259    }
260
261  } while( treject < gmax*G4UniformRand() );
262
263  //
264  // scattered Gamma angles. ( Z - axis along the parent positron)
265  //
266   
267  G4double cost = (epsil*gamap1-1.)/(epsil*sqg2m1);
268  G4double sint = std::sqrt((1.+cost)*(1.-cost));
269  G4double phi  = 0.;
270  G4double   beamTrans = std::sqrt(sqr(theBeamPolarization.p1()) + sqr(theBeamPolarization.p2()));
271  G4double targetTrans = std::sqrt(sqr(theTargetPolarization.p1()) + sqr(theTargetPolarization.p2()));
272
273  //  G4cout<<"phi dicing START"<<G4endl;
274  do{
275    phi  = twopi * G4UniformRand();
276    crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,2);
277
278    G4double gdiced =crossSectionCalculator->getVar(0);
279    gdiced += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
280    gdiced += 1.*(std::fabs(crossSectionCalculator->getVar(1)) 
281                  + std::fabs(crossSectionCalculator->getVar(2)))*beamTrans*targetTrans;
282    gdiced += 1.*std::fabs(crossSectionCalculator->getVar(4))
283      *(std::fabs(theBeamPolarization.p3())*targetTrans + std::fabs(theTargetPolarization.p3())*beamTrans);
284
285    G4double gdist = crossSectionCalculator->getVar(0);
286    gdist += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
287    gdist += crossSectionCalculator->getVar(1)*(std::cos(phi)*theBeamPolarization.p1() 
288                                                + std::sin(phi)*theBeamPolarization.p2())
289                                              *(std::cos(phi)*theTargetPolarization.p1() 
290                                                + std::sin(phi)*theTargetPolarization.p2());
291    gdist += crossSectionCalculator->getVar(2)*(std::cos(phi)*theBeamPolarization.p2() 
292                                                - std::sin(phi)*theBeamPolarization.p1())
293                                              *(std::cos(phi)*theTargetPolarization.p2() 
294                                                - std::sin(phi)*theTargetPolarization.p1());
295    gdist += crossSectionCalculator->getVar(4)
296      *(std::cos(phi)*theBeamPolarization.p3()*theTargetPolarization.p1()
297        + std::cos(phi)*theBeamPolarization.p1()*theTargetPolarization.p3() 
298        + std::sin(phi)*theBeamPolarization.p3()*theTargetPolarization.p2() 
299        + std::sin(phi)*theBeamPolarization.p2()*theTargetPolarization.p3());
300
301    treject = gdist/gdiced;
302    //G4cout<<" treject = "<<treject<<" at phi = "<<phi<<G4endl;
303     if (treject>1.+1.e-10 || treject<0){
304       G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
305             <<" phi rejection does not work properly: "<<treject<<G4endl;
306       G4cout<<" gdiced = "<<gdiced<<G4endl;
307       G4cout<<" gdist = "<<gdist<<G4endl;
308       G4cout<<" epsil = "<<epsil<<G4endl;
309     }
310     
311     if (treject<1.e-3) {
312       G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
313            <<" phi rejection does not work properly: "<<treject<<"\n";
314       G4cout<<" gdiced="<<gdiced<<"   gdist="<<gdist<<"\n";
315       G4cout<<" epsil = "<<epsil<<G4endl;
316     }
317
318  } while( treject < G4UniformRand() );
319  //  G4cout<<"phi dicing END"<<G4endl;
320
321  G4double dirx = sint*std::cos(phi) , diry = sint*std::sin(phi) , dirz = cost;
322
323  //
324  // kinematic of the created pair
325  //
326  G4double TotalAvailableEnergy = PositKinEnergy + 2*electron_mass_c2;
327  G4double Phot1Energy = epsil*TotalAvailableEnergy;
328  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
329
330  // *** prepare calculation of polarization transfer ***
331  G4ThreeVector Phot1Direction (dirx, diry, dirz);
332
333  // get interaction frame
334  G4ThreeVector  nInteractionFrame = 
335    G4PolarizationHelper::GetFrame(PositDirection,Phot1Direction);
336     
337  // define proper in-plane and out-of-plane component of initial spins
338  theBeamPolarization.InvRotateAz(nInteractionFrame,PositDirection);
339  theTargetPolarization.InvRotateAz(nInteractionFrame,PositDirection);
340
341  // calculate spin transfere matrix
342
343  crossSectionCalculator->Initialize(epsil,gama,phi,theBeamPolarization,theTargetPolarization,2);
344
345  // **********************************************************************
346
347  Phot1Direction.rotateUz(PositDirection);   
348  // create G4DynamicParticle object for the particle1 
349  G4DynamicParticle* aParticle1= new G4DynamicParticle (G4Gamma::Gamma(),
350                                                        Phot1Direction, Phot1Energy);
351  finalGamma1Polarization=crossSectionCalculator->GetPol2();
352  G4double n1=finalGamma1Polarization.mag2();
353  if (n1>1) {
354    G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "
355          <<epsil<<" is too large!!! \n"
356          <<"annihi pol1= "<<finalGamma1Polarization<<", ("<<n1<<")\n";
357    finalGamma1Polarization+=1./std::sqrt(n1);
358  }
359
360  // define polarization of first final state photon
361  finalGamma1Polarization.SetPhoton();
362  finalGamma1Polarization.RotateAz(nInteractionFrame,Phot1Direction);
363  aParticle1->SetPolarization(finalGamma1Polarization.p1(),
364                              finalGamma1Polarization.p2(),
365                              finalGamma1Polarization.p3());
366
367  fvect->push_back(aParticle1);
368
369
370  // **********************************************************************
371
372  G4double Eratio= Phot1Energy/Phot2Energy;
373  G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
374  G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio,
375                                (PositP-dirz*Phot1Energy)/Phot2Energy); 
376  Phot2Direction.rotateUz(PositDirection); 
377  // create G4DynamicParticle object for the particle2
378  G4DynamicParticle* aParticle2= new G4DynamicParticle (G4Gamma::Gamma(),
379                                                        Phot2Direction, Phot2Energy);
380
381  // define polarization of second final state photon
382  finalGamma2Polarization=crossSectionCalculator->GetPol3();
383  G4double n2=finalGamma2Polarization.mag2();
384  if (n2>1) {
385    G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<<epsil<<" is too large!!! \n";
386    G4cout<<"annihi pol2= "<<finalGamma2Polarization<<", ("<<n2<<")\n";
387   
388    finalGamma2Polarization+=1./std::sqrt(n2);
389  }
390  finalGamma2Polarization.SetPhoton();
391  finalGamma2Polarization.RotateAz(nInteractionFrame,Phot2Direction);
392  aParticle2->SetPolarization(finalGamma2Polarization.p1(),
393                              finalGamma2Polarization.p2(),
394                              finalGamma2Polarization.p3());
395
396  fvect->push_back(aParticle2);
397}
Note: See TracBrowser for help on using the repository browser.