Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PolarizedAnnihilationModel.cc
이 파일의 문서화 페이지로 가기
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 109176 2018-04-03 06:53:39Z 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 {
72 }
73 
75 {
77 }
78 
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83  const G4DataVector& dv)
84 {
86  if(gParticleChange) { return; }
88 }
89 
90 G4double
92 {
93  // cross section from base model
95 
99  if (polzz!=0 || poltt!=0) {
100  G4double xval,lasym,tasym;
101  ComputeAsymmetriesPerElectron(kinEnergy,xval,lasym,tasym);
102  xs*=(1.+polzz*lasym+poltt*tasym);
103  }
104 
105  return xs;
106 }
107 
109  G4double & valueX,
110  G4double & valueA,
111  G4double & valueT)
112 {
113  // *** calculate asymmetries
114  G4double gam = 1. + ene/electron_mass_c2;
127  G4double xsT=0.5*(xsT1+xsT2);
128 
129  valueX=xs0;
130  valueA=xsA/xs0-1.;
131  valueT=xsT/xs0-1.;
132  // G4cout<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
133  if ( (valueA < -1) || (1 < valueA)) {
134  G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
135  G4cout<< " something wrong in total cross section calculation (valueA)\n";
136  G4cout<< " LONG: "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
137  }
138  if ( (valueT < -1) || (1 < valueT)) {
139  G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
140  G4cout<< " something wrong in total cross section calculation (valueT)\n";
141  G4cout<< " TRAN: "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<" energy = "<<gam<<G4endl;
142  }
143 }
144 
145 void G4PolarizedAnnihilationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
146  const G4MaterialCutsCouple*,
147  const G4DynamicParticle* dp,
148  G4double, G4double)
149 {
150  const G4Track * aTrack = gParticleChange->GetCurrentTrack();
151 
152  // kill primary
155 
156  // V.Ivanchenko add protection against zero kin energy
157  G4double PositKinEnergy = dp->GetKineticEnergy();
158 
159  if(PositKinEnergy == 0.0) {
160 
161  G4double cosTeta = 2.*G4UniformRand()-1.;
162  G4double sinTeta = std::sqrt((1.0 - cosTeta)*(1.0 + cosTeta));
163  G4double phi = twopi * G4UniformRand();
164  G4ThreeVector dir(sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
165  fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(), dir, electron_mass_c2));
166  fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(),-dir, electron_mass_c2));
167  return;
168  }
169 
170  // *** obtain and save target and beam polarization ***
172 
173  // obtain polarization of the beam
175 
176  // obtain polarization of the media
177  G4VPhysicalVolume* aPVolume = aTrack->GetVolume();
178  G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
179  const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
180  theTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
181 
182  if (verboseLevel >= 1) {
183  G4cout << "G4PolarizedComptonModel::SampleSecondaries in "
184  << aLVolume->GetName() << G4endl;
185  }
186 
187  // transfer target electron polarization in frame of positron
188  if (targetIsPolarized)
190 
191  G4ParticleMomentum PositDirection = dp->GetMomentumDirection();
192 
193  // polar asymmetry:
195 
196  G4double gamam1 = PositKinEnergy/electron_mass_c2;
197  G4double gama = gamam1+1. , gamap1 = gamam1+2.;
198  G4double sqgrate = std::sqrt(gamam1/gamap1)/2. , sqg2m1 = std::sqrt(gamam1*gamap1);
199 
200  // limits of the energy sampling
201  G4double epsilmin = 0.5 - sqgrate , epsilmax = 0.5 + sqgrate;
202  G4double epsilqot = epsilmax/epsilmin;
203 
204  //
205  // sample the energy rate of the created gammas
206  // note: for polarized partices, the actual dicing strategy
207  // will depend on the energy, and the degree of polarization !!
208  //
209  G4double epsil;
210  G4double gmax=1. + std::fabs(polarization); // crude estimate
211 
212  //G4bool check_range=true;
213 
216  G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
217  <<"epsilmin DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
218  //check_range=false;
219  }
220 
223  G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
224  <<"epsilmax DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
225  //check_range=false;
226  }
227 
228  G4int ncount=0;
229  G4double trejectmax=0.;
230  G4double treject;
231 
232 
233  do {
234  //
235  epsil = epsilmin*std::pow(epsilqot,G4UniformRand());
236 
238 
239  treject = crossSectionCalculator->DiceEpsilon();
240  treject*=epsil;
241 
242  if (treject>gmax || treject<0.)
243  G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
244  <<" eps ("<<epsil<<") rejection does not work properly: "<<treject<<G4endl;
245  ++ncount;
246  if (treject>trejectmax) trejectmax=treject;
247  if (ncount>1000) {
248  G4cout<<"WARNING in PolarizedAnnihilationPS::PostStepDoIt\n"
249  <<"eps dicing very inefficient ="<<trejectmax/gmax
250  <<", "<<treject/gmax<<". For secondary energy = "<<epsil<<" "<<ncount<<G4endl;
251  break;
252  }
253 
254  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
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();
271 
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 
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  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
313  } while( treject < G4UniformRand() );
314  // G4cout<<"phi dicing END"<<G4endl;
315 
316  G4double dirx = sint*std::cos(phi) , diry = sint*std::sin(phi) , dirz = cost;
317 
318  //
319  // kinematic of the created pair
320  //
321  G4double TotalAvailableEnergy = PositKinEnergy + 2*electron_mass_c2;
322  G4double Phot1Energy = epsil*TotalAvailableEnergy;
323  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
324 
325  // *** prepare calculation of polarization transfer ***
326  G4ThreeVector Phot1Direction (dirx, diry, dirz);
327 
328  // get interaction frame
329  G4ThreeVector nInteractionFrame =
330  G4PolarizationHelper::GetFrame(PositDirection,Phot1Direction);
331 
332  // define proper in-plane and out-of-plane component of initial spins
333  theBeamPolarization.InvRotateAz(nInteractionFrame,PositDirection);
334  theTargetPolarization.InvRotateAz(nInteractionFrame,PositDirection);
335 
336  // calculate spin transfere matrix
337 
339 
340  // **********************************************************************
341 
342  Phot1Direction.rotateUz(PositDirection);
343  // create G4DynamicParticle object for the particle1
345  Phot1Direction, Phot1Energy);
348  if (n1>1) {
349  G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "
350  <<epsil<<" is too large!!! \n"
351  <<"annihi pol1= "<<finalGamma1Polarization<<", ("<<n1<<")\n";
352  finalGamma1Polarization*=1./std::sqrt(n1);
353  }
354 
355  // define polarization of first final state photon
357  finalGamma1Polarization.RotateAz(nInteractionFrame,Phot1Direction);
361 
362  fvect->push_back(aParticle1);
363 
364 
365  // **********************************************************************
366 
367  G4double Eratio= Phot1Energy/Phot2Energy;
368  G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
369  G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio,
370  (PositP-dirz*Phot1Energy)/Phot2Energy);
371  Phot2Direction.rotateUz(PositDirection);
372  // create G4DynamicParticle object for the particle2
374  Phot2Direction, Phot2Energy);
375 
376  // define polarization of second final state photon
379  if (n2>1) {
380  G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<<epsil<<" is too large!!! \n";
381  G4cout<<"annihi pol2= "<<finalGamma2Polarization<<", ("<<n2<<")\n";
382 
383  finalGamma2Polarization*=1./std::sqrt(n2);
384  }
386  finalGamma2Polarization.RotateAz(nInteractionFrame,Phot2Direction);
390 
391  fvect->push_back(aParticle2);
392 }
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4LogicalVolume * GetLogicalVolume() const
static const G4StokesVector P2
virtual void Initialize(G4double eps, G4double gamma, G4double phi, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0) override
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
const char * p
Definition: xmltok.h:285
G4double p3() const
G4PolarizedAnnihilationModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="Polarized-Annihilation")
const G4ThreeVector & GetPolarization() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double z() const
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
G4double p1() const
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy)
virtual G4double ComputeCrossSectionPerElectron(G4double kinEnergy) final
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
TString part[npart]
Definition: Style.C:32
static constexpr double electron_mass_c2
static const G4StokesVector P3
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
void ComputeAsymmetriesPerElectron(G4double gammaEnergy, G4double &valueX, G4double &valueA, G4double &valueT)
double mag2() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4Track * GetCurrentTrack() const
G4PolarizedAnnihilationCrossSection * crossSectionCalculator
static const G4StokesVector ZERO
G4ParticleChangeForGamma * gParticleChange
int G4int
Definition: G4Types.hh:78
TDirectory * dir
virtual G4double TotalXSection(G4double xmin, G4double xmax, G4double y, const G4StokesVector &pol0, const G4StokesVector &pol1) override
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) final
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
double x() const
static const G4StokesVector P1
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
G4double p2() const
T sqr(const T &x)
Definition: templates.hh:145
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:131
double y() const
static G4PolarizationManager * GetInstance()
void ProposeTrackStatus(G4TrackStatus status)
const G4String & GetName() const
G4VPhysicalVolume * GetVolume() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
bool IsPolarized(G4LogicalVolume *lVol) const