CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/SimMuon/CSCDigitizer/src/CSCGasCollisions.cc

Go to the documentation of this file.
00001 // This implements the CSC gas ionization simulation.
00002 // It requires the GEANT3-generated tables in the input file 
00003 // (default) MuonData/EndcapData/collisions.dat
00004 
00005 // Outstanding problems to fix 15-Oct-2003
00006 // - ework must be re-tuned once eion or ework are removed before delta e range estimation.
00007 // - logic of long-range delta e's must be revisited.
00008 // - The code here needs to have diagnostic output removed, or reduced. PARTICULARLY filling of std::vectors!
00009 // - 'gap' is a misnomer, when simhit entry and exit don't coincide with gap edges
00010 //  so the CSCCrossGap might better be named as a simhit-related thing. 
00011 // 22-Jan-2004 Corrected position of trap for 'infinite' loop while
00012 //    generating steps. Output files (debugV-flagged only) require SC flag.
00013 //    Increase deCut from 10 keV to 100 keV to accomodate protons!
00014 
00015 #include "FWCore/Utilities/interface/Exception.h"
00016 #include "DataFormats/GeometryVector/interface/LocalVector.h"
00017 #include "SimMuon/CSCDigitizer/src/CSCGasCollisions.h"
00018 #include "Utilities/General/interface/FileInPath.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 #include <iostream>
00021 #include <fstream>
00022 #include <cassert>
00023 
00024 // stdlib math functions
00025 // for 'abs':
00026 #include <cstdlib>
00027 // for usual math functions:
00028 #include <cmath>
00029 
00030 // stdlib container trickery
00031 // for 'for_each':
00032 #include <algorithm>
00033 // for 'greater' and 'less':
00034 #include <functional>
00035 // for 'accumulate':
00036 #include <numeric>
00037 #include <iterator>
00038 
00039 using namespace std;
00040 /* Gas mixture is Ar/CO2/CF4 = 40/50/10
00041 We'll use the ioinzation energies
00042 Ar    15.8 eV
00043 CO2   13.7 eV
00044 CF4    17.8
00045 to arrive at a weighted average of 14.95 */
00046 
00047 CSCGasCollisions::CSCGasCollisions() : me("CSCGasCollisions"),
00048   gasDensity( 2.1416e-03 ), 
00049   deCut( 1.e05 ), eion( 14.95 ), ework( 34.0 ), clusterExtent( 0.001 ),
00050   theGammaBins(N_GAMMA, 0.), theEnergyBins(N_ENERGY, 0.), 
00051   theCollisionTable(N_ENTRIES, 0.), theCrossGap( 0 ),
00052   theParticleDataTable(0),
00053   theRandFlat(0),
00054   theRandExponential(0),
00055   saveGasCollisions ( false )
00056 {
00057 
00058   edm::LogInfo(me) << "Constructing a " << me << ":"; 
00059   edm::LogInfo(me) << "gas density = " << gasDensity << " g/cm3";
00060   edm::LogInfo(me) << "max eloss per collision allowed = " << deCut/1000. 
00061        << " keV (for higher elosses, hits should have been simulated.)";
00062   edm::LogInfo(me) << "ionization threshold = " << eion << " eV";
00063   edm::LogInfo(me) << "effective work function = " << ework << " eV";
00064   edm::LogInfo(me) << "cluster extent = " << clusterExtent*1.e04 << " micrometres";
00065   edm::LogInfo(me) << "Save gas collision info to files when using debugV? " << saveGasCollisions;
00066 
00067   readCollisionTable();
00068 }
00069 
00070 CSCGasCollisions::~CSCGasCollisions() {
00071   edm::LogInfo(me) << "Destructing a " << me;
00072   delete theCrossGap;
00073   delete theRandFlat;
00074   delete theRandExponential;
00075 }
00076 
00077 void CSCGasCollisions::readCollisionTable() {
00078 
00079   // I'd prefer to allow comments in the data file which means
00080   // complications of reading line-by-line and then item-by-item.
00081 
00082   // Is float OK? Or do I need double?
00083 
00084   // Do I need any sort of error trapping?
00085 
00086   // We use the default CMSSW data file path
00087   // We use the default file name 'collisions.dat'
00088    
00089   // This can be reset in .orcarc by SimpleConfigurable
00090   //    Muon:Endcap:CollisionsFile
00091 
00092   string path( getenv( "CMSSW_SEARCH_PATH" ) );
00093   // TODO make configurable
00094   string colliFile = "SimMuon/CSCDigitizer/data/collisions.dat";
00095   FileInPath f1( path, colliFile );
00096   if (f1() == 0 ) {
00097     string errorMessage = "Input file " + colliFile + "not found";
00098     edm::LogError("CSCGasCollisions") << errorMessage << " in path " << path
00099           << "\nSet Muon:Endcap:CollisionsFile in .orcarc to the "
00100          " location of the file relative to ORCA_DATA_PATH." ;
00101     throw cms::Exception( " Endcap Muon gas collisions data file not found.");
00102   }
00103   else {
00104     edm::LogInfo(me) << ": reading " << f1.name();
00105   }
00106 
00107   ifstream & fin = *f1();
00108 
00109   if (fin == 0) {
00110     string errorMessage = "Cannot open input file " + path + colliFile;
00111     edm::LogError("CSCGasCollisions") << errorMessage;
00112     throw cms::Exception(errorMessage);
00113   }
00114 
00115   fin.clear( );                  // Clear eof read status
00116   fin.seekg( 0, ios::beg );      // Position at start of file
00117 
00118   // @@ We had better have the right sizes everywhere or all
00119   // hell will break loose. There's no trapping.
00120   LogTrace(me) << "Reading gamma bins";
00121   int j = 0;
00122   for(int i = 0; i<N_GAMMA; ++i) {
00123     fin >> theGammaBins[i];
00124     LogTrace(me) << ++j << "   " << theGammaBins[i];
00125   }
00126 
00127   j = 0;
00128   LogTrace(me) << "Reading energy bins \n";
00129   for(int i = 0; i<N_ENERGY; ++i) {
00130      fin >> theEnergyBins[i];
00131      LogTrace(me) << ++j << "   " << theEnergyBins[i];
00132   }
00133 
00134   j = 0;
00135   LogTrace(me) << "Reading collisions table \n";
00136 
00137   for(int i = 0; i<N_ENTRIES; ++i) {
00138      fin >> theCollisionTable[i];
00139      LogTrace(me) << ++j << "   " << theCollisionTable[i];
00140   }
00141 
00142   fin.close();
00143 }
00144 
00145 
00146 void CSCGasCollisions::setParticleDataTable(const ParticleDataTable * pdt)
00147 {
00148   theParticleDataTable = pdt;
00149 }
00150 
00151 
00152 void CSCGasCollisions::setRandomEngine(CLHEP::HepRandomEngine & engine)
00153 {
00154   theRandFlat = new CLHEP::RandFlat(engine);
00155   theRandExponential = new CLHEP::RandExponential(engine);
00156 }
00157 
00158 
00159 void CSCGasCollisions::simulate( const PSimHit& simHit, 
00160   std::vector<LocalPoint>& positions, std::vector<int>& electrons ) {
00161 
00162   const float epsilonL = 0.01;                     // Shortness of simhit 'length'
00163   //  const float max_gap_z = 1.5;                     // Gas gaps are 0.5 or 1.0 cm
00164 
00165   // Note that what I call the 'gap' may in fact be the 'length' of a PSimHit which
00166   // does not start and end on the gap edges. This confuses the nomenclature at least.
00167 
00168   double mom    = simHit.pabs();
00169   int iam       = simHit.particleType();           // PDG type
00170   delete theCrossGap;                              // before building new one
00171   assert(theParticleDataTable != 0);
00172   ParticleData const * particle = theParticleDataTable->particle( simHit.particleType() );
00173   double mass = 0.105658; // assume a muon
00174   if(particle == 0)
00175   {
00176      edm::LogError("CSCGasCollisions") << "Cannot find particle of type " << simHit.particleType()
00177                                        << " in the PDT";
00178   }
00179   else 
00180   {
00181     mass = particle->mass();
00182   }
00183 
00184   theCrossGap   = new CSCCrossGap( mass, mom, simHit.exitPoint() - simHit.entryPoint() );
00185   float gapSize = theCrossGap->length();
00186 
00187   // Test the simhit 'length' (beware of angular effects)
00188   //  if ( gapSize <= epsilonL || gapSize > max_gap_z ) {
00189   if ( gapSize <= epsilonL ) {
00190     LogTrace(me) << ": WARNING! simhit entry and exit are very close: \n"
00191       << "\n entry = " << simHit.entryPoint()
00192       << "\n exit  = " << simHit.exitPoint()
00193       << "\n particle type = " << iam << ", momentum = " << mom 
00194       << ", gap length = " << gapSize;
00195     return; //@@ Just skip this PSimHit
00196   }
00197 
00198   // Interpolate the table for current gamma value
00199   // Extract collisions binned by energy loss values, for this gamma
00200   std::vector<float> collisions(N_ENERGY);
00201   double loggam = theCrossGap->logGamma(); 
00202   fillCollisionsForThisGamma( static_cast<float>(loggam), collisions );
00203 
00204   double anmin = exp(collisions[N_ENERGY-1]);
00205   double anmax = exp(collisions[0]);
00206   double amu   = anmax - anmin;
00207 
00208   LogTrace(me) << "collisions extremes = " << collisions[N_ENERGY-1]
00209      << ", " << collisions[0] << "\n"
00210      << "anmin = " << anmin << ", anmax = " << anmax << "\n"
00211      << "amu = " << amu << "\n";
00212 
00213   float dedx       = 0.; // total energy loss
00214   double sum_steps = 0.; // total distance across gap (along simhit direction)
00215   int n_steps      = 0;  // no. of steps/primary collisions
00216   int n_try        = 0;  // no. of tries to generate steps
00217   double step      = -1.; // Sentinel for start
00218 
00219   LocalPoint layerLocalPoint( simHit.entryPoint() );
00220 
00221   // step/primary collision loop
00222   while ( sum_steps < gapSize) {
00223     ++n_try;
00224     if ( n_try > MAX_STEPS ) {
00225           LogTrace(me) << ": n_try=" << n_try << " is too large. Skip simhit."
00226            << "\n particle type=" << iam << ", momentum= " << mom
00227            << "\n gapSize=" << gapSize << ", last step=" << step 
00228            << ", sum_steps=" << sum_steps << ", n_steps=" << n_steps;
00229         break;
00230     }
00231     step = generateStep( amu );
00232     if ( sum_steps + step > gapSize ) break;
00233 
00234     float eloss = generateEnergyLoss( amu, anmin, anmax, collisions );
00235 
00236     // Is the eloss too large? (then GEANT should have produced hits!)
00237     if ( eloss > deCut ) {
00238       LogTrace(me) << "eloss > " << deCut << " = " << eloss;
00239       continue; // to generate another collision/step
00240     }
00241 
00242     dedx      += eloss; // the energy lost from the ionizing particle
00243     sum_steps += step;  // the position of the ionizing particle
00244     ++n_steps;          // the number of primary collisions
00245     
00246     if (n_steps > MAX_STEPS ) { // Extra-careful trap for bizarreness
00247         edm::LogInfo(me) << ": n_steps=" << n_steps << " is too large. Skip simhit.";
00248         edm::LogInfo(me) << "particle type=" << iam << ", momentum= " << mom;
00249         edm::LogInfo(me) << "gapSize=" << gapSize << ", last step=" << step 
00250              << ", sum_steps=" << sum_steps;
00251         break;
00252     }
00253     LogTrace(me) << "sum_steps = " << sum_steps << ", dedx = " << dedx;
00254     
00255     // Generate ionization. 
00256     // eion is the minimum energy at which ionization can occur in the gas
00257     if ( eloss > eion  ) {
00258       layerLocalPoint += step * theCrossGap->unitVector(); // local point where the collision occurs
00259       ionize( eloss, layerLocalPoint );
00260     }
00261     else {
00262       LogTrace(me) << "Energy available = " << eloss <<
00263         ", too low for ionization.";
00264     }
00265 
00266   } // step/collision loop
00267 
00268   //TODO port this
00269   //if ( debugV ) writeSummary( n_steps, sum_steps, dedx, simHit.energyLoss() );
00270 
00271   // Return values in two container arguments
00272   positions = theCrossGap->ionClusters();
00273   electrons = theCrossGap->electrons();
00274   return;
00275 }
00276 
00277 double CSCGasCollisions::generateStep( double avCollisions ) const
00278 {
00279 // Generate a m.f.p.  (1/avCollisions = cm/collision)
00280   double step = (theRandExponential->fire())/avCollisions;
00281 
00282 // Without using CLHEP: approx random exponential by...
00283 //    double da = double(rand())/double(RAND_MAX);
00284 //    double step = -log(1.-da)/avCollisions;
00285 
00286     LogTrace(me)  << " step = " << step;
00287     // Next line only used to fill a container of 'step's for later diagnostic dumps
00288     //if ( debugV ) theCrossGap->addStep( step );
00289     return step;
00290 }
00291 
00292 float CSCGasCollisions::generateEnergyLoss( double avCollisions,
00293    double anmin, double anmax, const std::vector<float>& collisions ) const
00294 {
00295 // Generate a no. of collisions between collisions[0] and [N_ENERGY-1]
00296     float lnColl = log(theRandFlat->fire(anmin, anmax));
00297 
00298 // Without using CLHEP: approx random between anmin and anmax
00299 //    double ra = double(rand())/double(RAND_MAX)*avCollisions;
00300 //    cout << "ra = " << ra << std::endl;
00301 //    float lnColl = static_cast<float>( log( ra ) );
00302 
00303     // Find energy loss for that number
00304     float lnE    = lnEnergyLoss( lnColl, collisions );
00305     float eloss  = exp(lnE);
00306     // Compensate if gamma was actually below 1.1
00307     if ( theCrossGap->gamma() < 1.1 ) eloss = eloss * 0.173554/theCrossGap->beta2();
00308     LogTrace(me) << "eloss = " << eloss;
00309     // Next line only used to fill container of eloss's for later diagnostic dumps
00310     // if ( debugV ) theCrossGap->addEloss( eloss );
00311     return eloss;
00312 }
00313 
00314 void  CSCGasCollisions::ionize( double energyAvailable, LocalPoint startHere ) const
00315 {
00316   while ( energyAvailable > eion ) {
00317     LogTrace(me) << "     NEW CLUSTER " << theCrossGap->noOfClusters() + 1 <<
00318                " AT " << startHere;
00319     LocalPoint newCluster( startHere );
00320     theCrossGap->addCluster(newCluster);
00321 
00322   //@@ I consider NOT subtracting eion before calculating range to be a bug.
00323   //@@ But this changes tuning of the algorithm so leave it until after the big rush to 7_5_0
00324   //@@  energyAvailable -= eion; 
00325 
00326   // Sauli CERN 77-09: delta e range with E in MeV (Sauli references Kobetich & Katz 1968,
00327   // but I cannot find this expression in that set of papers.)
00328   // Take HALF that range. //@@ Why? Why not...
00329     double range = 0.5 * (0.71/gasDensity)*pow( energyAvailable*1.E-6, 1.72);
00330     LogTrace(me) << " range = " << range;
00331     if ( range < clusterExtent ) {
00332 
00333       // short-range delta e
00334           // How many electrons can we make? Now use *average* energy for ionization (not *minimum*)
00335       int nelec = static_cast<int>(energyAvailable/ework);
00336       LogTrace(me) << "s-r delta energy in = " << energyAvailable;
00337       //energyAvailable -= nelec*(energyAvailable/ework);
00338       energyAvailable -= nelec*ework;
00339             // If still above eion (minimum, not average) add one more e
00340       if ( energyAvailable > eion ) {
00341         ++nelec;
00342         energyAvailable -= eion;
00343       }
00344       LogTrace(me) << "s-r delta energy out = " << energyAvailable << ", nelec = " << nelec;
00345       theCrossGap->addElectrons( nelec );
00346       break;
00347 
00348      }
00349      else {
00350       // long-range delta e
00351          LogTrace(me) << "l-r delta \n"
00352              << "no. of electrons in cluster now = " << theCrossGap->noOfElectrons();
00353          theCrossGap->addElectrons( 1 ); // Position is at startHere still
00354 
00355           bool new_range = false;
00356           while ( !new_range && (energyAvailable>ework) ) {
00357             energyAvailable -= ework;
00358             while ( energyAvailable > eion ) {
00359               double range2 = 0.5 * 0.71/gasDensity*pow( 1.E-6*energyAvailable, 1.72);
00360               double drange = range - range2;
00361               LogTrace(me) << "  energy left = " << energyAvailable << 
00362                         ", range2 = " << range2 << ", drange = " << drange;
00363               if ( drange < clusterExtent ) {
00364                 theCrossGap->addElectronToBack(); // increment last element
00365               }
00366               else {
00367                 startHere += drange*theCrossGap->unitVector(); // update delta e start position
00368                 range = range2;                       // update range
00369                 new_range = true;                     // Test range again
00370                 LogTrace(me) << "reset range to range2 and iterate";
00371               }
00372               break; // out of inner while energyAvailable>eion
00373 
00374             } // inner while energyAvailable>eion
00375 
00376           } // while !new_range && energyAvailable>ework
00377 
00378           // energyAvailable now less than ework, but still may be over eion...add an e
00379           if ( energyAvailable > eion ) {
00380             energyAvailable -= ework; // yes, it may go negative
00381             theCrossGap->addElectronToBack(); // add one more e
00382           }
00383 
00384         } // if range
00385 
00386       } // outer while energyAvailable>eion
00387 }
00388 
00389 void CSCGasCollisions::writeSummary( int n_steps, double sum_steps, float dedx, float simHiteloss ) const
00390 {
00391   std::vector<LocalPoint> ion_clusters = theCrossGap->ionClusters();
00392   std::vector<int> electrons             = theCrossGap->electrons();
00393   std::vector<float> elosses             = theCrossGap->eLossPerStep();
00394   std::vector<double> steps              = theCrossGap->stepLengths();
00395 
00396     cout << "------------------" << std::endl;
00397     cout << "AFTER CROSSING GAP" << std::endl;
00398     cout << "No. of steps = " << n_steps << std::endl;
00399     cout << "Check:  stored steps = " << theCrossGap->noOfSteps() << std::endl;
00400 
00401     cout << "Lengths of steps: " << std::endl;
00402     std::copy( steps.begin(), steps.end(), std::ostream_iterator<float>(cout,"\n"));
00403     cout << std::endl;
00404 
00405     if ( saveGasCollisions ) {
00406       ofstream of0("osteplen.dat",ios::app);
00407       std::copy( steps.begin(), steps.end(), std::ostream_iterator<float>(of0,"\n"));
00408     }
00409 
00410     cout << "Total sum of steps = " << sum_steps << std::endl;
00411     if ( n_steps > 0 ) cout << "Average step length = " << 
00412        sum_steps/float(n_steps) << std::endl;
00413     cout << std::endl;
00414     
00415     cout << "Energy loss per collision:" << std::endl;
00416     std::copy( elosses.begin(), elosses.end(), std::ostream_iterator<float>(cout,"\n"));
00417     cout << std::endl;
00418 
00419     if ( saveGasCollisions ) {
00420       ofstream of1("olperc.dat",ios::app);
00421       std::copy( elosses.begin(), elosses.end(), std::ostream_iterator<float>(of1,"\n"));
00422     }
00423 
00424     cout << "Total energy loss across gap = " << dedx << " eV = " <<
00425        dedx/1000. << " keV" << std::endl;
00426     int n_ic = count_if( elosses.begin(), elosses.end(),
00427                      bind2nd(greater<float>(), eion) );
00428     cout << "No. of primary ionizing collisions across gap = " << n_ic << std::endl;
00429     if ( n_steps > 0 ) cout << "Average energy loss/collision = " << 
00430        dedx/float(n_steps) << " eV" << std::endl;
00431     cout << std::endl;
00432 
00433     cout << "No. of ion clusters = " << ion_clusters.size() << std::endl;
00434     cout << "Positions of clusters:" << std::endl;
00435     std::copy( ion_clusters.begin(), ion_clusters.end(), 
00436          std::ostream_iterator<LocalPoint>(cout,"\n"));
00437     cout << std::endl;
00438 
00439     if ( saveGasCollisions ) {
00440       ofstream of2("oclpos.dat",ios::app);
00441       std::copy( ion_clusters.begin(), ion_clusters.end(), 
00442          std::ostream_iterator<LocalPoint>(of2,"\n"));
00443     }
00444 
00445     cout << "No. of electrons per cluster:" << std::endl;
00446     std::copy( electrons.begin(), electrons.end(), std::ostream_iterator<int>(cout,"\n"));
00447     cout << std::endl;
00448 
00449     if ( saveGasCollisions ) {
00450       ofstream of3("oepercl.dat",ios::app);
00451       std::copy( electrons.begin(), electrons.end(), std::ostream_iterator<int>(of3,"\n"));
00452     }
00453 
00454     // Check for zero-e clusters
00455     std::vector<int>::const_iterator bigger = find(electrons.begin(),
00456            electrons.end(), 0 );
00457     if ( bigger != electrons.end() ) {
00458        cout << "Error! There is a cluster with 0 electrons." << std::endl;
00459     }
00460     int n_e = accumulate(electrons.begin(), electrons.end(), 0 );
00461     if ( n_steps > 0 ) {
00462       cout << "#        cm       cm             keV               eV          eV       eV     keV" << std::endl;
00463       cout << " " << n_steps << "  " << sum_steps << " " << sum_steps/float(n_steps) << "    " <<
00464          ion_clusters.size() << "    " <<
00465          dedx/1000. << "    " << n_ic << "    " << dedx/float(n_steps) << "  " << n_e << "  " <<
00466          dedx/float(n_e) << " " << float(n_e)/float(ion_clusters.size()) << " " << simHiteloss*1.E6 << std::endl;
00467     }
00468 }
00469 
00470 float CSCGasCollisions::lnEnergyLoss( float lnCollisions, 
00471    const std::vector<float>& collisions ) const {
00472 
00473   float lnE = -1.;
00474 
00475    // Find collision[] bin in which lnCollisions falls
00476     std::vector<float>::const_iterator it = find(collisions.begin(),
00477       collisions.end(), lnCollisions );
00478 
00479     if ( it != collisions.end() ) {
00480       // found the value
00481       std::vector<float>::difference_type ihi = it - collisions.begin();
00482       LogTrace(me) << ": using one energy bin " << ihi << " = " 
00483                          << theEnergyBins[ihi]
00484                          << " for lnCollisions = " << lnCollisions;
00485       lnE = theEnergyBins[ihi];
00486     }
00487     else {
00488       // interpolate the value
00489       std::vector<float>::const_iterator loside = find_if(collisions.begin(),
00490         collisions.end(), bind2nd(less<float>(), lnCollisions));
00491       std::vector<float>::difference_type ilo = loside - collisions.begin();
00492       if ( ilo > 0 ) {
00493         LogTrace(me) << ": using energy bin " 
00494                            << ilo-1 << " and " << ilo;
00495         lnE = theEnergyBins[ilo-1] + (lnCollisions-collisions[ilo-1])*
00496           (theEnergyBins[ilo]-theEnergyBins[ilo-1]) /
00497              (collisions[ilo]-collisions[ilo-1]);
00498       }
00499       else {
00500         LogTrace(me) << ": using one energy bin 0 = " 
00501                            << theEnergyBins[0]
00502                            << " for lnCollisions = " << lnCollisions;
00503         lnE = theEnergyBins[0]; //@@ WHAT ELSE TO DO?
00504       }
00505     }
00506 
00507     return lnE;
00508 }
00509 
00510 void CSCGasCollisions::fillCollisionsForThisGamma( float logGamma,
00511    std::vector<float>& collisions ) const
00512 {
00513   std::vector<float>::const_iterator bigger = find_if(theGammaBins.begin(),
00514     theGammaBins.end(), bind2nd(greater<float>(), logGamma));
00515 
00516   if ( bigger == theGammaBins.end() ) {
00517     // use highest bin
00518     LogTrace(me) << ": using highest gamma bin" 
00519                       << " for logGamma = " << logGamma;
00520     for (int i=0; i<N_ENERGY; ++i)
00521       collisions[i] = theCollisionTable[i*N_GAMMA];
00522   }
00523   else {
00524     // use bigger and its lower neighbour
00525     std::vector<float>::difference_type ihi = bigger - theGammaBins.begin();
00526     if ( ihi > 0 ) {
00527       double dlg2 = *bigger--; // and decrement after deref
00528       //LogTrace(me) << ": using gamma bins " 
00529       //                   << ihi-1 << " and " << ihi;
00530       double dlg1 = *bigger;   // now the preceding element
00531       double dlg = (logGamma-dlg1)/(dlg2-dlg1);
00532       double omdlg = 1. - dlg;
00533       for (int i=0; i<N_ENERGY; ++i)
00534         collisions[i] = theCollisionTable[i*N_GAMMA+ihi-1]*omdlg +
00535           theCollisionTable[i*N_GAMMA+ihi]*dlg;
00536     }
00537     else {
00538       // bigger has no lower neighbour
00539        LogTrace(me) << ": using lowest gamma bin" 
00540                          << " for logGamma = " << logGamma;
00541 
00542       for (int i=0; i<N_ENERGY; ++i)
00543         collisions[i] = theCollisionTable[i*N_GAMMA];
00544     }
00545   }
00546 }