Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/src/G4LivermorePolarizedRayleighModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // Author: Sebastien Incerti
 28 //         30 October 2008
 29 //         on base of G4LowEnergyPolarizedRayleigh developed by R. Capra
 30 //
 31 // History:
 32 // --------
 33 // 02 May 2009   S Incerti as V. Ivanchenko proposed in G4LivermoreRayleighModel.cc
 34 //
 35 // Cleanup initialisation and generation of secondaries:
 36 //                  - apply internal high-energy limit only in constructor 
 37 //                  - do not apply low-energy limit (default is 0)
 38 //                  - remove GetMeanFreePath method and table
 39 //                  - remove initialisation of element selector 
 40 //                  - use G4ElementSelector
 41 
 42 #include "G4LivermorePolarizedRayleighModel.hh"
 43 #include "G4PhysicalConstants.hh"
 44 #include "G4SystemOfUnits.hh"
 45 #include "G4LogLogInterpolation.hh"
 46 #include "G4CompositeEMDataSet.hh"
 47 #include "G4AutoLock.hh"
 48 
 49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 50 
 51 using namespace std;
 52 namespace { G4Mutex LivermorePolarizedRayleighModelMutex = G4MUTEX_INITIALIZER; }
 53 
 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 55 
 56 G4PhysicsFreeVector* G4LivermorePolarizedRayleighModel::dataCS[] = {nullptr};
 57 G4PhysicsFreeVector* G4LivermorePolarizedRayleighModel::formFactorData[] = {nullptr};
 58 
 59 G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel(const G4ParticleDefinition*,
 60                    const G4String& nam)
 61   :G4VEmModel(nam),fParticleChange(nullptr),isInitialised(false)
 62 {
 63   lowEnergyLimit = 250 * CLHEP::eV;
 64   //
 65   verboseLevel= 0;
 66   // Verbosity scale:
 67   // 0 = nothing 
 68   // 1 = warning for energy non-conservation 
 69   // 2 = details of energy budget
 70   // 3 = calculation of cross sections, file openings, sampling of atoms
 71   // 4 = entering in methods
 72 
 73   if(verboseLevel > 0) {
 74     G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl
 75            << "Energy range: " << LowEnergyLimit() / CLHEP::eV << " eV - "
 76      << HighEnergyLimit() / CLHEP::GeV << " GeV"
 77            << G4endl;
 78   }
 79 }
 80 
 81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 82 
 83 G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel()
 84 {  
 85   if(IsMaster()) {
 86     for(G4int i=0; i<maxZ; ++i) {
 87       if(dataCS[i]) { 
 88   delete dataCS[i];
 89   dataCS[i] = nullptr;
 90   delete formFactorData[i];
 91   formFactorData[i] = nullptr; 
 92       }
 93     }
 94   }
 95 }
 96 
 97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 98 
 99 void G4LivermorePolarizedRayleighModel::Initialise(const G4ParticleDefinition* particle,
100                                        const G4DataVector& cuts)
101 {
102   // Rayleigh process:                      The Quantum Theory of Radiation
103   //                                        W. Heitler,       Oxford at the Clarendon Press, Oxford (1954)       
104   // Scattering function:                   A simple model of photon transport
105   //                                        D.E. Cullen,      Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510               
106   // Polarization of the outcoming photon:  Beam test of a prototype detector array for the PoGO astronomical hard 
107   //                                        X-ray/soft gamma-ray polarimeter
108   //                                        T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168 
109   if (verboseLevel > 3)
110     G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
111 
112   if(IsMaster()) {
113     
114     // Initialise element selector
115     InitialiseElementSelectors(particle, cuts);
116     
117     // Access to elements
118     const char* path = G4FindDataDir("G4LEDATA");
119     auto elmTable = G4Element::GetElementTable();
120     for (auto const & elm : *elmTable) {
121       G4int Z = std::min(elm->GetZasInt(), maxZ);
122       if( nullptr == dataCS[Z] ) { ReadData(Z, path); }
123     }
124   }
125   
126   if(isInitialised) { return; }
127   fParticleChange = GetParticleChangeForGamma();
128   isInitialised = true;
129 }
130 
131 
132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133 
134 void G4LivermorePolarizedRayleighModel::InitialiseLocal(
135                 const G4ParticleDefinition*, G4VEmModel* masterModel)
136 {
137   SetElementSelectors(masterModel->GetElementSelectors());
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
142 void G4LivermorePolarizedRayleighModel::ReadData(std::size_t Z, const char* path)
143 {
144   if (verboseLevel > 1) {
145     G4cout << "Calling ReadData() of G4LivermoreRayleighModel" 
146      << G4endl;
147   }
148   
149   if(nullptr != dataCS[Z]) { return; }
150   
151   const char* datadir = path;
152   
153   if(nullptr == datadir) 
154     {
155       datadir = G4FindDataDir("G4LEDATA");
156       if(nullptr == datadir) 
157   {
158     G4Exception("G4LivermoreRayleighModelModel::ReadData()","em0006",
159           FatalException,
160           "Environment variable G4LEDATA not defined");
161     return;
162   }
163     }
164   dataCS[Z] = new G4PhysicsFreeVector();
165   formFactorData[Z] = new G4PhysicsFreeVector();
166   
167   std::ostringstream ostCS;
168   ostCS << datadir << "/livermore/rayl/re-cs-" << Z <<".dat";
169   std::ifstream finCS(ostCS.str().c_str());
170   
171   if( !finCS .is_open() ) {
172     G4ExceptionDescription ed;
173     ed << "G4LivermorePolarizedRayleighModel data file <" << ostCS.str().c_str()
174        << "> is not opened!" << G4endl;
175     G4Exception("G4LivermorePolarizedRayleighModel::ReadData()","em0003",
176     FatalException,
177     ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
178     return;
179   } else {
180     if(verboseLevel > 3) { 
181       G4cout << "File " << ostCS.str() 
182        << " is opened by G4LivermoreRayleighModel" << G4endl;
183     }
184     dataCS[Z]->Retrieve(finCS, true);
185   }
186 
187   std::ostringstream ostFF;
188   ostFF << datadir << "/livermore/rayl/re-ff-" << Z <<".dat";
189   std::ifstream finFF(ostFF.str().c_str());
190   
191   if( !finFF.is_open() ) {
192     G4ExceptionDescription ed;
193     ed << "G4LivermorePolarizedRayleighModel data file <" << ostFF.str().c_str()
194        << "> is not opened!" << G4endl;
195     G4Exception("G4LivermorePolarizedRayleighModel::ReadData()","em0003",
196     FatalException,
197     ed,"G4LEDATA version should be G4EMLOW8.0 or later.");
198     return;
199   } else {
200     if(verboseLevel > 3) { 
201       G4cout << "File " << ostFF.str() 
202                << " is opened by G4LivermoreRayleighModel" << G4endl;
203     }
204     formFactorData[Z]->Retrieve(finFF, true);
205   } 
206 }
207  
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
209  
210 G4double G4LivermorePolarizedRayleighModel::ComputeCrossSectionPerAtom(
211                                         const G4ParticleDefinition*,
212           G4double GammaEnergy,
213           G4double Z, G4double,
214           G4double, G4double)
215 {
216   if (verboseLevel > 1) {
217     G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()" 
218      << G4endl;
219   }
220  
221   if(GammaEnergy < lowEnergyLimit) { return 0.0; }
222    
223   G4double xs = 0.0;
224    
225   G4int intZ = G4lrint(Z);
226   if(intZ < 1 || intZ > maxZ) { return xs; }
227    
228   G4PhysicsFreeVector* pv = dataCS[intZ];
229  
230   // if element was not initialised
231   // do initialisation safely for MT mode
232   if(nullptr == pv) { 
233     InitialiseForElement(0, intZ);
234     pv = dataCS[intZ];
235     if(nullptr == pv) { return xs; }
236   }
237  
238   G4int n = G4int(pv->GetVectorLength() - 1);
239   G4double e = GammaEnergy/MeV;
240   if(e >= pv->Energy(n)) {
241     xs = (*pv)[n]/(e*e);  
242   } else if(e >= pv->Energy(0)) {
243     xs = pv->Value(e)/(e*e);  
244   }
245   return xs;
246 }
247 
248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
249 
250 void G4LivermorePolarizedRayleighModel::SampleSecondaries(
251                 std::vector<G4DynamicParticle*>*,
252     const G4MaterialCutsCouple* couple,
253     const G4DynamicParticle* aDynamicGamma,
254     G4double, G4double)
255 {
256   if (verboseLevel > 3)
257     G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
258 
259   G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
260   
261   if (photonEnergy0 <= lowEnergyLimit)
262   {
263     fParticleChange->ProposeTrackStatus(fStopAndKill);
264     fParticleChange->SetProposedKineticEnergy(0.);
265     fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
266     return;
267   }
268 
269   G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
270 
271   // Select randomly one element in the current material
272   const G4ParticleDefinition* particle =  aDynamicGamma->GetDefinition();
273   const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
274   G4int Z = elm->GetZasInt();
275 
276   G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
277   G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
278   G4double beta = GeneratePolarizationAngle();
279  
280   // incomingPhoton reference frame:
281   // z = versor parallel to the incomingPhotonDirection
282   // x = versor parallel to the incomingPhotonPolarization
283   // y = defined as z^x
284  
285   // outgoingPhoton reference frame:
286   // z' = versor parallel to the outgoingPhotonDirection
287   // x' = defined as x-x*z'z' normalized
288   // y' = defined as z'^x' 
289   G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit()); 
290   G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
291   G4ThreeVector y(z.cross(x));
292  
293   // z' = std::cos(phi)*std::sin(theta) 
294   // x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
295   G4double xDir;
296   G4double yDir;
297   G4double zDir;
298   zDir=outcomingPhotonCosTheta;
299   xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
300   yDir=xDir;
301   xDir*=std::cos(outcomingPhotonPhi);
302   yDir*=std::sin(outcomingPhotonPhi);
303  
304   G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
305   G4ThreeVector xPrime(x.perpPart(zPrime).unit());
306   G4ThreeVector yPrime(zPrime.cross(xPrime));
307  
308   // outgoingPhotonPolarization is directed as 
309   // x' std::cos(beta) + y' std::sin(beta)
310   G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
311  
312   fParticleChange->ProposeMomentumDirection(zPrime);
313   fParticleChange->ProposePolarization(outcomingPhotonPolarization);
314   fParticleChange->SetProposedKineticEnergy(photonEnergy0); 
315 }
316 
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
318 
319 G4double G4LivermorePolarizedRayleighModel::GenerateCosTheta(G4double incomingPhotonEnergy, G4int Z) const
320 {
321   //  d sigma                                                                    k0
322   // --------- =  r0^2 * pi * F^2(x, Z) * ( 2 - sin^2 theta) * std::sin (theta), x = ---- std::sin(theta/2)
323   //  d theta                                                                    hc
324  
325   //  d sigma                                             k0          1 - y
326   // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + y^2), x = ---- std::sqrt ( ------- ), y = std::cos(theta)
327   //    d y                                               hc            2
328 
329   //              Z
330   // F(x, Z) ~ --------
331   //            a + bx
332   //
333   // The time to exit from the outer loop grows as ~ k0
334   // On pcgeant2 the time is ~ 1 s for k0 ~ 1 MeV on the oxygen element. A 100 GeV
335   // event will take ~ 10 hours.
336   //
337   // On the avarage the inner loop does 1.5 iterations before exiting
338   const G4double xxfact = CLHEP::cm/(CLHEP::h_Planck*CLHEP::c_light); 
339   const G4double xFactor = incomingPhotonEnergy*xxfact;
340 
341   G4double cosTheta;
342   G4double fCosTheta;
343   G4double x;
344   G4double fValue;
345 
346   if (incomingPhotonEnergy > 5.*CLHEP::MeV)
347   {
348     cosTheta = 1.;
349   }
350   else
351   {
352     do
353     {
354       do
355   {
356     cosTheta = 2.*G4UniformRand()-1.;
357     fCosTheta = (1.+cosTheta*cosTheta)/2.;
358   }
359       while (fCosTheta < G4UniformRand());
360   
361       x = xFactor*std::sqrt((1.-cosTheta)/2.);
362   
363       if (x > 1.e+005)
364   fValue = formFactorData[Z]->Value(x);
365       else
366   fValue = formFactorData[Z]->Value(0.);
367    
368       fValue /= Z;
369       fValue *= fValue;
370     }
371     while(fValue < G4UniformRand());
372   }
373 
374   return cosTheta;
375 }
376 
377 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
378 
379 G4double G4LivermorePolarizedRayleighModel::GeneratePhi(G4double cosTheta) const
380 {
381   //  d sigma
382   // --------- = alpha * ( 1 - sin^2 (theta) * cos^2 (phi) )
383   //   d phi
384  
385   // On the average the loop takes no more than 2 iterations before exiting 
386 
387   G4double phi;
388   G4double cosPhi;
389   G4double phiProbability;
390   G4double sin2Theta;
391  
392   sin2Theta=1.-cosTheta*cosTheta;
393  
394   do
395     {
396       phi = CLHEP::twopi * G4UniformRand();
397       cosPhi = std::cos(phi);
398       phiProbability= 1. - sin2Theta*cosPhi*cosPhi;
399     }
400   while (phiProbability < G4UniformRand());
401  
402   return phi;
403 }
404 
405 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
406 
407 G4double G4LivermorePolarizedRayleighModel::GeneratePolarizationAngle(void) const
408 {
409   // Rayleigh polarization is always on the x' direction
410   return 0;
411 }
412 
413 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
414 
415 G4ThreeVector G4LivermorePolarizedRayleighModel::GetPhotonPolarization(const G4DynamicParticle&  photon)
416 {
417   // From G4VLowEnergyDiscretePhotonProcess.cc
418   G4ThreeVector photonMomentumDirection;
419   G4ThreeVector photonPolarization;
420 
421   photonPolarization = photon.GetPolarization(); 
422   photonMomentumDirection = photon.GetMomentumDirection();
423 
424   if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
425     {
426       // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
427       // then polarization is choosen randomly. 
428       G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
429       G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
430   
431       G4double angle(G4UniformRand() * CLHEP::twopi);
432   
433       e1*=std::cos(angle);
434       e2*=std::sin(angle);
435   
436       photonPolarization=e1+e2;
437     }
438   else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
439     {
440       // if |photonPolarization * photonDirection0| != 0.
441       // then polarization is made orthonormal;  
442       photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
443     }
444  
445   return photonPolarization.unit();
446 }
447 
448 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
449  
450 void G4LivermorePolarizedRayleighModel::InitialiseForElement(
451                 const G4ParticleDefinition*, G4int Z)
452 {
453   G4AutoLock l(&LivermorePolarizedRayleighModelMutex);
454   if(nullptr == dataCS[Z]) { ReadData(Z); }
455   l.unlock();
456 }
457