Geant4  10.02.p02
G4PolarizedAnnihilationModel.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: G4PolarizedAnnihilationModel.cc 91742 2015-08-04 11:48:51Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4PolarizedAnnihilationModel
34 //
35 // Author: Andreas Schaelicke
36 //
37 // Creation date: 01.05.2005
38 //
39 // Modifications:
40 // 18-07-06 use newly calculated cross sections (P. Starovoitov)
41 // 21-08-06 update interface (A. Schaelicke)
42 // 17-11-06 add protection agaist e+ zero energy PostStep (V.Ivanchenko)
43 // 10-07-07 copied Initialise() method from G4eeToTwoGammaModel to provide a
44 // local ParticleChangeForGamma object and reduce overhead
45 // in SampleSecondaries() (A. Schaelicke)
46 //
47 //
48 // Class Description:
49 //
50 // Implementation of polarized gamma Annihilation scattering on free electron
51 //
52 
53 // -------------------------------------------------------------------
55 #include "G4PhysicalConstants.hh"
56 #include "G4PolarizationManager.hh"
57 #include "G4PolarizationHelper.hh"
58 #include "G4StokesVector.hh"
61 #include "G4TrackStatus.hh"
62 #include "G4Gamma.hh"
63 
65  const G4String& nam)
66  : G4eeToTwoGammaModel(p,nam),crossSectionCalculator(0),verboseLevel(0),gParticleChange(0),
67  gIsInitialised(false)
68 {
70 }
71 
73 {
75 }
76 
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79 
81  const G4DataVector&)
82 {
83  // G4eeToTwoGammaModel::Initialise(part,dv);
84  if(gIsInitialised) return;
86  gIsInitialised = true;
87 }
88 
90  const G4ParticleDefinition* pd,
91  G4double kinEnergy,
92  G4double cut,
93  G4double emax)
94 {
96  cut,emax);
97 
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 
111  G4double & valueX,
112  G4double & valueA,
113  G4double & valueT)
114 {
115  // *** calculate asymmetries
116  G4double gam = 1. + ene/electron_mass_c2;
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 
148 void 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
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 ***
178 
179  // obtain polarization of the beam
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)
191 
192  G4ParticleMomentum PositDirection = dp->GetMomentumDirection();
193 
194  // polar asymmetry:
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 
217  G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
218  <<"epsilmin DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
219  //check_range=false;
220  }
221 
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 
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  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
256  } while( treject < gmax*G4UniformRand() );
257 
258  //
259  // scattered Gamma angles. ( Z - axis along the parent positron)
260  //
261 
262  G4double cost = (epsil*gamap1-1.)/(epsil*sqg2m1);
263  G4double sint = std::sqrt((1.+cost)*(1.-cost));
264  G4double phi = 0.;
265  G4double beamTrans = std::sqrt(sqr(theBeamPolarization.p1()) + sqr(theBeamPolarization.p2()));
266  G4double targetTrans = std::sqrt(sqr(theTargetPolarization.p1()) + sqr(theTargetPolarization.p2()));
267 
268  // G4cout<<"phi dicing START"<<G4endl;
269  do{
270  phi = twopi * G4UniformRand();
272 
275  gdiced += 1.*(std::fabs(crossSectionCalculator->getVar(1))
276  + std::fabs(crossSectionCalculator->getVar(2)))*beamTrans*targetTrans;
277  gdiced += 1.*std::fabs(crossSectionCalculator->getVar(4))
278  *(std::fabs(theBeamPolarization.p3())*targetTrans + std::fabs(theTargetPolarization.p3())*beamTrans);
279 
282  gdist += crossSectionCalculator->getVar(1)*(std::cos(phi)*theBeamPolarization.p1()
283  + std::sin(phi)*theBeamPolarization.p2())
284  *(std::cos(phi)*theTargetPolarization.p1()
285  + std::sin(phi)*theTargetPolarization.p2());
286  gdist += crossSectionCalculator->getVar(2)*(std::cos(phi)*theBeamPolarization.p2()
287  - std::sin(phi)*theBeamPolarization.p1())
288  *(std::cos(phi)*theTargetPolarization.p2()
289  - std::sin(phi)*theTargetPolarization.p1());
290  gdist += crossSectionCalculator->getVar(4)
291  *(std::cos(phi)*theBeamPolarization.p3()*theTargetPolarization.p1()
292  + std::cos(phi)*theBeamPolarization.p1()*theTargetPolarization.p3()
293  + std::sin(phi)*theBeamPolarization.p3()*theTargetPolarization.p2()
294  + std::sin(phi)*theBeamPolarization.p2()*theTargetPolarization.p3());
295 
296  treject = gdist/gdiced;
297  //G4cout<<" treject = "<<treject<<" at phi = "<<phi<<G4endl;
298  if (treject>1.+1.e-10 || treject<0){
299  G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
300  <<" phi rejection does not work properly: "<<treject<<G4endl;
301  G4cout<<" gdiced = "<<gdiced<<G4endl;
302  G4cout<<" gdist = "<<gdist<<G4endl;
303  G4cout<<" epsil = "<<epsil<<G4endl;
304  }
305 
306  if (treject<1.e-3) {
307  G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
308  <<" phi rejection does not work properly: "<<treject<<"\n";
309  G4cout<<" gdiced="<<gdiced<<" gdist="<<gdist<<"\n";
310  G4cout<<" epsil = "<<epsil<<G4endl;
311  }
312 
313  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
314  } while( treject < G4UniformRand() );
315  // G4cout<<"phi dicing END"<<G4endl;
316 
317  G4double dirx = sint*std::cos(phi) , diry = sint*std::sin(phi) , dirz = cost;
318 
319  //
320  // kinematic of the created pair
321  //
322  G4double TotalAvailableEnergy = PositKinEnergy + 2*electron_mass_c2;
323  G4double Phot1Energy = epsil*TotalAvailableEnergy;
324  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
325 
326  // *** prepare calculation of polarization transfer ***
327  G4ThreeVector Phot1Direction (dirx, diry, dirz);
328 
329  // get interaction frame
330  G4ThreeVector nInteractionFrame =
331  G4PolarizationHelper::GetFrame(PositDirection,Phot1Direction);
332 
333  // define proper in-plane and out-of-plane component of initial spins
334  theBeamPolarization.InvRotateAz(nInteractionFrame,PositDirection);
335  theTargetPolarization.InvRotateAz(nInteractionFrame,PositDirection);
336 
337  // calculate spin transfere matrix
338 
340 
341  // **********************************************************************
342 
343  Phot1Direction.rotateUz(PositDirection);
344  // create G4DynamicParticle object for the particle1
346  Phot1Direction, Phot1Energy);
349  if (n1>1) {
350  G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "
351  <<epsil<<" is too large!!! \n"
352  <<"annihi pol1= "<<finalGamma1Polarization<<", ("<<n1<<")\n";
353  finalGamma1Polarization*=1./std::sqrt(n1);
354  }
355 
356  // define polarization of first final state photon
358  finalGamma1Polarization.RotateAz(nInteractionFrame,Phot1Direction);
362 
363  fvect->push_back(aParticle1);
364 
365 
366  // **********************************************************************
367 
368  G4double Eratio= Phot1Energy/Phot2Energy;
369  G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
370  G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio,
371  (PositP-dirz*Phot1Energy)/Phot2Energy);
372  Phot2Direction.rotateUz(PositDirection);
373  // create G4DynamicParticle object for the particle2
375  Phot2Direction, Phot2Energy);
376 
377  // define polarization of second final state photon
380  if (n2>1) {
381  G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<<epsil<<" is too large!!! \n";
382  G4cout<<"annihi pol2= "<<finalGamma2Polarization<<", ("<<n2<<")\n";
383 
384  finalGamma2Polarization*=1./std::sqrt(n2);
385  }
387  finalGamma2Polarization.RotateAz(nInteractionFrame,Phot2Direction);
391 
392  fvect->push_back(aParticle2);
393 }
const G4ThreeVector & GetPolarization() const
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double p2() const
static G4PolarizationManager * GetInstance()
virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1)
int G4int
Definition: G4Types.hh:78
G4double p3() const
static const G4StokesVector P3
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static const G4StokesVector P2
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cutEnergy=0., G4double maxEnergy=DBL_MAX)
G4ParticleChangeForGamma * gParticleChange
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax)
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
G4PolarizedAnnihilationCrossSection * crossSectionCalculator
static const double twopi
Definition: G4SIunits.hh:75
static const G4StokesVector P1
G4PolarizedAnnihilationModel(const G4ParticleDefinition *p=0, const G4String &nam="Polarized-Annihilation")
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4double p1() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
virtual void Initialize(G4double eps, G4double gamma, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
void ComputeAsymmetriesPerElectron(G4double gammaEnergy, G4double &valueX, G4double &valueA, G4double &valueT)
static const G4double emax
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4LogicalVolume * GetLogicalVolume() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
#define DBL_MIN
Definition: templates.hh:75
G4VPhysicalVolume * GetVolume() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ThreeVector G4ParticleMomentum
static const G4StokesVector ZERO
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134