CMS 3D CMS Logo

SiG4UniversalFluctuation.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * DISCLAIMER                                                       *
00004 // *                                                                  *
00005 // * The following disclaimer summarizes all the specific disclaimers *
00006 // * of contributors to this software. The specific disclaimers,which *
00007 // * govern, are listed with their locations in:                      *
00008 // *   https://cern.ch/geant4/license                                  *
00009 // *                                                                  *
00010 // * Neither the authors of this software system, nor their employing *
00011 // * institutes,nor the agencies providing financial support for this *
00012 // * work  make  any representation or  warranty, express or implied, *
00013 // * regarding  this  software system or assume any liability for its *
00014 // * use.                                                             *
00015 // *                                                                  *
00016 // * This  code  implementation is the  intellectual property  of the *
00017 // * GEANT4 collaboration.                                            *
00018 // * By copying,  distributing  or modifying the Program (or any work *
00019 // * based  on  the Program)  you indicate  your  acceptance of  this *
00020 // * statement, and all its terms.                                    *
00021 // ********************************************************************
00022 //
00023 // $Id: SiG4UniversalFluctuation.cc,v 1.6 2007/05/10 08:55:16 fambrogl Exp $
00024 // GEANT4 tag $Name: V02-00-00 $
00025 //
00026 // -------------------------------------------------------------------
00027 //
00028 // GEANT4 Class file
00029 //
00030 //
00031 // File name:     G4UniversalFluctuation
00032 //
00033 // Author:        Vladimir Ivanchenko 
00034 // 
00035 // Creation date: 03.01.2002
00036 //
00037 // Modifications: 
00038 //
00039 // 28-12-02 add method Dispersion (V.Ivanchenko)
00040 // 07-02-03 change signature (V.Ivanchenko)
00041 // 13-02-03 Add name (V.Ivanchenko)
00042 // 16-10-03 Changed interface to Initialisation (V.Ivanchenko)
00043 // 07-11-03 Fix problem of rounding of double in G4UniversalFluctuations
00044 // 06-02-04 Add control on big sigma > 2*meanLoss (V.Ivanchenko)
00045 // 26-04-04 Comment out the case of very small step (V.Ivanchenko)
00046 // 07-02-05 define problim = 5.e-3 (mma)
00047 // 03-05-05 conditions of Gaussian fluctuation changed (bugfix)
00048 //          + smearing for very small loss (L.Urban)
00049 //  
00050 // Modified for standalone use in CMSSW. danek k. 2/06        
00051 
00052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00054 
00055 #include "SimTracker/Common/interface/SiG4UniversalFluctuation.h"
00056 #include "CLHEP/Units/SystemOfUnits.h"
00057 #include "CLHEP/Units/PhysicalConstants.h"
00058 #include "CLHEP/Random/RandGaussQ.h"
00059 #include "CLHEP/Random/RandPoisson.h"
00060 #include "CLHEP/Random/RandFlat.h"
00061 #include <math.h>
00062 //#include "G4UniversalFluctuation.hh"
00063 //#include "Randomize.hh"
00064 //#include "G4Poisson.hh"
00065 //#include "G4Step.hh"
00066 //#include "G4Material.hh"
00067 //#include "G4DynamicParticle.hh"
00068 //#include "G4ParticleDefinition.hh"
00069 
00070 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00071 
00072 using namespace std;
00073 
00074 SiG4UniversalFluctuation::SiG4UniversalFluctuation(CLHEP::HepRandomEngine& eng)
00075   :rndEngine(eng),
00076    gaussQDistribution(0),
00077    poissonDistribution(0),
00078    flatDistribution(0),
00079    minNumberInteractionsBohr(10.0),
00080    theBohrBeta2(50.0*keV/proton_mass_c2),
00081    minLoss(10.*eV),
00082    problim(5.e-3),  
00083    alim(10.),
00084    nmaxCont1(4.),
00085    nmaxCont2(16.)
00086 {
00087   sumalim = -log(problim);
00088   //lastMaterial = 0;
00089   
00090   // Add these definitions d.k.
00091   chargeSquare   = 1.;  //Assume all particles have charge 1
00092   // Taken from Geant4 printout, HARDWIRED for Silicon.
00093   ipotFluct = 0.0001736; //material->GetIonisation()->GetMeanExcitationEnergy();  
00094   electronDensity = 6.797E+20; // material->GetElectronDensity();
00095   f1Fluct      = 0.8571; // material->GetIonisation()->GetF1fluct();
00096   f2Fluct      = 0.1429; //material->GetIonisation()->GetF2fluct();
00097   e1Fluct      = 0.000116;// material->GetIonisation()->GetEnergy1fluct();
00098   e2Fluct      = 0.00196; //material->GetIonisation()->GetEnergy2fluct();
00099   e1LogFluct   = -9.063;  //material->GetIonisation()->GetLogEnergy1fluct();
00100   e2LogFluct   = -6.235;  //material->GetIonisation()->GetLogEnergy2fluct();
00101   rateFluct    = 0.4;     //material->GetIonisation()->GetRateionexcfluct();
00102   ipotLogFluct = -8.659;  //material->GetIonisation()->GetLogMeanExcEnergy();
00103   e0 = 1.E-5;             //material->GetIonisation()->GetEnergy0fluct();
00104   
00105   gaussQDistribution = new CLHEP::RandGaussQ(rndEngine);
00106   poissonDistribution = new CLHEP::RandPoisson(rndEngine);
00107   flatDistribution = new CLHEP::RandFlat(rndEngine);
00108 
00109   //cout << " init new fluct +++++++++++++++++++++++++++++++++++++++++"<<endl;
00110 }
00111 
00112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00113 // The main dedx fluctuation routine.
00114 // Arguments: momentum in MeV/c, mass in MeV, delta ray cut (tmax) in
00115 // MeV, silicon thickness in mm, mean eloss in MeV.
00116 
00117 SiG4UniversalFluctuation::~SiG4UniversalFluctuation()
00118 {
00119   delete gaussQDistribution;
00120   delete poissonDistribution;
00121   delete flatDistribution;
00122 
00123 }
00124 
00125 
00126 double SiG4UniversalFluctuation::SampleFluctuations(const double momentum,
00127                                                     const double mass,
00128                                                     double& tmax,
00129                                                     const double length,
00130                                                     const double meanLoss)
00131 {
00132 // Calculate actual loss from the mean loss.
00133 // The model used to get the fluctuations is essentially the same
00134 // as in Glandz in Geant3 (Cern program library W5013, phys332).
00135 // L. Urban et al. NIM A362, p.416 (1995) and Geant4 Physics Reference Manual
00136 
00137   // shortcut for very very small loss (out of validity of the model)
00138   //
00139   if (meanLoss < minLoss) return meanLoss;
00140 
00141   //if(!particle) InitialiseMe(dp->GetDefinition());
00142   //G4double tau   = dp->GetKineticEnergy()/particleMass;
00143   //G4double gam   = tau + 1.0;
00144   //G4double gam2  = gam*gam;
00145   //G4double beta2 = tau*(tau + 2.0)/gam2;
00146 
00147   particleMass   = mass; // dp->GetMass();
00148   double gam2   = (momentum*momentum)/(particleMass*particleMass) + 1.0;
00149   double beta2 = 1.0 - 1.0/gam2;
00150   double gam = sqrt(gam2);
00151 
00152   double loss(0.), siga(0.);
00153   
00154   // Gaussian regime
00155   // for heavy particles only and conditions
00156   // for Gauusian fluct. has been changed 
00157   //
00158   if ((particleMass > electron_mass_c2) &&
00159       (meanLoss >= minNumberInteractionsBohr*tmax))
00160   {
00161     double massrate = electron_mass_c2/particleMass ;
00162     double tmaxkine = 2.*electron_mass_c2*beta2*gam2/
00163       (1.+massrate*(2.*gam+massrate)) ;
00164     if (tmaxkine <= 2.*tmax)   
00165     {
00166       //electronDensity = material->GetElectronDensity();
00167       siga  = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
00168                                 * electronDensity * chargeSquare;
00169       siga = sqrt(siga);
00170       double twomeanLoss = meanLoss + meanLoss;
00171       if (twomeanLoss < siga) {
00172         double x;
00173         do {
00174           loss = twomeanLoss*flatDistribution->fire();
00175           x = (loss - meanLoss)/siga;
00176         } while (1.0 - 0.5*x*x < flatDistribution->fire());
00177       } else {
00178         do {
00179           loss = gaussQDistribution->fire(meanLoss,siga);
00180         } while (loss < 0. || loss > twomeanLoss);
00181       }
00182       return loss;
00183     }
00184   }
00185 
00186   // Glandz regime : initialisation
00187   //
00188 //   if (material != lastMaterial) {
00189 //     f1Fluct      = material->GetIonisation()->GetF1fluct();
00190 //     f2Fluct      = material->GetIonisation()->GetF2fluct();
00191 //     e1Fluct      = material->GetIonisation()->GetEnergy1fluct();
00192 //     e2Fluct      = material->GetIonisation()->GetEnergy2fluct();
00193 //     e1LogFluct   = material->GetIonisation()->GetLogEnergy1fluct();
00194 //     e2LogFluct   = material->GetIonisation()->GetLogEnergy2fluct();
00195 //     rateFluct    = material->GetIonisation()->GetRateionexcfluct();
00196 //     ipotFluct    = material->GetIonisation()->GetMeanExcitationEnergy();
00197 //     ipotLogFluct = material->GetIonisation()->GetLogMeanExcEnergy();
00198 //     lastMaterial = material;
00199 //   }
00200 
00201   double a1 = 0. , a2 = 0., a3 = 0. ;
00202   double p1,p2,p3;
00203   double rate = rateFluct ;
00204 
00205   double w1 = tmax/ipotFluct;
00206   double w2 = log(2.*electron_mass_c2*beta2*gam2)-beta2;
00207 
00208   if(w2 > ipotLogFluct)
00209   {
00210     double C = meanLoss*(1.-rateFluct)/(w2-ipotLogFluct);
00211     a1 = C*f1Fluct*(w2-e1LogFluct)/e1Fluct;
00212     a2 = C*f2Fluct*(w2-e2LogFluct)/e2Fluct;
00213     if(a2 < 0.)
00214     {
00215       a1 = 0. ;
00216       a2 = 0. ;
00217       rate = 1. ;  
00218     }
00219   }
00220   else
00221   {
00222     rate = 1. ;
00223   }
00224 
00225   a3 = rate*meanLoss*(tmax-ipotFluct)/(ipotFluct*tmax*log(w1));
00226 
00227   double suma = a1+a2+a3;
00228   
00229   // Glandz regime
00230   //
00231   if (suma > sumalim)
00232   {
00233     p1 = 0., p2 = 0 ;
00234     if((a1+a2) > 0.)
00235     {
00236       // excitation type 1
00237       if (a1>alim) {
00238         siga=sqrt(a1) ;
00239         p1 = max(0.,gaussQDistribution->fire(a1,siga)+0.5);
00240       } else {
00241         p1 = double(poissonDistribution->fire(a1));
00242       }
00243     
00244       // excitation type 2
00245       if (a2>alim) {
00246         siga=sqrt(a2) ;
00247         p2 = max(0.,gaussQDistribution->fire(a2,siga)+0.5);
00248       } else {
00249         p2 = double(poissonDistribution->fire(a2));
00250       }
00251     
00252       loss = p1*e1Fluct+p2*e2Fluct;
00253  
00254       // smearing to avoid unphysical peaks
00255       if (p2 > 0.)
00256         loss += (1.-2.*flatDistribution->fire())*e2Fluct;
00257       else if (loss>0.)
00258         loss += (1.-2.*flatDistribution->fire())*e1Fluct;   
00259       if (loss < 0.) loss = 0.0;
00260     }
00261 
00262     // ionisation
00263     if (a3 > 0.) {
00264       if (a3>alim) {
00265         siga=sqrt(a3) ;
00266         p3 = max(0.,gaussQDistribution->fire(a3,siga)+0.5);
00267       } else {
00268         p3 = double(poissonDistribution->fire(a3));
00269       }
00270       double lossc = 0.;
00271       if (p3 > 0) {
00272         double na = 0.; 
00273         double alfa = 1.;
00274         if (p3 > nmaxCont2) {
00275           double rfac   = p3/(nmaxCont2+p3);
00276           double namean = p3*rfac;
00277           double sa     = nmaxCont1*rfac;
00278           na              = gaussQDistribution->fire(namean,sa);
00279           if (na > 0.) {
00280             alfa   = w1*(nmaxCont2+p3)/(w1*nmaxCont2+p3);
00281             double alfa1  = alfa*log(alfa)/(alfa-1.);
00282             double ea     = na*ipotFluct*alfa1;
00283             double sea    = ipotFluct*sqrt(na*(alfa-alfa1*alfa1));
00284             lossc += gaussQDistribution->fire(ea,sea);
00285           }
00286         }
00287 
00288         if (p3 > na) {
00289           w2 = alfa*ipotFluct;
00290           double w  = (tmax-w2)/tmax;
00291           int    nb = int(p3-na);
00292           for (int k=0; k<nb; k++) lossc += w2/(1.-w*flatDistribution->fire());
00293         }
00294       }        
00295       loss += lossc;  
00296     }
00297     return loss;
00298   }
00299   
00300   // suma < sumalim;  very small energy loss;  
00301   //
00302   //double e0 = material->GetIonisation()->GetEnergy0fluct();
00303 
00304   a3 = meanLoss*(tmax-e0)/(tmax*e0*log(tmax/e0));
00305   if (a3 > alim)
00306   {
00307     siga=sqrt(a3);
00308     p3 = max(0.,gaussQDistribution->fire(a3,siga)+0.5);
00309   } else {
00310     p3 = double(poissonDistribution->fire(a3));
00311   }
00312   if (p3 > 0.) {
00313     double w = (tmax-e0)/tmax;
00314     double corrfac = 1.;
00315     if (p3 > nmaxCont2) {
00316       corrfac = p3/nmaxCont2;
00317       p3 = nmaxCont2;
00318     } 
00319     int ip3 = (int)p3;
00320     for (int i=0; i<ip3; i++) loss += 1./(1.-w*flatDistribution->fire());
00321     loss *= e0*corrfac;
00322     // smearing for losses near to e0
00323     if(p3 <= 2.)
00324       loss += e0*(1.-2.*flatDistribution->fire()) ;
00325    }
00326     
00327    return loss;
00328 }
00329 
00330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00331 // G4double SiG4UniversalFluctuation::Dispersion(
00332 //                           const G4Material* material,
00333 //                           const G4DynamicParticle* dp,
00334 //                              G4double& tmax,
00335 //                              G4double& length)
00336 // {
00337 //   if(!particle) InitialiseMe(dp->GetDefinition());
00338 
00339 //   electronDensity = material->GetElectronDensity();
00340 
00341 //   G4double gam   = (dp->GetKineticEnergy())/particleMass + 1.0;
00342 //   G4double beta2 = 1.0 - 1.0/(gam*gam);
00343 
00344 //   G4double siga  = (1.0/beta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
00345 //                  * electronDensity * chargeSquare;
00346 
00347 //   return siga;
00348 // }
00349 
00350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

Generated on Tue Jun 9 17:47:52 2009 for CMSSW by  doxygen 1.5.4