CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FastSimulation/MaterialEffects/src/NuclearInteractionSimulator.cc

Go to the documentation of this file.
00001 //Framework Headers
00002 #include "FWCore/Utilities/interface/Exception.h"
00003 #include "FWCore/ParameterSet/interface/FileInPath.h"
00004 
00005 #include "FastSimulation/MaterialEffects/interface/NuclearInteractionSimulator.h"
00006 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00007 
00008 #include "FastSimDataFormats/NuclearInteractions/interface/NUEvent.h"
00009 
00010 //#include "DQMServices/Core/interface/DaqMonitorBEInterface.h"
00011 //#include "FWCore/ServiceRegistry/interface/Service.h"
00012 
00013 #include <iostream>
00014 #include <sys/stat.h>
00015 #include <cmath>
00016 #include "TFile.h"
00017 #include "TTree.h"
00018 #include "TROOT.h"
00019 
00020 NuclearInteractionSimulator::NuclearInteractionSimulator(
00021   std::vector<double>& pionEnergies,
00022   std::vector<int>& pionTypes,
00023   std::vector<std::string>& pionNames,
00024   std::vector<double>& pionMasses,
00025   std::vector<double>& pionPMin,
00026   double pionEnergy,
00027   std::vector<double>& lengthRatio,
00028   std::vector< std::vector<double> >& ratios,
00029   std::map<int,int >& idMap,
00030   std::string inputFile,
00031   unsigned int distAlgo,
00032   double distCut,
00033   const RandomEngine* engine) 
00034   :
00035   MaterialEffectsSimulator(engine),
00036   thePionEN(pionEnergies),
00037   thePionID(pionTypes),
00038   thePionNA(pionNames),
00039   thePionMA(pionMasses),
00040   thePionPMin(pionPMin),
00041   thePionEnergy(pionEnergy),
00042   theLengthRatio(lengthRatio),
00043   theRatios(ratios),
00044   theIDMap(idMap),
00045   theDistAlgo(distAlgo),
00046   theDistCut(distCut)
00047 
00048 {
00049 
00050   gROOT->cd();
00051 
00052   std::string fullPath;
00053 
00054   // Prepare the map of files
00055   // Loop over the particle names
00056   std::vector<TFile*> aVFile(thePionEN.size(),static_cast<TFile*>(0));
00057   std::vector<TTree*> aVTree(thePionEN.size(),static_cast<TTree*>(0));
00058   std::vector<TBranch*> aVBranch(thePionEN.size(),static_cast<TBranch*>(0));
00059   std::vector<NUEvent*> aVNUEvents(thePionEN.size(),static_cast<NUEvent*>(0));
00060   std::vector<unsigned> aVCurrentEntry(thePionEN.size(),static_cast<unsigned>(0));
00061   std::vector<unsigned> aVCurrentInteraction(thePionEN.size(),static_cast<unsigned>(0));
00062   std::vector<unsigned> aVNumberOfEntries(thePionEN.size(),static_cast<unsigned>(0));
00063   std::vector<unsigned> aVNumberOfInteractions(thePionEN.size(),static_cast<unsigned>(0));
00064   std::vector<std::string> aVFileName(thePionEN.size(),static_cast<std::string>(""));
00065   std::vector<double> aVPionCM(thePionEN.size(),static_cast<double>(0));
00066   theFiles.resize(thePionNA.size());
00067   theTrees.resize(thePionNA.size());
00068   theBranches.resize(thePionNA.size());
00069   theNUEvents.resize(thePionNA.size());
00070   theCurrentEntry.resize(thePionNA.size());
00071   theCurrentInteraction.resize(thePionNA.size());
00072   theNumberOfEntries.resize(thePionNA.size());
00073   theNumberOfInteractions.resize(thePionNA.size());
00074   theFileNames.resize(thePionNA.size());
00075   thePionCM.resize(thePionNA.size());
00076   for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) { 
00077     theFiles[iname] = aVFile;
00078     theTrees[iname] = aVTree;
00079     theBranches[iname] = aVBranch;
00080     theNUEvents[iname] = aVNUEvents;
00081     theCurrentEntry[iname] = aVCurrentEntry;
00082     theCurrentInteraction[iname] = aVCurrentInteraction;
00083     theNumberOfEntries[iname] = aVNumberOfEntries;
00084     theNumberOfInteractions[iname] = aVNumberOfInteractions;
00085     theFileNames[iname] = aVFileName;
00086     thePionCM[iname] = aVPionCM;
00087   } 
00088 
00089   // Read the information from a previous run (to keep reproducibility)
00090   bool input = this->read(inputFile);
00091   if ( input ) 
00092     std::cout << "***WARNING*** You are reading nuclear-interaction information from the file "
00093               << inputFile << " created in an earlier run."
00094               << std::endl;
00095 
00096   // Open the file for saving the information of the current run
00097   myOutputFile.open ("NuclearInteractionOutputFile.txt");
00098   myOutputBuffer = 0;
00099 
00100 
00101   // Open the root files
00102   //  for ( unsigned file=0; file<theFileNames.size(); ++file ) {
00103   unsigned fileNb = 0;
00104   for ( unsigned iname=0; iname<thePionNA.size(); ++iname ) {
00105     for ( unsigned iene=0; iene<thePionEN.size(); ++iene ) {
00106       //std::cout << "iname/iene " << iname << " " << iene << std::endl; 
00107       std::ostringstream filename;
00108       double theEne = thePionEN[iene];
00109       filename << "NuclearInteractionsVal_" << thePionNA[iname] << "_E"<< theEne << ".root";
00110       theFileNames[iname][iene] = filename.str();
00111       //std::cout << "thePid/theEne " << thePionID[iname] << " " << theEne << std::endl; 
00112 
00113       edm::FileInPath myDataFile("FastSimulation/MaterialEffects/data/"+theFileNames[iname][iene]);
00114       fullPath = myDataFile.fullPath();
00115       //    theFiles[file] = TFile::Open(theFileNames[file].c_str());
00116       theFiles[iname][iene] = TFile::Open(fullPath.c_str());
00117       if ( !theFiles[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects") 
00118         << "File " << theFileNames[iname][iene] << " " << fullPath <<  " not found ";
00119       ++fileNb;
00120       //
00121       theTrees[iname][iene] = (TTree*) theFiles[iname][iene]->Get("NuclearInteractions"); 
00122       if ( !theTrees[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects") 
00123         << "Tree with name NuclearInteractions not found in " << theFileNames[iname][iene];
00124       //
00125       theBranches[iname][iene] = theTrees[iname][iene]->GetBranch("nuEvent");
00126       //std::cout << "The branch = " << theBranches[iname][iene] << std::endl;
00127       if ( !theBranches[iname][iene] ) throw cms::Exception("FastSimulation/MaterialEffects") 
00128         << "Branch with name nuEvent not found in " << theFileNames[iname][iene];
00129       //
00130       theNUEvents[iname][iene] = new NUEvent();
00131       //std::cout << "The branch = " << theBranches[iname][iene] << std::endl;
00132       theBranches[iname][iene]->SetAddress(&theNUEvents[iname][iene]);
00133       //
00134       theNumberOfEntries[iname][iene] = theTrees[iname][iene]->GetEntries();
00135 
00136       // Add some randomness (if there was no input file)
00137       if ( !input )
00138         theCurrentEntry[iname][iene] = (unsigned) (theNumberOfEntries[iname][iene] * random->flatShoot());
00139 
00140       theTrees[iname][iene]->GetEntry(theCurrentEntry[iname][iene]);
00141       unsigned NInteractions = theNUEvents[iname][iene]->nInteractions();
00142       theNumberOfInteractions[iname][iene] = NInteractions;
00143       // Add some randomness (if there was no input file)
00144       if ( !input )
00145         theCurrentInteraction[iname][iene] = (unsigned) (theNumberOfInteractions[iname][iene] * random->flatShoot());
00146 
00147       /*
00148       std::cout << "File " << theFileNames[iname][iene]
00149                 << " is opened with " << theNumberOfEntries[iname][iene] 
00150                 << " entries and will be read from Entry/Interaction "
00151                 << theCurrentEntry[iname][iene] << "/" << theCurrentInteraction[iname][iene] 
00152                 << std::endl;
00153       */
00154       
00155       //
00156       // Compute the corresponding cm energies of the nuclear interactions
00157       XYZTLorentzVector Proton(0.,0.,0.,0.986);
00158       XYZTLorentzVector 
00159         Reference(0.,
00160                   0.,
00161                   std::sqrt(thePionEN[iene]*thePionEN[iene]
00162                             -thePionMA[iname]*thePionMA[iname]),
00163                   thePionEN[iene]);
00164       thePionCM[iname][iene] = (Reference+Proton).M();
00165 
00166     }
00167 
00168   }
00169 
00170   // Find the index for which EN = 4. (or thereabout)
00171   ien4 = 0;
00172   while ( thePionEN[ien4] < 4.0 ) ++ien4;
00173 
00174   // Return Loot in the same state as it was when entering. 
00175   gROOT->cd();
00176 
00177   // Information (Should be on LogInfo)
00178 //  std::cout << " ---> A total of " << fileNb 
00179 //          << " nuclear-interaction files was sucessfully open" << std::endl;
00180 
00181   //  dbe = edm::Service<DaqMonitorBEInterface>().operator->();
00182   //  htot = dbe->book1D("Total", "All particles",150,0.,150.);
00183   //  helas = dbe->book1D("Elastic", "Elastic interactions",150,0.,150.);
00184   //  hinel = dbe->book1D("Inelastic", "Inelastic interactions",150,0.,150.);
00185   //  hscatter = dbe->book1D("Scattering","Elastic Scattering angle",200,0.,2.); 
00186   //  hscatter2 = dbe->book2D("Scattering2","Elastic Scattering angle vs p",100,0.,10.,200,0.,2.); 
00187   //  hAfter = dbe->book1D("eAfter","Energy after collision",200,0.,4.); 
00188   //  hAfter2 = dbe->book2D("eAfter2","Energy after collision",100,-2.5,2.5,100,0.,4); 
00189   //  hAfter3 = dbe->book2D("eAfter3","Energy after collision",100,0.,1000.,100,0.,4); 
00190 
00191 }
00192 
00193 NuclearInteractionSimulator::~NuclearInteractionSimulator() {
00194 
00195   // Close all local files
00196   // Among other things, this allows the TROOT destructor to end up 
00197   // without crashing, while trying to close these files from outside
00198   for ( unsigned ifile=0; ifile<theFiles.size(); ++ifile ) { 
00199     for ( unsigned iene=0; iene<theFiles[ifile].size(); ++iene ) {
00200       // std::cout << "Closing file " << iene << " with name " << theFileNames[ifile][iene] << std::endl;
00201       theFiles[ifile][iene]->Close(); 
00202     }
00203   }
00204 
00205   // Close the output file
00206   myOutputFile.close();
00207 
00208   // And return Loot in the same state as it was when entering. 
00209   gROOT->cd();
00210 
00211   //  dbe->save("test.root");
00212 
00213 }
00214 
00215 void NuclearInteractionSimulator::compute(ParticlePropagator& Particle)
00216 {
00217 
00218   // Read a Nuclear Interaction in a random manner
00219 
00220   double pHadron = std::sqrt(Particle.Vect().Mag2()); 
00221   //  htot->Fill(pHadron);
00222 
00223   // The hadron has enough momentum to create some relevant final state
00224   if ( pHadron > thePionEnergy ) { 
00225 
00226     // The particle type
00227     std::map<int,int>::const_iterator thePit = theIDMap.find(Particle.pid());
00228     
00229     int thePid = thePit != theIDMap.end() ? thePit->second : Particle.pid(); 
00230 
00231     // Is this particle type foreseen?
00232     unsigned fPid = abs(thePid);
00233     if ( fPid != 211 && fPid != 130 && fPid != 321 && fPid != 2112 && fPid != 2212 ) { 
00234       return;
00235       //std::cout << "Unknown particle type = " << thePid << std::endl;
00236       //thePid = 211;
00237     }
00238 
00239     // The inelastic interaction length at p(pion) = 5 GeV/c
00240     unsigned thePidIndex = index(thePid);
00241     double theInelasticLength = radLengths * theLengthRatio[thePidIndex];
00242 
00243     // The elastic interaction length
00244     // The baroque parameterization is a fit to Fig. 40.13 of the PDG
00245     double ee = pHadron > 0.6 ? 
00246       exp(-std::sqrt((pHadron-0.6)/1.122)) : exp(std::sqrt((0.6-pHadron)/1.122));
00247     double theElasticLength = ( 0.8753 * ee + 0.15 )
00248     //    double theElasticLength = ( 0.15 + 0.195 / log(pHadron/0.4) )
00249     //    double theElasticLength = ( 0.15 + 0.305 / log(pHadron/0.35) )
00250                             * theInelasticLength;
00251 
00252     // The total interaction length
00253     double theTotalInteractionLength = theInelasticLength + theElasticLength;
00254 
00255     // Probability to interact is dl/L0 (maximum for 4 GeV pion)
00256     double aNuclInteraction = -std::log(random->flatShoot());
00257     if ( aNuclInteraction < theTotalInteractionLength ) { 
00258 
00259       // The elastic part
00260       double elastic = random->flatShoot();
00261       if ( elastic < theElasticLength/theTotalInteractionLength ) {
00262 
00263         //              helas->Fill(pHadron);
00264         
00265         // Characteristic scattering angle for the elastic part
00266         double theta0 = std::sqrt(3.)/ std::pow(theA(),1./3.) * Particle.mass()/pHadron; 
00267         
00268         // Draw an angle with theta/theta0*exp[(-theta/2theta0)**2] shape 
00269         double theta = theta0 * std::sqrt(-2.*std::log(random->flatShoot()));
00270         double phi = 2. * 3.14159265358979323 * random->flatShoot();
00271         
00272         // Rotate the particle accordingly
00273         RawParticle::Rotation rotation1(orthogonal(Particle.Vect()),theta);
00274         RawParticle::Rotation rotation2(Particle.Vect(),phi);
00275         Particle.rotate(rotation1);
00276         Particle.rotate(rotation2);
00277         
00278         // Distance 
00279         double distance = std::sin(theta);
00280 
00281         // Create a daughter if the kink is large engough 
00282         if ( distance > theDistCut ) { 
00283           _theUpdatedState.resize(1);
00284           _theUpdatedState[0].SetXYZT(Particle.Px(), Particle.Py(),
00285                                       Particle.Pz(), Particle.E());
00286           _theUpdatedState[0].setID(Particle.pid());
00287         }
00288 
00289         //      hscatter->Fill(myTheta);
00290         //      hscatter2->Fill(pHadron,myTheta);
00291         
00292       } 
00293 
00294       // The inelastic part
00295       else {
00296     
00297         // Avoid multiple map access
00298         const std::vector<double>& aPionCM = thePionCM[thePidIndex];
00299         const std::vector<double>& aRatios = theRatios[thePidIndex];
00300         // Find the file with the closest c.m energy
00301         // The target nucleon
00302         XYZTLorentzVector Proton(0.,0.,0.,0.939);
00303         // The current particle
00304         const XYZTLorentzVector& Hadron = (const XYZTLorentzVector&)Particle;
00305         // The smallest momentum for inelastic interactions
00306         double pMin = thePionPMin[thePidIndex]; 
00307         // The correspong smallest four vector
00308         XYZTLorentzVector Hadron0(0.,0.,pMin,std::sqrt(pMin*pMin+Particle.M2()));
00309 
00310         // The current centre-of-mass energy
00311         double ecm = (Proton+Hadron).M();
00312         //std::cout << "Proton = " << Proton << std::endl;
00313         //std::cout << "Hadron = " << Hadron << std::endl;
00314         //std::cout << "ecm = " << ecm << std::endl;
00315         // Get the files of interest (closest c.m. energies)
00316         unsigned ene1=0;
00317         unsigned ene2=0;
00318         // The smallest centre-of-mass energy
00319         //      double ecm1=1.63; 
00320         double ecm1= (Proton+Hadron0).M();
00321         //std::cout << "ecm1 = " << ecm1 << std::endl;
00322         //std::cout << "ecm[0] = " << aPionCM[0] << std::endl;
00323         //std::cout << "ecm[11] = " << aPionCM [ aPionCM.size()-1 ] << std::endl;
00324         double ecm2=aPionCM[0]; 
00325         double ratio1=0.;
00326         double ratio2=aRatios[0];
00327         if ( ecm > aPionCM[0] && ecm < aPionCM [ aPionCM.size()-1 ] ) {
00328           for ( unsigned ene=1; 
00329                 ene < aPionCM.size() && ecm > aPionCM[ene-1]; 
00330                 ++ene ) {
00331             if ( ecm<aPionCM[ene] ) { 
00332               ene2 = ene;
00333               ene1 = ene2-1;
00334               ecm1 = aPionCM[ene1];
00335               ecm2 = aPionCM[ene2];
00336               ratio1 = aRatios[ene1];
00337               ratio2 = aRatios[ene2];
00338             } 
00339           }
00340         } else if ( ecm > aPionCM[ aPionCM.size()-1 ] ) { 
00341           ene1 = aPionCM.size()-1;
00342           ene2 = aPionCM.size()-2;
00343           ecm1 = aPionCM[ene1];
00344           ecm2 = aPionCM[ene2];
00345           ratio1 = aRatios[ene2];
00346           ratio2 = aRatios[ene2];
00347         } 
00348 
00349         
00350         // The inelastic part of the cross section depends cm energy
00351         double slope = (std::log10(ecm )-std::log10(ecm1)) 
00352                      / (std::log10(ecm2)-std::log10(ecm1));
00353         double inelastic = ratio1 + (ratio2-ratio1) * slope;
00354         double inelastic4 =  pHadron < 4. ? aRatios[ien4] : 1.;  
00355 
00356         //std::cout << "Inelastic = " << ratio1 << " " << ratio2 << " " << inelastic << std::endl;
00357         //      std::cout << "Energy/Inelastic : " 
00358         //              << Hadron.e() << " " << inelastic << std::endl;
00359         
00360         // Simulate an inelastic interaction
00361         if ( elastic > 1.- (inelastic*theInelasticLength)
00362                            /theTotalInteractionLength ) { 
00363 
00364           // Avoid mutliple map access
00365           std::vector<unsigned>& aCurrentInteraction = theCurrentInteraction[thePidIndex];
00366           std::vector<unsigned>& aNumberOfInteractions = theNumberOfInteractions[thePidIndex];
00367           std::vector<NUEvent*>& aNUEvents = theNUEvents[thePidIndex];
00368           //      hinel->Fill(pHadron);
00369           //      std::cout << "INELASTIC INTERACTION ! " 
00370           //                << pHadron << " " << theInelasticLength << " "
00371           //                << inelastic * theInelasticLength << std::endl;
00372           // Choice of the file to read according the the log10(ecm) distance
00373           // and protection against low momentum proton and neutron that never interacts 
00374           // (i.e., empty files)
00375           unsigned ene;
00376           if ( random->flatShoot() < slope || aNumberOfInteractions[ene1] == 0 )  
00377             ene = ene2;
00378           else
00379             ene = ene1;
00380 
00381           //std::cout << "Ecm1/2 = " << ecm1 << " " << ecm2 << std::endl;
00382           //std::cout << "Ratio1/2 = " << ratio1 << " " << ratio2 << std::endl;
00383           //std::cout << "Ene = " << ene << " slope = " << slope << std::endl;
00384 
00385           //std::cout << "Pion energy = " << Hadron.E() 
00386           //        << "File chosen " << theFileNames[thePidIndex][ene]
00387           //        << std::endl;
00388           
00389           // The boost characteristics
00390           XYZTLorentzVector theBoost = Proton + Hadron;
00391           theBoost /= theBoost.e();
00392           
00393           // std::cout << "File chosen : " << thePid << "/" << ene 
00394           //        << " Current interaction = " << aCurrentInteraction[ene] 
00395           //        << " Total interactions = " << aNumberOfInteractions[ene]
00396           //        << std::endl;
00397           //      theFiles[thePidIndex][ene]->cd();
00398           //      gDirectory->ls();
00399 
00400           // Check we are not either at the end of an interaction bunch 
00401           // or at the end of a file
00402           if ( aCurrentInteraction[ene] == aNumberOfInteractions[ene] ) {
00403             // std::cout << "End of interaction bunch ! ";
00404             std::vector<unsigned>& aCurrentEntry = theCurrentEntry[thePidIndex];
00405             std::vector<unsigned>& aNumberOfEntries = theNumberOfEntries[thePidIndex];
00406             std::vector<TTree*>& aTrees = theTrees[thePidIndex];
00407             ++aCurrentEntry[ene];
00408             // std::cerr << "Read the next entry " 
00409             //           << aCurrentEntry[ene] << std::endl;
00410             aCurrentInteraction[ene] = 0;
00411             if ( aCurrentEntry[ene] == aNumberOfEntries[ene] ) { 
00412               aCurrentEntry[ene] = 0;
00413               //  std::cout << "End of file - Rewind! " << std::endl;
00414             }
00415             unsigned myEntry = aCurrentEntry[ene];
00416             // std::cout << "The new entry " << myEntry 
00417             //           << " is read ... in TTree " << aTrees[ene] << " "; 
00418             aTrees[ene]->GetEntry(myEntry);
00419             // std::cout << "The number of interactions in the new entry is ... "; 
00420             aNumberOfInteractions[ene] = aNUEvents[ene]->nInteractions();
00421             // std::cout << aNumberOfInteractions[ene] << std::endl;
00422           }
00423           
00424           // Read the interaction
00425           NUEvent::NUInteraction anInteraction 
00426             = aNUEvents[ene]->theNUInteractions()[aCurrentInteraction[ene]];
00427           
00428           unsigned firstTrack = anInteraction.first; 
00429           unsigned lastTrack = anInteraction.last;
00430           //      std::cout << "First and last tracks are " << firstTrack << " " << lastTrack << std::endl;
00431           
00432           _theUpdatedState.resize(lastTrack-firstTrack+1);
00433 
00434           double distMin = 1E99;
00435 
00436           // Some rotation around the boost axis, for more randomness
00437           XYZVector theAxis = theBoost.Vect().Unit();
00438           double theAngle = random->flatShoot() * 2. * 3.14159265358979323;
00439           RawParticle::Rotation axisRotation(theAxis,theAngle);
00440           RawParticle::Boost axisBoost(theBoost.x(),theBoost.y(),theBoost.z());
00441 
00442           // A rotation to bring the particles back to the pion direction
00443           XYZVector zAxis(0.,0.,1.); 
00444           XYZVector orthAxis = (zAxis.Cross(theBoost.Vect())).Unit(); 
00445           double orthAngle = acos(theBoost.Vect().Unit().Z());
00446           RawParticle::Rotation orthRotation(orthAxis,orthAngle);
00447 
00448           // A few checks
00449           // double eAfter = 0.;
00450           
00451           // Loop on the nuclear interaction products
00452           for ( unsigned iTrack=firstTrack; iTrack<=lastTrack; ++iTrack ) {
00453             
00454             unsigned idaugh = iTrack - firstTrack;
00455             NUEvent::NUParticle aParticle = aNUEvents[ene]->theNUParticles()[iTrack];
00456             //  std::cout << "Track " << iTrack 
00457             //            << " id/px/py/pz/mass "
00458             //            << aParticle.id << " " 
00459             //            << aParticle.px << " " 
00460             //            << aParticle.py << " " 
00461             //            << aParticle.pz << " " 
00462             //            << aParticle.mass << " " << endl; 
00463             
00464             // Add a RawParticle with the proper energy in the c.m frame of 
00465             // the nuclear interaction
00466             double energy = std::sqrt( aParticle.px*aParticle.px
00467                                      + aParticle.py*aParticle.py
00468                                      + aParticle.pz*aParticle.pz
00469                                      + aParticle.mass*aParticle.mass/(ecm*ecm) );
00470 
00471             RawParticle& aDaughter = _theUpdatedState[idaugh]; 
00472             aDaughter.SetXYZT(aParticle.px*ecm,aParticle.py*ecm,
00473                               aParticle.pz*ecm,energy*ecm);         
00474             aDaughter.setID(aParticle.id);
00475 
00476             // Rotate to the collision axis
00477             aDaughter.rotate(orthRotation);
00478 
00479             // Rotate around the boost axis for more randomness
00480             aDaughter.rotate(axisRotation);
00481             
00482             // Boost it in the lab frame
00483             aDaughter.boost(axisBoost);
00484 
00485             // Store the closest daughter index (for later tracking purposes, so charged particles only) 
00486             double distance = distanceToPrimary(Particle,aDaughter);
00487             // Find the closest daughter, if closer than a given upper limit.
00488             if ( distance < distMin && distance < theDistCut ) {
00489               distMin = distance;
00490               theClosestChargedDaughterId = idaugh;
00491             }
00492 
00493             // eAfter += aDaughter.E();
00494 
00495           }
00496 
00497           /*
00498           double eBefore = Particle.E();
00499           double rapid = Particle.momentum().Eta();
00500           if ( eBefore > 0. ) { 
00501             hAfter->Fill(eAfter/eBefore);
00502             hAfter2->Fill(rapid,eAfter/eBefore);
00503             hAfter3->Fill(eBefore,eAfter/eBefore);
00504           }
00505           */
00506 
00507           // Increment for next time
00508           ++aCurrentInteraction[ene];
00509           
00510         // Simulate a stopping hadron (low momentum)
00511         } else if ( pHadron < 4. &&  
00512                     elastic > 1.- (inelastic4*theInelasticLength)
00513                                   /theTotalInteractionLength ) { 
00514           // A fake particle with 0 momentum as a daughter!
00515           _theUpdatedState.resize(1);
00516           _theUpdatedState[0].SetXYZT(0.,0.,0.,0.);
00517           _theUpdatedState[0].setID(22);
00518         }
00519 
00520       }
00521 
00522     }
00523 
00524   }
00525 
00526 }
00527 
00528 double 
00529 NuclearInteractionSimulator::distanceToPrimary(const RawParticle& Particle,
00530                                                const RawParticle& aDaughter) const {
00531 
00532   double distance = 2E99;
00533 
00534   // Compute the distance only for charged primaries
00535   if ( fabs(Particle.charge()) > 1E-12 ) { 
00536 
00537     // The secondary must have the same charge
00538     double chargeDiff = fabs(aDaughter.charge()-Particle.charge());
00539     if ( fabs(chargeDiff) < 1E-12 ) {
00540 
00541       // Here are two distance definitions * to be tuned *
00542       switch ( theDistAlgo ) { 
00543         
00544       case 1:
00545         // sin(theta12)
00546         distance = (aDaughter.Vect().Unit().Cross(Particle.Vect().Unit())).R();
00547         break;
00548         
00549       case 2: 
00550         // sin(theta12) * p1/p2
00551         distance = (aDaughter.Vect().Cross(Particle.Vect())).R()
00552           /aDaughter.Vect().Mag2();
00553         break;
00554         
00555       default:
00556         // Should not happen
00557         distance = 2E99;
00558         break;
00559         
00560       }
00561 
00562     }
00563 
00564   }
00565       
00566   return distance;
00567 
00568 }
00569 
00570 void
00571 NuclearInteractionSimulator::save() {
00572 
00573   // Size of buffer
00574   ++myOutputBuffer;
00575 
00576   // Periodically close the current file and open a new one
00577   if ( myOutputBuffer/1000*1000 == myOutputBuffer ) { 
00578     myOutputFile.close();
00579     myOutputFile.open ("NuclearInteractionOutputFile.txt");
00580     //    myOutputFile.seekp(0); // No need to rewind in that case
00581   }
00582   //
00583   unsigned size1 = 
00584     theCurrentEntry.size()*
00585     theCurrentEntry.begin()->size();
00586   std::vector<unsigned> theCurrentEntries;
00587   theCurrentEntries.resize(size1);
00588   size1*=sizeof(unsigned);
00589   //
00590   unsigned size2 = 
00591     theCurrentInteraction.size()*
00592     theCurrentInteraction.begin()->size();
00593   std::vector<unsigned> theCurrentInteractions;
00594   theCurrentInteractions.resize(size2);
00595   size2 *= sizeof(unsigned);
00596 
00597   // Save the current entries 
00598   std::vector< std::vector<unsigned> >::const_iterator aCurrentEntry = theCurrentEntry.begin();
00599   std::vector< std::vector<unsigned> >::const_iterator lastCurrentEntry = theCurrentEntry.end();
00600   unsigned allEntries=0;
00601   for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) { 
00602     unsigned size = aCurrentEntry->size();
00603     for ( unsigned iene=0; iene<size; ++iene )
00604       theCurrentEntries[allEntries++] = (*aCurrentEntry)[iene]; 
00605   }
00606   
00607   // Save the current interactions
00608   std::vector< std::vector<unsigned> >::const_iterator aCurrentInteraction = theCurrentInteraction.begin();
00609   std::vector< std::vector<unsigned> >::const_iterator lastCurrentInteraction = theCurrentInteraction.end();
00610   unsigned allInteractions=0;
00611   for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) { 
00612     unsigned size = aCurrentInteraction->size();
00613     for ( unsigned iene=0; iene<size; ++iene )
00614       theCurrentInteractions[allInteractions++] = (*aCurrentInteraction)[iene]; 
00615   }
00616   // 
00617   myOutputFile.write((const char*)(&theCurrentEntries.front()), size1);
00618   myOutputFile.write((const char*)(&theCurrentInteractions.front()), size2);  
00619   myOutputFile.flush();
00620 
00621 }
00622 
00623 bool
00624 NuclearInteractionSimulator::read(std::string inputFile) {
00625 
00626   ifstream myInputFile;
00627   struct stat results;
00628   //
00629   unsigned size1 = 
00630     theCurrentEntry.size()*
00631     theCurrentEntry.begin()->size();
00632   std::vector<unsigned> theCurrentEntries;
00633   theCurrentEntries.resize(size1);
00634   size1*=sizeof(unsigned);
00635   //
00636   unsigned size2 = 
00637     theCurrentInteraction.size()*
00638     theCurrentInteraction.begin()->size();
00639   std::vector<unsigned> theCurrentInteractions;
00640   theCurrentInteractions.resize(size2);
00641   size2 *= sizeof(unsigned);
00642   //
00643   unsigned size = 0;
00644 
00645 
00646   // Open the file (if any)
00647   myInputFile.open (inputFile.c_str());
00648   if ( myInputFile.is_open() ) { 
00649 
00650     // Get the size of the file
00651     if ( stat(inputFile.c_str(), &results) == 0 ) size = results.st_size;
00652     else return false; // Something is wrong with that file !
00653 
00654     // Position the pointer just before the last record
00655     myInputFile.seekg(size-size1-size2);
00656     myInputFile.read((char*)(&theCurrentEntries.front()),size1);
00657     myInputFile.read((char*)(&theCurrentInteractions.front()),size2);
00658     myInputFile.close();
00659 
00660     // Read the current entries
00661     std::vector< std::vector<unsigned> >::iterator aCurrentEntry = theCurrentEntry.begin();
00662     std::vector< std::vector<unsigned> >::iterator lastCurrentEntry = theCurrentEntry.end();
00663     unsigned allEntries=0;
00664     for ( ; aCurrentEntry!=lastCurrentEntry; ++aCurrentEntry ) { 
00665       unsigned size = aCurrentEntry->size();
00666       for ( unsigned iene=0; iene<size; ++iene )
00667         (*aCurrentEntry)[iene] = theCurrentEntries[allEntries++]; 
00668     }
00669 
00670     // Read the current interactions
00671     std::vector< std::vector<unsigned> >::iterator aCurrentInteraction = theCurrentInteraction.begin();
00672     std::vector< std::vector<unsigned> >::iterator lastCurrentInteraction = theCurrentInteraction.end();
00673     unsigned allInteractions=0;
00674     for ( ; aCurrentInteraction!=lastCurrentInteraction; ++aCurrentInteraction ) { 
00675       unsigned size = aCurrentInteraction->size();
00676       for ( unsigned iene=0; iene<size; ++iene )
00677         (*aCurrentInteraction)[iene] =  theCurrentInteractions[allInteractions++];
00678     }
00679     
00680     return true;
00681   }
00682  
00683   return false;
00684 
00685 }
00686 
00687 unsigned 
00688 NuclearInteractionSimulator::index(int thePid) { 
00689   
00690   unsigned myIndex=0;
00691   while ( thePid != thePionID[myIndex] ) ++myIndex; 
00692   //  std::cout << "pid/index = " << thePid << " " << myIndex << std::endl;
00693   return myIndex;
00694 
00695 }