Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/eRosita/physics/src/G4LowEnergyRayleigh.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 ]

Diff markup

Differences between /examples/advanced/eRosita/physics/src/G4LowEnergyRayleigh.cc (Version 11.3.0) and /examples/advanced/eRosita/physics/src/G4LowEnergyRayleigh.cc (Version 10.4)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
 26 // -------------------------------------------     26 // --------------------------------------------------------------------
 27 //                                                 27 //
                                                   >>  28 // $Id: G4LowEnergyRayleigh.cc 107396 2017-11-10 08:28:08Z gcosmo $
                                                   >>  29 // GEANT4 tag $Name: geant4-09-01-ref-00 $
 28 //                                                 30 //
 29 // Author: A. Forti                                31 // Author: A. Forti
 30 //         Maria Grazia Pia (Maria.Grazia.Pia@     32 //         Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
 31 //                                                 33 //
 32 // History:                                        34 // History:
 33 // --------                                        35 // --------
 34 // Added Livermore data table construction met     36 // Added Livermore data table construction methods A. Forti
 35 // Added BuildMeanFreePath A. Forti                37 // Added BuildMeanFreePath A. Forti
 36 // Added PostStepDoIt A. Forti                     38 // Added PostStepDoIt A. Forti
 37 // Added SelectRandomAtom A. Forti                 39 // Added SelectRandomAtom A. Forti
 38 // Added map of the elements  A.Forti              40 // Added map of the elements  A.Forti
 39 // 24.04.01 V.Ivanchenko remove RogueWave          41 // 24.04.01 V.Ivanchenko remove RogueWave
 40 // 11.08.2001 MGP - Major revision according t     42 // 11.08.2001 MGP - Major revision according to a design iteration
 41 // 06.10.2001 MGP - Added strategy to test ran     43 // 06.10.2001 MGP - Added strategy to test range for secondary generation
 42 // 05.06.2002 F.Longo and G.Depaola  - bug fix     44 // 05.06.2002 F.Longo and G.Depaola  - bug fixed in angular distribution
 43 // 20.10.2002 G. Depaola   - Change sampling m     45 // 20.10.2002 G. Depaola   - Change sampling method of theta
 44 // 22.01.2003 V.Ivanchenko - Cut per region        46 // 22.01.2003 V.Ivanchenko - Cut per region
 45 // 24.04.2003 V.Ivanchenko - Cut per region mf     47 // 24.04.2003 V.Ivanchenko - Cut per region mfpt
 46 //                                                 48 //
 47 // -------------------------------------------     49 // --------------------------------------------------------------------
 48                                                    50 
 49 #include "G4LowEnergyRayleigh.hh"                  51 #include "G4LowEnergyRayleigh.hh"
 50 #include "Randomize.hh"                            52 #include "Randomize.hh"
 51 #include "G4PhysicalConstants.hh"                  53 #include "G4PhysicalConstants.hh"
 52 #include "G4SystemOfUnits.hh"                      54 #include "G4SystemOfUnits.hh"
 53 #include "G4ParticleDefinition.hh"                 55 #include "G4ParticleDefinition.hh"
 54 #include "G4Track.hh"                              56 #include "G4Track.hh"
 55 #include "G4Step.hh"                               57 #include "G4Step.hh"
 56 #include "G4ForceCondition.hh"                     58 #include "G4ForceCondition.hh"
 57 #include "G4Gamma.hh"                              59 #include "G4Gamma.hh"
 58 #include "G4Electron.hh"                           60 #include "G4Electron.hh"
 59 #include "G4DynamicParticle.hh"                    61 #include "G4DynamicParticle.hh"
 60 #include "G4VParticleChange.hh"                    62 #include "G4VParticleChange.hh"
 61 #include "G4ThreeVector.hh"                        63 #include "G4ThreeVector.hh"
 62 #include "G4RDVCrossSectionHandler.hh"             64 #include "G4RDVCrossSectionHandler.hh"
 63 #include "G4RDCrossSectionHandler.hh"              65 #include "G4RDCrossSectionHandler.hh"
 64 #include "G4RDVEMDataSet.hh"                       66 #include "G4RDVEMDataSet.hh"
 65 #include "G4RDCompositeEMDataSet.hh"               67 #include "G4RDCompositeEMDataSet.hh"
 66 #include "G4RDVDataSetAlgorithm.hh"                68 #include "G4RDVDataSetAlgorithm.hh"
 67 #include "G4RDLogLogInterpolation.hh"              69 #include "G4RDLogLogInterpolation.hh"
 68                                                    70 
 69 #include "G4MaterialCutsCouple.hh"                 71 #include "G4MaterialCutsCouple.hh"
 70                                                    72 
 71 G4LowEnergyRayleigh::G4LowEnergyRayleigh(const     73 G4LowEnergyRayleigh::G4LowEnergyRayleigh(const G4String& processName)
 72   : G4VDiscreteProcess(processName),               74   : G4VDiscreteProcess(processName),
 73     lowEnergyLimit(250*eV),                        75     lowEnergyLimit(250*eV),
 74     highEnergyLimit(100*GeV),                      76     highEnergyLimit(100*GeV),
 75     intrinsicLowEnergyLimit(10*eV),                77     intrinsicLowEnergyLimit(10*eV),
 76     intrinsicHighEnergyLimit(100*GeV)              78     intrinsicHighEnergyLimit(100*GeV)
 77 {                                                  79 {
 78   if (lowEnergyLimit < intrinsicLowEnergyLimit     80   if (lowEnergyLimit < intrinsicLowEnergyLimit ||
 79       highEnergyLimit > intrinsicHighEnergyLim     81       highEnergyLimit > intrinsicHighEnergyLimit)
 80     {                                              82     {
 81       G4Exception("G4LowEnergyRayleigh::G4LowE     83       G4Exception("G4LowEnergyRayleigh::G4LowEnergyRayleigh()",
 82                   "OutOfRange", FatalException     84                   "OutOfRange", FatalException,
 83                   "Energy limit outside intrin     85                   "Energy limit outside intrinsic process validity range!");
 84     }                                              86     }
 85                                                    87 
 86   crossSectionHandler = new G4RDCrossSectionHa     88   crossSectionHandler = new G4RDCrossSectionHandler();
 87                                                    89 
 88   G4RDVDataSetAlgorithm* ffInterpolation = new     90   G4RDVDataSetAlgorithm* ffInterpolation = new G4RDLogLogInterpolation;
 89   G4String formFactorFile = "rayl/re-ff-";         91   G4String formFactorFile = "rayl/re-ff-";
 90   formFactorData = new G4RDCompositeEMDataSet(     92   formFactorData = new G4RDCompositeEMDataSet(ffInterpolation,1.,1.);
 91   formFactorData->LoadData(formFactorFile);        93   formFactorData->LoadData(formFactorFile);
 92                                                    94 
 93   meanFreePathTable = 0;                           95   meanFreePathTable = 0;
 94                                                    96 
 95    if (verboseLevel > 0)                           97    if (verboseLevel > 0)
 96      {                                             98      {
 97        G4cout << GetProcessName() << " is crea     99        G4cout << GetProcessName() << " is created " << G4endl
 98         << "Energy range: "                       100         << "Energy range: "
 99         << lowEnergyLimit / keV << " keV - "      101         << lowEnergyLimit / keV << " keV - "
100         << highEnergyLimit / GeV << " GeV"        102         << highEnergyLimit / GeV << " GeV"
101         << G4endl;                                103         << G4endl;
102      }                                            104      }
103 }                                                 105 }
104                                                   106 
105 G4LowEnergyRayleigh::~G4LowEnergyRayleigh()       107 G4LowEnergyRayleigh::~G4LowEnergyRayleigh()
106 {                                                 108 {
107   delete meanFreePathTable;                       109   delete meanFreePathTable;
108   delete crossSectionHandler;                     110   delete crossSectionHandler;
109   delete formFactorData;                          111   delete formFactorData;
110 }                                                 112 }
111                                                   113 
112 void G4LowEnergyRayleigh::BuildPhysicsTable(co    114 void G4LowEnergyRayleigh::BuildPhysicsTable(const G4ParticleDefinition& )
113 {                                                 115 {
114                                                   116 
115   crossSectionHandler->Clear();                   117   crossSectionHandler->Clear();
116   G4String crossSectionFile = "rayl/re-cs-";      118   G4String crossSectionFile = "rayl/re-cs-";
117   crossSectionHandler->LoadData(crossSectionFi    119   crossSectionHandler->LoadData(crossSectionFile);
118                                                   120 
119   delete meanFreePathTable;                       121   delete meanFreePathTable;
120   meanFreePathTable = crossSectionHandler->Bui    122   meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
121 }                                                 123 }
122                                                   124 
123 G4VParticleChange* G4LowEnergyRayleigh::PostSt    125 G4VParticleChange* G4LowEnergyRayleigh::PostStepDoIt(const G4Track& aTrack,
124                  const G4Step& aStep)             126                  const G4Step& aStep)
125 {                                                 127 {
126                                                   128 
127   aParticleChange.Initialize(aTrack);             129   aParticleChange.Initialize(aTrack);
128                                                   130 
129   const G4DynamicParticle* incidentPhoton = aT    131   const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
130   G4double photonEnergy0 = incidentPhoton->Get    132   G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
131                                                   133 
132   if (photonEnergy0 <= lowEnergyLimit)            134   if (photonEnergy0 <= lowEnergyLimit)
133     {                                             135     {
134       aParticleChange.ProposeTrackStatus(fStop    136       aParticleChange.ProposeTrackStatus(fStopAndKill);
135       aParticleChange.ProposeEnergy(0.);          137       aParticleChange.ProposeEnergy(0.);
136       aParticleChange.ProposeLocalEnergyDeposi    138       aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
137       return G4VDiscreteProcess::PostStepDoIt(    139       return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
138     }                                             140     }
139                                                   141 
140   //  G4double e0m = photonEnergy0 / electron_    142   //  G4double e0m = photonEnergy0 / electron_mass_c2 ;
141   G4ParticleMomentum photonDirection0 = incide    143   G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
142                                                   144 
143   // Select randomly one element in the curren    145   // Select randomly one element in the current material
144   const G4MaterialCutsCouple* couple = aTrack.    146   const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
145   G4int Z = crossSectionHandler->SelectRandomA    147   G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
146                                                   148 
147   // Sample the angle of the scattered photon     149   // Sample the angle of the scattered photon
148                                                   150 
149   G4double wlPhoton = h_Planck*c_light/photonE    151   G4double wlPhoton = h_Planck*c_light/photonEnergy0;
150                                                   152 
151   G4double gReject,x,dataFormFactor;              153   G4double gReject,x,dataFormFactor;
152   G4double randomFormFactor;                      154   G4double randomFormFactor;
153   G4double cosTheta;                              155   G4double cosTheta;
154   G4double sinTheta;                              156   G4double sinTheta;
155   G4double fcostheta;                             157   G4double fcostheta;
156                                                   158 
157   do                                              159   do
158     {                                             160     {
159       do                                          161       do
160       {                                           162       {
161       cosTheta = 2. * G4UniformRand() - 1.;       163       cosTheta = 2. * G4UniformRand() - 1.;
162       fcostheta = ( 1. + cosTheta*cosTheta)/2.    164       fcostheta = ( 1. + cosTheta*cosTheta)/2.;
163       } while (fcostheta < G4UniformRand());      165       } while (fcostheta < G4UniformRand());
164                                                   166 
165       G4double sinThetaHalf = std::sqrt((1. -     167       G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.);
166       x = sinThetaHalf / (wlPhoton/cm);           168       x = sinThetaHalf / (wlPhoton/cm);
167       if (x > 1.e+005)                            169       if (x > 1.e+005)
168          dataFormFactor = formFactorData->Find    170          dataFormFactor = formFactorData->FindValue(x,Z-1);
169       else                                        171       else
170          dataFormFactor = formFactorData->Find    172          dataFormFactor = formFactorData->FindValue(0.,Z-1);
171       randomFormFactor = G4UniformRand() * Z *    173       randomFormFactor = G4UniformRand() * Z * Z;
172       sinTheta = std::sqrt(1. - cosTheta*cosTh    174       sinTheta = std::sqrt(1. - cosTheta*cosTheta);
173       gReject = dataFormFactor * dataFormFacto    175       gReject = dataFormFactor * dataFormFactor;
174                                                   176 
175     } while( gReject < randomFormFactor);         177     } while( gReject < randomFormFactor);
176                                                   178 
177   // Scattered photon angles. ( Z - axis along    179   // Scattered photon angles. ( Z - axis along the parent photon)
178   G4double phi = twopi * G4UniformRand() ;        180   G4double phi = twopi * G4UniformRand() ;
179   G4double dirX = sinTheta*std::cos(phi);         181   G4double dirX = sinTheta*std::cos(phi);
180   G4double dirY = sinTheta*std::sin(phi);         182   G4double dirY = sinTheta*std::sin(phi);
181   G4double dirZ = cosTheta;                       183   G4double dirZ = cosTheta;
182                                                   184 
183   // Update G4VParticleChange for the scattere    185   // Update G4VParticleChange for the scattered photon
184   G4ThreeVector photonDirection1(dirX, dirY, d    186   G4ThreeVector photonDirection1(dirX, dirY, dirZ);
185                                                   187 
186   photonDirection1.rotateUz(photonDirection0);    188   photonDirection1.rotateUz(photonDirection0);
187   aParticleChange.ProposeEnergy(photonEnergy0)    189   aParticleChange.ProposeEnergy(photonEnergy0);
188   aParticleChange.ProposeMomentumDirection(pho    190   aParticleChange.ProposeMomentumDirection(photonDirection1);
189                                                   191 
190   aParticleChange.SetNumberOfSecondaries(0);      192   aParticleChange.SetNumberOfSecondaries(0);
191                                                   193 
192   return G4VDiscreteProcess::PostStepDoIt(aTra    194   return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
193 }                                                 195 }
194                                                   196 
195 G4bool G4LowEnergyRayleigh::IsApplicable(const    197 G4bool G4LowEnergyRayleigh::IsApplicable(const G4ParticleDefinition& particle)
196 {                                                 198 {
197   return ( &particle == G4Gamma::Gamma() );       199   return ( &particle == G4Gamma::Gamma() ); 
198 }                                                 200 }
199                                                   201 
200 G4double G4LowEnergyRayleigh::GetMeanFreePath(    202 G4double G4LowEnergyRayleigh::GetMeanFreePath(const G4Track& track, 
201                 G4double, // previousStepSize     203                 G4double, // previousStepSize
202                 G4ForceCondition*)                204                 G4ForceCondition*)
203 {                                                 205 {
204   const G4DynamicParticle* photon = track.GetD    206   const G4DynamicParticle* photon = track.GetDynamicParticle();
205   G4double energy = photon->GetKineticEnergy()    207   G4double energy = photon->GetKineticEnergy();
206   const G4MaterialCutsCouple* couple = track.G    208   const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
207   size_t materialIndex = couple->GetIndex();      209   size_t materialIndex = couple->GetIndex();
208                                                   210 
209   G4double meanFreePath;                          211   G4double meanFreePath;
210   if (energy > highEnergyLimit) meanFreePath =    212   if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
211   else if (energy < lowEnergyLimit) meanFreePa    213   else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
212   else meanFreePath = meanFreePathTable->FindV    214   else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
213   return meanFreePath;                            215   return meanFreePath;
214 }                                                 216 }
215                                                   217