processes/electromagnetic/lowenergy/src/G4LowEnergyRayleigh.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 // --------------------------------------------------------------------
00027 //
00028 // $Id: G4LowEnergyRayleigh.cc,v 1.37 2006/06/29 19:40:29 gunter Exp $
00029 // GEANT4 tag $Name: geant4-09-01-patch-01 $
00030 //
00031 // Author: A. Forti
00032 //         Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
00033 //
00034 // History:
00035 // --------
00036 // Added Livermore data table construction methods A. Forti
00037 // Added BuildMeanFreePath A. Forti
00038 // Added PostStepDoIt A. Forti
00039 // Added SelectRandomAtom A. Forti
00040 // Added map of the elements  A.Forti
00041 // 24.04.01 V.Ivanchenko remove RogueWave
00042 // 11.08.2001 MGP - Major revision according to a design iteration
00043 // 06.10.2001 MGP - Added strategy to test range for secondary generation
00044 // 05.06.2002 F.Longo and G.Depaola  - bug fixed in angular distribution
00045 // 20.10.2002 G. Depaola   - Change sampling method of theta
00046 // 22.01.2003 V.Ivanchenko - Cut per region
00047 // 24.04.2003 V.Ivanchenko - Cut per region mfpt
00048 //
00049 // --------------------------------------------------------------------
00050 
00051 #include "G4LowEnergyRayleigh.hh"
00052 #include "Randomize.hh"
00053 #include "G4ParticleDefinition.hh"
00054 #include "G4Track.hh"
00055 #include "G4Step.hh"
00056 #include "G4ForceCondition.hh"
00057 #include "G4Gamma.hh"
00058 #include "G4Electron.hh"
00059 #include "G4DynamicParticle.hh"
00060 #include "G4VParticleChange.hh"
00061 #include "G4ThreeVector.hh"
00062 #include "G4VCrossSectionHandler.hh"
00063 #include "G4CrossSectionHandler.hh"
00064 #include "G4VEMDataSet.hh"
00065 #include "G4CompositeEMDataSet.hh"
00066 #include "G4VDataSetAlgorithm.hh"
00067 #include "G4LogLogInterpolation.hh"
00068 
00069 #include "G4MaterialCutsCouple.hh"
00070 
00071 G4LowEnergyRayleigh::G4LowEnergyRayleigh(const G4String& processName)
00072   : G4VDiscreteProcess(processName),
00073     lowEnergyLimit(250*eV),
00074     highEnergyLimit(100*GeV),
00075     intrinsicLowEnergyLimit(10*eV),
00076     intrinsicHighEnergyLimit(100*GeV)
00077 {
00078   if (lowEnergyLimit < intrinsicLowEnergyLimit ||
00079       highEnergyLimit > intrinsicHighEnergyLimit)
00080     {
00081       G4Exception("G4LowEnergyRayleigh::G4LowEnergyRayleigh - energy limit outside intrinsic process validity range");
00082     }
00083 
00084   crossSectionHandler = new G4CrossSectionHandler();
00085 
00086   G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
00087   G4String formFactorFile = "rayl/re-ff-";
00088   formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
00089   formFactorData->LoadData(formFactorFile);
00090 
00091   meanFreePathTable = 0;
00092 
00093    if (verboseLevel > 0)
00094      {
00095        G4cout << GetProcessName() << " is created " << G4endl
00096               << "Energy range: "
00097               << lowEnergyLimit / keV << " keV - "
00098               << highEnergyLimit / GeV << " GeV"
00099               << G4endl;
00100      }
00101 }
00102 
00103 G4LowEnergyRayleigh::~G4LowEnergyRayleigh()
00104 {
00105   delete meanFreePathTable;
00106   delete crossSectionHandler;
00107   delete formFactorData;
00108 }
00109 
00110 void G4LowEnergyRayleigh::BuildPhysicsTable(const G4ParticleDefinition& )
00111 {
00112 
00113   crossSectionHandler->Clear();
00114   G4String crossSectionFile = "rayl/re-cs-";
00115   crossSectionHandler->LoadData(crossSectionFile);
00116 
00117   delete meanFreePathTable;
00118   meanFreePathTable = crossSectionHandler->BuildMeanFreePathForMaterials();
00119 }
00120 
00121 G4VParticleChange* G4LowEnergyRayleigh::PostStepDoIt(const G4Track& aTrack,
00122                                                      const G4Step& aStep)
00123 {
00124 
00125   aParticleChange.Initialize(aTrack);
00126 
00127   const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
00128   G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
00129 
00130   if (photonEnergy0 <= lowEnergyLimit)
00131     {
00132       aParticleChange.ProposeTrackStatus(fStopAndKill);
00133       aParticleChange.ProposeEnergy(0.);
00134       aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
00135       return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
00136     }
00137 
00138   //  G4double e0m = photonEnergy0 / electron_mass_c2 ;
00139   G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
00140 
00141   // Select randomly one element in the current material
00142   const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
00143   G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
00144 
00145   // Sample the angle of the scattered photon
00146 
00147   G4double wlPhoton = h_Planck*c_light/photonEnergy0;
00148 
00149   G4double gReject,x,dataFormFactor;
00150   G4double randomFormFactor;
00151   G4double cosTheta;
00152   G4double sinTheta;
00153   G4double fcostheta;
00154 
00155   do
00156     {
00157       do
00158       {
00159       cosTheta = 2. * G4UniformRand() - 1.;
00160       fcostheta = ( 1. + cosTheta*cosTheta)/2.;
00161       } while (fcostheta < G4UniformRand());
00162 
00163       G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.);
00164       x = sinThetaHalf / (wlPhoton/cm);
00165       if (x > 1.e+005)
00166          dataFormFactor = formFactorData->FindValue(x,Z-1);
00167       else
00168          dataFormFactor = formFactorData->FindValue(0.,Z-1);
00169       randomFormFactor = G4UniformRand() * Z * Z;
00170       sinTheta = std::sqrt(1. - cosTheta*cosTheta);
00171       gReject = dataFormFactor * dataFormFactor;
00172 
00173     } while( gReject < randomFormFactor);
00174 
00175   // Scattered photon angles. ( Z - axis along the parent photon)
00176   G4double phi = twopi * G4UniformRand() ;
00177   G4double dirX = sinTheta*std::cos(phi);
00178   G4double dirY = sinTheta*std::sin(phi);
00179   G4double dirZ = cosTheta;
00180 
00181   // Update G4VParticleChange for the scattered photon
00182   G4ThreeVector photonDirection1(dirX, dirY, dirZ);
00183 
00184   photonDirection1.rotateUz(photonDirection0);
00185   aParticleChange.ProposeEnergy(photonEnergy0);
00186   aParticleChange.ProposeMomentumDirection(photonDirection1);
00187 
00188   aParticleChange.SetNumberOfSecondaries(0);
00189 
00190   return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
00191 }
00192 
00193 G4bool G4LowEnergyRayleigh::IsApplicable(const G4ParticleDefinition& particle)
00194 {
00195   return ( &particle == G4Gamma::Gamma() ); 
00196 }
00197 
00198 G4double G4LowEnergyRayleigh::GetMeanFreePath(const G4Track& track, 
00199                                               G4double, // previousStepSize
00200                                               G4ForceCondition*)
00201 {
00202   const G4DynamicParticle* photon = track.GetDynamicParticle();
00203   G4double energy = photon->GetKineticEnergy();
00204   const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
00205   size_t materialIndex = couple->GetIndex();
00206 
00207   G4double meanFreePath;
00208   if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
00209   else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
00210   else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
00211   return meanFreePath;
00212 }

Generated on Fri Apr 11 17:10:08 2008 for Geant4 by  doxygen 1.4.7