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