G4TransitionRadiation.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // G4TransitionRadiation class -- implementation file
00029 
00030 // GEANT 4 class implementation file --- Copyright CERN 1995
00031 // CERN Geneva Switzerland
00032 
00033 // For information related to this code, please, contact
00034 // CERN, CN Division, ASD Group
00035 // History:
00036 // 1st version 11.09.97 V. Grichine (Vladimir.Grichine@cern.ch )
00037 // 2nd version 16.12.97 V. Grichine
00038 // 3rd version 28.07.05, P.Gumplinger add G4ProcessType to constructor
00039 
00040 
00041 #include <cmath>
00042 
00043 #include "G4TransitionRadiation.hh"
00044 #include "G4Material.hh"
00045 #include "G4EmProcessSubType.hh"
00046 
00047 // Local constants
00048 
00049 const G4int   G4TransitionRadiation::fSympsonNumber = 100 ;
00050 const G4int   G4TransitionRadiation::fGammaNumber = 15 ;
00051 const G4int   G4TransitionRadiation::fPointNumber = 100 ;
00052 
00053 
00055 //
00056 // Constructor for selected couple of materials
00057 //
00058 
00059 G4TransitionRadiation::
00060 G4TransitionRadiation( const G4String& processName, G4ProcessType type )
00061   : G4VDiscreteProcess(processName, type)
00062 {
00063   SetProcessSubType(fTransitionRadiation);
00064   fMatIndex1 = fMatIndex2 = 0;
00065 
00066   fGamma = fEnergy = fVarAngle = fMinEnergy = fMaxEnergy = fMaxTheta = fSigma1 = fSigma2 = 0.0;
00067 }
00068 
00070 //
00071 // Destructor
00072 //
00073 
00074 G4TransitionRadiation::~G4TransitionRadiation()
00075 {}
00076 
00077 G4bool 
00078 G4TransitionRadiation::IsApplicable(const G4ParticleDefinition& aParticleType)
00079 {
00080   return ( aParticleType.GetPDGCharge() != 0.0 );
00081 }
00082 
00083 G4double G4TransitionRadiation::GetMeanFreePath(const G4Track&,
00084                                                 G4double,
00085                                                 G4ForceCondition* condition)
00086 {
00087   *condition = Forced;
00088   return DBL_MAX;      // so TR doesn't limit mean free path
00089 }
00090 
00091 G4VParticleChange* G4TransitionRadiation::PostStepDoIt(const G4Track&,
00092                                                        const G4Step&)
00093 {
00094   ClearNumberOfInteractionLengthLeft();
00095   return &aParticleChange;
00096 }
00097 
00099 //
00100 // Sympson integral of TR spectral-angle density over energy between
00101 // the limits energy 1 and energy2 at fixed varAngle = 1 - std::cos(Theta)
00102 
00103 G4double
00104 G4TransitionRadiation::IntegralOverEnergy( G4double energy1,
00105                                            G4double energy2,
00106                                            G4double varAngle     )  const
00107 {
00108   G4int i ;
00109   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00110   h = 0.5*(energy2 - energy1)/fSympsonNumber ;
00111   for(i=1;i<fSympsonNumber;i++)
00112   {
00113     sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle)  ;
00114     sumOdd  += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
00115   }
00116   sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
00117   return h*(  SpectralAngleTRdensity(energy1,varAngle)
00118             + SpectralAngleTRdensity(energy2,varAngle)
00119             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00120 }
00121 
00122 
00123 
00125 //
00126 // Sympson integral of TR spectral-angle density over energy between
00127 // the limits varAngle1 and varAngle2 at fixed energy
00128 
00129 G4double
00130 G4TransitionRadiation::IntegralOverAngle( G4double energy,
00131                                           G4double varAngle1,
00132                                           G4double varAngle2     ) const
00133 {
00134   G4int i ;
00135   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00136   h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
00137   for(i=1;i<fSympsonNumber;i++)
00138   {
00139     sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h)  ;
00140     sumOdd  += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
00141   }
00142   sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
00143 
00144   return h*(  SpectralAngleTRdensity(energy,varAngle1)
00145             + SpectralAngleTRdensity(energy,varAngle2)
00146             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00147 }
00148 
00150 //
00151 // The number of transition radiation photons generated in the
00152 // angle interval between varAngle1 and varAngle2
00153 //
00154 
00155 G4double G4TransitionRadiation::
00156 AngleIntegralDistribution( G4double varAngle1,
00157                            G4double varAngle2     )   const
00158 {
00159   G4int i ;
00160   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00161   h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
00162   for(i=1;i<fSympsonNumber;i++)
00163   {
00164    sumEven += IntegralOverEnergy(fMinEnergy,
00165                                  fMinEnergy +0.3*(fMaxEnergy-fMinEnergy),
00166                                  varAngle1 + 2*i*h)
00167             + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00168                                  fMaxEnergy,
00169                                  varAngle1 + 2*i*h);
00170    sumOdd  += IntegralOverEnergy(fMinEnergy,
00171                                  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00172                                  varAngle1 + (2*i - 1)*h)
00173             + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00174                                  fMaxEnergy,
00175                                  varAngle1 + (2*i - 1)*h) ;
00176   }
00177   sumOdd += IntegralOverEnergy(fMinEnergy,
00178                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00179                                varAngle1 + (2*fSympsonNumber - 1)*h)
00180           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00181                                fMaxEnergy,
00182                                varAngle1 + (2*fSympsonNumber - 1)*h) ;
00183 
00184   return h*(IntegralOverEnergy(fMinEnergy,
00185                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00186                                varAngle1)
00187           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00188                                fMaxEnergy,
00189                                varAngle1)
00190           + IntegralOverEnergy(fMinEnergy,
00191                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00192                                varAngle2)
00193           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00194                                fMaxEnergy,
00195                                varAngle2)
00196             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00197 }
00198 
00200 //
00201 // The number of transition radiation photons, generated in the
00202 // energy interval between energy1 and energy2
00203 //
00204 
00205 G4double G4TransitionRadiation::
00206 EnergyIntegralDistribution( G4double energy1,
00207                             G4double energy2     )  const
00208 {
00209   G4int i ;
00210   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00211   h = 0.5*(energy2 - energy1)/fSympsonNumber ;
00212   for(i=1;i<fSympsonNumber;i++)
00213   {
00214    sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta )
00215             + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta);
00216    sumOdd  += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta)
00217             + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ;
00218   }
00219   sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
00220                               0.0,0.01*fMaxTheta)
00221           + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
00222                               0.01*fMaxTheta,fMaxTheta) ;
00223 
00224   return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta)
00225           + IntegralOverAngle(energy1,0.01*fMaxTheta,fMaxTheta)
00226           + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta)
00227           + IntegralOverAngle(energy2,0.01*fMaxTheta,fMaxTheta)
00228             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00229 }
00230 
00231 
00232 
00233 
00234 // end of G4TransitionRadiation implementation file --------------------------

Generated on Mon May 27 17:50:02 2013 for Geant4 by  doxygen 1.4.7