Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/xrays/src/G4TransitionRadiation.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 // G4TransitionRadiation class -- implementation file
 27 
 28 // GEANT 4 class implementation file --- Copyright CERN 1995
 29 
 30 // For information related to this code, please, contact
 31 // CERN, CN Division, ASD Group
 32 // History:
 33 // 1st version 11.09.97 V. Grichine (Vladimir.Grichine@cern.ch )
 34 // 2nd version 16.12.97 V. Grichine
 35 // 3rd version 28.07.05, P.Gumplinger add G4ProcessType to constructor
 36 
 37 //#include <cmath>
 38 
 39 #include "G4TransitionRadiation.hh"
 40 
 41 #include "G4EmProcessSubType.hh"
 42 
 43 ///////////////////////////////////////////////////////////////////////
 44 // Constructor for selected couple of materials
 45 G4TransitionRadiation::G4TransitionRadiation(const G4String& processName,
 46                                              G4ProcessType type)
 47   : G4VDiscreteProcess(processName, type)
 48 {
 49   SetProcessSubType(fTransitionRadiation);
 50   fMatIndex1 = fMatIndex2 = 0;
 51 
 52   fGamma = fEnergy = fVarAngle = fMinEnergy = fMaxEnergy = fMaxTheta = 0.0;
 53   fSigma1 = fSigma2 = 0.0;
 54 }
 55 
 56 //////////////////////////////////////////////////////////////////////
 57 // Destructor
 58 G4TransitionRadiation::~G4TransitionRadiation() = default;
 59 
 60 void G4TransitionRadiation::ProcessDescription(std::ostream& out) const
 61 {
 62   out << "Base class for simulation of x-ray transition radiation.\n";
 63 }
 64 
 65 G4bool G4TransitionRadiation::IsApplicable(
 66   const G4ParticleDefinition& aParticleType)
 67 {
 68   return (aParticleType.GetPDGCharge() != 0.0);
 69 }
 70 
 71 G4double G4TransitionRadiation::GetMeanFreePath(const G4Track&, G4double,
 72                                                 G4ForceCondition* condition)
 73 {
 74   *condition = Forced;
 75   return DBL_MAX;  // so TR doesn't limit mean free path
 76 }
 77 
 78 G4VParticleChange* G4TransitionRadiation::PostStepDoIt(const G4Track&,
 79                                                        const G4Step&)
 80 {
 81   ClearNumberOfInteractionLengthLeft();
 82   return &aParticleChange;
 83 }
 84 
 85 ///////////////////////////////////////////////////////////////////
 86 // Sympson integral of TR spectral-angle density over energy between
 87 // the limits energy 1 and energy2 at fixed varAngle = 1 - std::cos(Theta)
 88 G4double G4TransitionRadiation::IntegralOverEnergy(G4double energy1,
 89                                                    G4double energy2,
 90                                                    G4double varAngle) const
 91 {
 92   G4int i;
 93   G4double h, sumEven = 0.0, sumOdd = 0.0;
 94   h = 0.5 * (energy2 - energy1) / fSympsonNumber;
 95   for(i = 1; i < fSympsonNumber; i++)
 96   {
 97     sumEven += SpectralAngleTRdensity(energy1 + 2 * i * h, varAngle);
 98     sumOdd += SpectralAngleTRdensity(energy1 + (2 * i - 1) * h, varAngle);
 99   }
100   sumOdd +=
101     SpectralAngleTRdensity(energy1 + (2 * fSympsonNumber - 1) * h, varAngle);
102   return h *
103          (SpectralAngleTRdensity(energy1, varAngle) +
104           SpectralAngleTRdensity(energy2, varAngle) + 4.0 * sumOdd +
105           2.0 * sumEven) /
106          3.0;
107 }
108 
109 ///////////////////////////////////////////////////////////////////
110 // Sympson integral of TR spectral-angle density over energy between
111 // the limits varAngle1 and varAngle2 at fixed energy
112 G4double G4TransitionRadiation::IntegralOverAngle(G4double energy,
113                                                   G4double varAngle1,
114                                                   G4double varAngle2) const
115 {
116   G4int i;
117   G4double h, sumEven = 0.0, sumOdd = 0.0;
118   h = 0.5 * (varAngle2 - varAngle1) / fSympsonNumber;
119   for(i = 1; i < fSympsonNumber; ++i)
120   {
121     sumEven += SpectralAngleTRdensity(energy, varAngle1 + 2 * i * h);
122     sumOdd += SpectralAngleTRdensity(energy, varAngle1 + (2 * i - 1) * h);
123   }
124   sumOdd +=
125     SpectralAngleTRdensity(energy, varAngle1 + (2 * fSympsonNumber - 1) * h);
126 
127   return h *
128          (SpectralAngleTRdensity(energy, varAngle1) +
129           SpectralAngleTRdensity(energy, varAngle2) + 4.0 * sumOdd +
130           2.0 * sumEven) /
131          3.0;
132 }
133 
134 ///////////////////////////////////////////////////////////////////
135 // The number of transition radiation photons generated in the
136 // angle interval between varAngle1 and varAngle2
137 G4double G4TransitionRadiation::AngleIntegralDistribution(
138   G4double varAngle1, G4double varAngle2) const
139 {
140   G4int i;
141   G4double h, sumEven = 0.0, sumOdd = 0.0;
142   h = 0.5 * (varAngle2 - varAngle1) / fSympsonNumber;
143   for(i = 1; i < fSympsonNumber; ++i)
144   {
145     sumEven += IntegralOverEnergy(fMinEnergy,
146                                   fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
147                                   varAngle1 + 2 * i * h) +
148                IntegralOverEnergy(fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
149                                   fMaxEnergy, varAngle1 + 2 * i * h);
150     sumOdd += IntegralOverEnergy(fMinEnergy,
151                                  fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
152                                  varAngle1 + (2 * i - 1) * h) +
153               IntegralOverEnergy(fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
154                                  fMaxEnergy, varAngle1 + (2 * i - 1) * h);
155   }
156   sumOdd +=
157     IntegralOverEnergy(fMinEnergy, fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
158                        varAngle1 + (2 * fSympsonNumber - 1) * h) +
159     IntegralOverEnergy(fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy), fMaxEnergy,
160                        varAngle1 + (2 * fSympsonNumber - 1) * h);
161 
162   return h *
163          (IntegralOverEnergy(fMinEnergy,
164                              fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
165                              varAngle1) +
166           IntegralOverEnergy(fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
167                              fMaxEnergy, varAngle1) +
168           IntegralOverEnergy(fMinEnergy,
169                              fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
170                              varAngle2) +
171           IntegralOverEnergy(fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
172                              fMaxEnergy, varAngle2) +
173           4.0 * sumOdd + 2.0 * sumEven) /
174          3.0;
175 }
176 
177 ///////////////////////////////////////////////////////////////////
178 // The number of transition radiation photons, generated in the
179 // energy interval between energy1 and energy2
180 G4double G4TransitionRadiation::EnergyIntegralDistribution(
181   G4double energy1, G4double energy2) const
182 {
183   G4int i;
184   G4double h, sumEven = 0.0, sumOdd = 0.0;
185   h = 0.5 * (energy2 - energy1) / fSympsonNumber;
186   for(i = 1; i < fSympsonNumber; ++i)
187   {
188     sumEven +=
189       IntegralOverAngle(energy1 + 2 * i * h, 0.0, 0.01 * fMaxTheta) +
190       IntegralOverAngle(energy1 + 2 * i * h, 0.01 * fMaxTheta, fMaxTheta);
191     sumOdd +=
192       IntegralOverAngle(energy1 + (2 * i - 1) * h, 0.0, 0.01 * fMaxTheta) +
193       IntegralOverAngle(energy1 + (2 * i - 1) * h, 0.01 * fMaxTheta, fMaxTheta);
194   }
195   sumOdd += IntegralOverAngle(energy1 + (2 * fSympsonNumber - 1) * h, 0.0,
196                               0.01 * fMaxTheta) +
197             IntegralOverAngle(energy1 + (2 * fSympsonNumber - 1) * h,
198                               0.01 * fMaxTheta, fMaxTheta);
199 
200   return h *
201          (IntegralOverAngle(energy1, 0.0, 0.01 * fMaxTheta) +
202           IntegralOverAngle(energy1, 0.01 * fMaxTheta, fMaxTheta) +
203           IntegralOverAngle(energy2, 0.0, 0.01 * fMaxTheta) +
204           IntegralOverAngle(energy2, 0.01 * fMaxTheta, fMaxTheta) +
205           4.0 * sumOdd + 2.0 * sumEven) /
206          3.0;
207 }
208