G4eeToHadronsModel.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: G4eeToHadronsModel.cc 66996 2013-01-29 14:50:52Z gcosmo $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class header file
00031 //
00032 //
00033 // File name:     G4eeToHadronsModel
00034 //
00035 // Author:        Vladimir Ivanchenko
00036 //
00037 // Creation date: 12.08.2003
00038 //
00039 // Modifications:
00040 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
00041 // 18-05-05 Use optimized interfaces (V.Ivantchenko)
00042 //
00043 //
00044 // -------------------------------------------------------------------
00045 //
00046 
00047 
00048 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00050 
00051 #include "G4eeToHadronsModel.hh"
00052 #include "Randomize.hh"
00053 #include "G4PhysicalConstants.hh"
00054 #include "G4SystemOfUnits.hh"
00055 #include "G4Electron.hh"
00056 #include "G4Gamma.hh"
00057 #include "G4Positron.hh"
00058 #include "G4PionPlus.hh"
00059 #include "Randomize.hh"
00060 #include "G4Vee2hadrons.hh"
00061 #include "G4PhysicsVector.hh"
00062 #include "G4PhysicsLogVector.hh"
00063 
00064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00065 
00066 using namespace std;
00067 
00068 G4eeToHadronsModel::G4eeToHadronsModel(G4Vee2hadrons* mod, G4int ver,
00069                                        const G4String& nam)
00070   : G4VEmModel(nam),
00071     model(mod),
00072     crossPerElectron(0),
00073     crossBornPerElectron(0),
00074     isInitialised(false),
00075     nbins(100),
00076     verbose(ver)
00077 {
00078   theGamma = G4Gamma::Gamma();
00079   highKinEnergy = HighEnergyLimit();
00080   lowKinEnergy  = LowEnergyLimit();
00081   emin = lowKinEnergy;
00082   emax = highKinEnergy;
00083   peakKinEnergy = highKinEnergy;
00084   epeak = emax;
00085 }
00086 
00087 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00088 
00089 G4eeToHadronsModel::~G4eeToHadronsModel()
00090 {
00091   delete model;
00092   delete crossPerElectron;
00093   delete crossBornPerElectron;
00094 }
00095 
00096 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00097 
00098 void G4eeToHadronsModel::Initialise(const G4ParticleDefinition*,
00099                                     const G4DataVector&)
00100 {
00101   if(isInitialised) { return; }
00102   isInitialised  = true;
00103 
00104   // Lab system
00105   highKinEnergy = HighEnergyLimit();
00106   lowKinEnergy  = LowEnergyLimit();
00107 
00108   // CM system
00109   emin  = model->LowEnergy();
00110   emax  = model->HighEnergy();
00111 
00112   G4double emin0 = 
00113     2.0*electron_mass_c2*sqrt(1.0 + 0.5*lowKinEnergy/electron_mass_c2);
00114   G4double emax0 = 
00115     2.0*electron_mass_c2*sqrt(1.0 + 0.5*highKinEnergy/electron_mass_c2);
00116 
00117   // recompute low energy
00118   if(emin0 > emax) {
00119     emin0 = emax;
00120     model->SetLowEnergy(emin0);
00121   }
00122   if(emin > emin0) {
00123     emin0 = emin;
00124     lowKinEnergy  = 0.5*emin*emin/electron_mass_c2 - 2.0*electron_mass_c2;
00125     SetLowEnergyLimit(lowKinEnergy);
00126   }
00127 
00128   // recompute high energy
00129   if(emax < emax0) {
00130     emax0 = emax;
00131     highKinEnergy = 0.5*emax*emax/electron_mass_c2 - 2.0*electron_mass_c2;
00132     SetHighEnergyLimit(highKinEnergy);
00133   }
00134 
00135   // peak energy
00136   epeak = std::min(model->PeakEnergy(), emax);
00137   peakKinEnergy  = 0.5*epeak*epeak/electron_mass_c2 - 2.0*electron_mass_c2;
00138 
00139   if(verbose>0) {
00140     G4cout << "G4eeToHadronsModel::Initialise: " << G4endl;
00141     G4cout << "LabSystem: emin(GeV)= " << lowKinEnergy/GeV
00142            << " epeak(GeV)= " << peakKinEnergy/GeV
00143            << " emax(GeV)= " << highKinEnergy/GeV
00144            << G4endl;
00145     G4cout << "SM System: emin(MeV)= " << emin/MeV
00146            << " epeak(MeV)= " << epeak/MeV
00147            << " emax(MeV)= " << emax/MeV
00148            << G4endl;
00149   }
00150 
00151   if(lowKinEnergy < peakKinEnergy) {
00152     crossBornPerElectron = model->PhysicsVector(emin, emax);
00153     crossPerElectron     = model->PhysicsVector(emin, emax);
00154     nbins = crossPerElectron->GetVectorLength();
00155     for(G4int i=0; i<nbins; i++) {
00156       G4double e  = crossPerElectron->GetLowEdgeEnergy(i);
00157       G4double cs = model->ComputeCrossSection(e);
00158       crossBornPerElectron->PutValue(i, cs);
00159     }
00160     ComputeCMCrossSectionPerElectron();
00161   }
00162   if(verbose>1) {
00163     G4cout << "G4eeToHadronsModel: Cross secsions per electron"
00164            << " nbins= " << nbins
00165            << " emin(MeV)= " << emin/MeV
00166            << " emax(MeV)= " << emax/MeV
00167            << G4endl;
00168     G4bool b;
00169     for(G4int i=0; i<nbins; i++) {
00170       G4double e  = crossPerElectron->GetLowEdgeEnergy(i);
00171       G4double s1 = crossPerElectron->GetValue(e, b);
00172       G4double s2 = crossBornPerElectron->GetValue(e, b);
00173       G4cout << "E(MeV)= " << e/MeV
00174              << "  cross(nb)= " << s1/nanobarn
00175              << "  crossBorn(nb)= " << s2/nanobarn
00176              << G4endl;
00177     }
00178   }
00179 }
00180 
00181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00182 
00183 G4double G4eeToHadronsModel::CrossSectionPerVolume(
00184                                       const G4Material* mat,
00185                                       const G4ParticleDefinition* p,
00186                                       G4double kineticEnergy,
00187                                       G4double, G4double)
00188 {
00189   return mat->GetElectronDensity()*
00190     ComputeCrossSectionPerElectron(p, kineticEnergy);
00191 }
00192 
00193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00194 
00195 G4double G4eeToHadronsModel::ComputeCrossSectionPerAtom(
00196                                       const G4ParticleDefinition* p,
00197                                       G4double kineticEnergy,
00198                                       G4double Z, G4double,
00199                                       G4double, G4double)
00200 {
00201   return Z*ComputeCrossSectionPerElectron(p, kineticEnergy);
00202 }
00203 
00204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00205 
00206 G4double G4eeToHadronsModel::ComputeCrossSectionPerElectron(
00207                                           const G4ParticleDefinition*,
00208                                                 G4double kineticEnergy,
00209                                                 G4double, G4double)
00210 {
00211   G4double cross = 0.0;
00212   if(crossPerElectron) {
00213     G4bool b;
00214     G4double e = 2.0*electron_mass_c2*
00215                  sqrt(1.0 + 0.5*kineticEnergy/electron_mass_c2);
00216     cross = crossPerElectron->GetValue(e, b);
00217   }
00218   //  G4cout << "e= " << kineticEnergy << " cross= " << cross << G4endl;
00219   return cross;
00220 }
00221 
00222 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00223 
00224 void G4eeToHadronsModel::SampleSecondaries(std::vector<G4DynamicParticle*>* newp,
00225                                            const G4MaterialCutsCouple*,
00226                                            const G4DynamicParticle* dParticle,
00227                                            G4double,
00228                                            G4double)
00229 {
00230   if(crossPerElectron) {
00231     G4double t = dParticle->GetKineticEnergy();
00232     G4double e = 2.0*electron_mass_c2*sqrt(1.0 + 0.5*t/electron_mass_c2);
00233     G4LorentzVector inlv = dParticle->Get4Momentum();
00234     G4ThreeVector inBoost = inlv.boostVector();
00235     if(e > emin) {
00236       G4DynamicParticle* gamma = GenerateCMPhoton(e);
00237       G4LorentzVector gLv = gamma->Get4Momentum();
00238       G4LorentzVector lv(0.0,0.0,0.0,e);
00239       lv -= gLv;
00240       G4double mass = lv.m();
00241       G4ThreeVector boost = lv.boostVector();
00242       const G4ThreeVector dir = gamma->GetMomentumDirection();
00243       model->SampleSecondaries(newp, mass, dir);
00244       G4int np = newp->size();
00245       for(G4int j=0; j<np; j++) {
00246         G4DynamicParticle* dp = (*newp)[j];
00247         G4LorentzVector v = dp->Get4Momentum();
00248         v.boost(boost);
00249         v.boost(inBoost);
00250         dp->Set4Momentum(v);
00251       }
00252       gLv.boost(inBoost);
00253       gamma->Set4Momentum(gLv);
00254       newp->push_back(gamma);
00255     }
00256   }
00257 }
00258 
00259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00260 
00261 void G4eeToHadronsModel::ComputeCMCrossSectionPerElectron()
00262 {
00263   G4bool b;
00264   for(G4int i=0; i<nbins; i++) {
00265     G4double e  = crossPerElectron->GetLowEdgeEnergy(i);
00266     G4double cs = 0.0;
00267     if(i > 0) {
00268       G4double L   = 2.0*log(e/electron_mass_c2);
00269       G4double bt  = 2.0*fine_structure_const*(L - 1.0)/pi;
00270       G4double btm1= bt - 1.0;
00271       G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
00272       G4double s1  = crossBornPerElectron->GetValue(e, b);
00273       G4double e1  = crossPerElectron->GetLowEdgeEnergy(i-1);
00274       G4double x1  = 1. - e1/e;
00275       cs += s1*(del*pow(x1,bt) - bt*(x1 - 0.25*x1*x1));
00276       if(i > 1) {
00277         G4double e2  = e1;
00278         G4double x2  = x1;
00279         G4double s2  = crossBornPerElectron->GetValue(e2, b);
00280         G4double w2  = bt*(del*pow(x2,btm1) - 1.0 + 0.5*x2);
00281         G4double w1;      
00282 
00283         for(G4int j=i-2; j>=0; j--) {
00284           e1  = crossPerElectron->GetLowEdgeEnergy(j);
00285           x1  = 1. - e1/e;
00286           s1  = crossBornPerElectron->GetValue(e1, b);
00287           w1  = bt*(del*pow(x1,btm1) - 1.0 + 0.5*x1);
00288           cs += 0.5*(x1 - x2)*(w2*s2 + w1*s1);
00289           e2 = e1;
00290           x2 = x1;
00291           s2 = s1;
00292           w2 = w1;
00293         }
00294       }
00295     }
00296     crossPerElectron->PutValue(i, cs);
00297     //    G4cout << "e= " << e << "  cs= " << cs << G4endl;
00298   }
00299 }
00300 
00301 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00302 
00303 G4DynamicParticle* G4eeToHadronsModel::GenerateCMPhoton(G4double e)
00304 {
00305   G4bool b;   
00306   G4double x;
00307   G4DynamicParticle* gamma = 0;
00308   G4double L   = 2.0*log(e/electron_mass_c2);
00309   G4double bt  = 2.0*fine_structure_const*(L - 1.)/pi;
00310   G4double btm1= bt - 1.0;
00311   G4double del = 1. + fine_structure_const*(1.5*L + pi*pi/3. -2.)/pi;
00312 
00313   G4double s0 = crossBornPerElectron->GetValue(e, b);
00314   G4double de = (emax - emin)/(G4double)nbins;
00315   G4double x0 = min(de,e - emin)/e;
00316   G4double ds = crossBornPerElectron->GetValue(e, b)
00317               *(del*pow(x0,bt) - bt*(x0 - 0.25*x0*x0));
00318   G4double e1 = e*(1. - x0);
00319 
00320   if(e1 < emax && s0*G4UniformRand()<ds) { 
00321     x = x0*pow(G4UniformRand(),1./bt);
00322   } else {    
00323 
00324     x  = 1. - e1/e;
00325     G4double s1 = crossBornPerElectron->GetValue(e1, b);
00326     G4double w1 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
00327     G4double grej = s1*w1;
00328     G4double f;
00329     //    G4cout << "e= " << e/GeV << " epeak= " << epeak/GeV 
00330     //       << " s1= " << s1 << " w1= " << w1 
00331     //       << " grej= " << grej << G4endl;
00332     // Above emax cross section is 0
00333     if(e1 > emax) {
00334       x  = 1. - emax/e;
00335       G4double s2 = crossBornPerElectron->GetValue(emax, b);
00336       G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
00337       grej = s2*w2;
00338       //  G4cout << "emax= " << emax << " s2= " << s2 << " w2= " << w2 
00339       //   << " grej= " << grej << G4endl;
00340     }
00341 
00342     if(e1 > epeak) {
00343       x  = 1. - epeak/e;
00344       G4double s2 = crossBornPerElectron->GetValue(epeak, b);
00345       G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
00346       grej = max(grej,s2*w2);
00347       //G4cout << "epeak= " << epeak << " s2= " << s2 << " w2= " << w2 
00348       //     << " grej= " << grej << G4endl;
00349     }
00350     G4double xmin = 1. - e1/e;
00351     if(e1 > emax) xmin = 1. - emax/e;
00352     G4double xmax = 1. - emin/e;
00353     do {
00354       x = xmin + G4UniformRand()*(xmax - xmin);
00355       G4double s2 = crossBornPerElectron->GetValue((1.0 - x)*e, b);
00356       G4double w2 = bt*(del*pow(x,btm1) - 1.0 + 0.5*x);
00357       //G4cout << "x= " << x << " xmin= " << xmin << " xmax= " << xmax
00358       //     << " s2= " << s2 << " w2= " << w2 
00359       //           << G4endl;
00360       f = s2*w2;
00361       if(f > grej) {
00362         G4cout << "G4DynamicParticle* G4eeToHadronsModel:WARNING "
00363                << f << " > " << grej << " majorant is`small!" 
00364                << G4endl; 
00365       }
00366     } while (f < grej*G4UniformRand());
00367   }
00368 
00369   G4ThreeVector dir(0.0,0.0,1.0);
00370   gamma = new G4DynamicParticle(theGamma,dir,x*e);
00371   return gamma;
00372 }
00373 
00374 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00375 

Generated on Mon May 27 17:48:05 2013 for Geant4 by  doxygen 1.4.7