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

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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