CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/IOMC/Input/src/HepMCFileReader.cc

Go to the documentation of this file.
00001 // $Id: HepMCFileReader.cc,v 1.11 2009/12/01 19:23:11 fabstoec Exp $
00002 
00012 #include <iostream>
00013 #include <iomanip>
00014 #include <string>
00015 
00016 #include "HepMC/GenEvent.h"
00017 #include "HepMC/GenParticle.h"
00018 #include "HepMC/IO_GenEvent.h"
00019 
00020 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00021 
00022 #include "IOMC/Input/interface/HepMCFileReader.h"
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 #include "FWCore/Utilities/interface/Exception.h"
00025 
00026 using namespace std;
00027 
00028 
00029 HepMCFileReader *HepMCFileReader::instance_=0;
00030 
00031 
00032 //-------------------------------------------------------------------------
00033 HepMCFileReader *HepMCFileReader::instance()
00034 {
00035   // Implement a HepMCFileReader singleton.
00036 
00037   if (instance_== 0) {
00038     instance_ = new HepMCFileReader();
00039   }
00040   return instance_;
00041 }
00042 
00043 
00044 //-------------------------------------------------------------------------
00045 HepMCFileReader::HepMCFileReader() :
00046   evt_(0), input_(0)
00047 { 
00048   // Default constructor.
00049   if (instance_ == 0) {
00050     instance_ = this;  
00051   } else {
00052     edm::LogError("HepMCFileReader") << "Constructing a new instance";
00053   }
00054 } 
00055 
00056 
00057 //-------------------------------------------------------------------------
00058 HepMCFileReader::~HepMCFileReader()
00059 {
00060   edm::LogInfo("HepMCFileReader") << "Destructing HepMCFileReader";
00061   
00062   instance_=0;    
00063   delete input_;
00064 }
00065 
00066 
00067 //-------------------------------------------------------------------------
00068 void HepMCFileReader::initialize(const string &filename)
00069 {
00070   if (isInitialized()) {
00071     edm::LogError("HepMCFileReader") << "Was already initialized... reinitializing";
00072     delete input_;
00073   }
00074 
00075   edm::LogInfo("HepMCFileReader") << "Opening file" << filename << "using HepMC::IO_GenEvent";
00076   input_ = new HepMC::IO_GenEvent(filename.c_str(), std::ios::in);
00077 
00078   if (rdstate() == std::ios::failbit) {
00079     throw cms::Exception("FileNotFound", "HepMCFileReader::initialize()")
00080       << "File " << filename << " was not found.\n";
00081   }
00082 }
00083 
00084 
00085 //-------------------------------------------------------------------------
00086 int HepMCFileReader::rdstate() const
00087 {
00088   // work around a HepMC IO_ inheritence shortfall
00089 
00090   HepMC::IO_GenEvent *p = dynamic_cast<HepMC::IO_GenEvent*>(input_);
00091   if (p) return p->rdstate();
00092 
00093   return std::ios::failbit;
00094 }
00095 
00096 
00097 //-------------------------------------------------------------------------
00098 bool HepMCFileReader::readCurrentEvent()
00099 {
00100   evt_ = input_->read_next_event();
00101   if (evt_) { 
00102     edm::LogInfo("HepMCFileReader") << "| --- Event Nr. "  << evt_->event_number()
00103          << " with " << evt_->particles_size() << " particles --- !";  
00104     ReadStats();
00105     //  printHepMcEvent();
00106     //          printEvent();
00107   } else {
00108     edm::LogInfo("HepMCFileReader") << "Got no event" <<endl;
00109   }
00110 
00111   return evt_ != 0;
00112 }
00113 
00114 
00115 //-------------------------------------------------------------------------
00116 bool HepMCFileReader::setEvent(int event)
00117 {
00118   return true;
00119 }
00120 
00121 
00122 //-------------------------------------------------------------------------
00123 bool HepMCFileReader::printHepMcEvent() const
00124 {
00125   if (evt_ != 0) evt_->print(); 
00126   return true;
00127 }
00128 
00129 
00130 //-------------------------------------------------------------------------
00131 HepMC::GenEvent * HepMCFileReader::fillCurrentEventData()
00132 {
00133   readCurrentEvent();
00134   return evt_;
00135 }
00136 
00137 
00138 //-------------------------------------------------------------------------
00139 // Print out in old CMKIN style for comparisons
00140 void HepMCFileReader::printEvent() const {
00141   int mo1=0,mo2=0,da1=0,da2=0,status=0,pid=0;   
00142   if (evt_ != 0) {   
00143     cout << "---#-------pid--st---Mo1---Mo2---Da1---Da2------px------py------pz-------E-";
00144     cout << "------m---------x---------y---------z---------t-";
00145     cout << endl;
00146     cout.setf(ios::right, ios::adjustfield);
00147     for(int n=1; n<=evt_->particles_size(); n++) {  
00148       HepMC::GenParticle *g = index_to_particle[n];  
00149       getStatsFromTuple( mo1,mo2,da1,da2,status,pid,n);
00150       cout << setw(4) << n
00151       << setw(8) << pid
00152       << setw(5) << status
00153       << setw(6) << mo1    
00154       << setw(6) << mo2
00155       << setw(6) << da1
00156       << setw(6) << da2;
00157       cout.setf(ios::fixed, ios::floatfield);
00158       cout.setf(ios::right, ios::adjustfield);
00159       cout << setw(10) << setprecision(2) << g->momentum().x();
00160       cout << setw(8) << setprecision(2) << g->momentum().y();
00161       cout << setw(10) << setprecision(2) << g->momentum().z();
00162       cout << setw(8) << setprecision(2) << g->momentum().t();
00163       cout << setw(8) << setprecision(2) << g->generatedMass();       
00164       // tau=L/(gamma*beta*c) 
00165       if (g->production_vertex() != 0 && g->end_vertex() != 0 && status == 2) {
00166         cout << setw(10) << setprecision(2) <<g->production_vertex()->position().x();
00167         cout << setw(10) << setprecision(2) <<g->production_vertex()->position().y();
00168         cout << setw(10) << setprecision(2) <<g->production_vertex()->position().z();
00169         
00170         double xm = g->production_vertex()->position().x();
00171         double ym = g->production_vertex()->position().y();
00172         double zm = g->production_vertex()->position().z();
00173         double xd = g->end_vertex()->position().x();
00174         double yd = g->end_vertex()->position().y();
00175         double zd = g->end_vertex()->position().z();
00176         double decl = sqrt((xd-xm)*(xd-xm)+(yd-ym)*(yd-ym)+(zd-zm)*(zd-zm));
00177         double labTime = decl/c_light;
00178         // convert lab time to proper time
00179         double properTime =  labTime/g->momentum().rho()*(g->generatedMass() );
00180         // set the proper time in nanoseconds
00181         cout << setw(8) << setprecision(2) <<properTime;
00182       }
00183       else{
00184         cout << setw(10) << setprecision(2) << 0.0;
00185         cout << setw(10) << setprecision(2) << 0.0;
00186         cout << setw(10) << setprecision(2) << 0.0;
00187         cout << setw(8) << setprecision(2) << 0.0;  
00188       }
00189       cout << endl; 
00190     }
00191   } else {
00192     cout << " HepMCFileReader: No event available !" << endl;
00193   }
00194 }
00195 
00196 
00197 //-------------------------------------------------------------------------
00198 void HepMCFileReader::ReadStats()
00199 {
00200   unsigned int particle_counter=0;
00201   index_to_particle.reserve(evt_->particles_size()+1); 
00202   index_to_particle[0] = 0; 
00203   for (HepMC::GenEvent::vertex_const_iterator v = evt_->vertices_begin();
00204        v != evt_->vertices_end(); ++v ) {
00205     // making a list of incoming particles of the vertices
00206     // so that the mother indices in HEPEVT can be filled properly
00207     for (HepMC::GenVertex::particles_in_const_iterator p1 = (*v)->particles_in_const_begin();
00208          p1 != (*v)->particles_in_const_end(); ++p1 ) {
00209       ++particle_counter;
00210       //particle_counter can be very large for heavy ions
00211       if(particle_counter >= index_to_particle.size() ) {
00212         //make it large enough to hold up to this index
00213          index_to_particle.resize(particle_counter+1);
00214       }                 
00215       index_to_particle[particle_counter] = *p1;
00216       particle_to_index[*p1] = particle_counter;
00217     }   
00218     // daughters are entered only if they aren't a mother of
00219     // another vertex
00220     for (HepMC::GenVertex::particles_out_const_iterator p2 = (*v)->particles_out_const_begin();
00221          p2 != (*v)->particles_out_const_end(); ++p2) {
00222       if (!(*p2)->end_vertex()) {
00223         ++particle_counter;
00224         //particle_counter can be very large for heavy ions
00225         if(particle_counter  >= index_to_particle.size() ) {        
00226           //make it large enough to hold up to this index
00227           index_to_particle.resize(particle_counter+1);
00228         }                   
00229         index_to_particle[particle_counter] = *p2;
00230         particle_to_index[*p2] = particle_counter;
00231       }
00232     }
00233   } 
00234 }
00235 
00236 
00237 //-------------------------------------------------------------------------
00238 void HepMCFileReader::getStatsFromTuple(int &mo1, int &mo2, int &da1, int &da2,
00239                                         int &status, int &pid, int j) const
00240 {
00241   if (!evt_) {
00242     cout <<  "HepMCFileReader: Got no event :-(  Game over already  ?" <<endl; 
00243   } else {
00244     status =  index_to_particle[j]->status();
00245     pid = index_to_particle[j]->pdg_id();
00246     if ( index_to_particle[j]->production_vertex() ) {
00247       
00248           //HepLorentzVector p = index_to_particle[j]->
00249           //production_vertex()->position();
00250       
00251       int num_mothers = index_to_particle[j]->production_vertex()->
00252       particles_in_size();
00253       int first_mother = find_in_map( particle_to_index,
00254       *(index_to_particle[j]->
00255       production_vertex()->
00256       particles_in_const_begin()));
00257       int last_mother = first_mother + num_mothers - 1;
00258       if ( first_mother == 0 ) last_mother = 0;
00259       mo1=first_mother;
00260       mo2=last_mother;
00261     } else {
00262       mo1 =0;
00263       mo2 =0;
00264     }   
00265     if (!index_to_particle[j]->end_vertex()) {
00266       //find # of 1. daughter
00267       int first_daughter = find_in_map( particle_to_index,
00268       *(index_to_particle[j]->end_vertex()->particles_begin(HepMC::children)));
00269       //cout <<"first_daughter "<< first_daughter <<  "num_daughters " << num_daughters << endl; 
00270       HepMC::GenVertex::particle_iterator ic;
00271       int last_daughter=0;
00272       //find # of last daughter
00273       for (ic = index_to_particle[j]->end_vertex()->particles_begin(HepMC::children);
00274       ic != index_to_particle[j]->end_vertex()->particles_end(HepMC::children); ++ic) 
00275       last_daughter= find_in_map( particle_to_index,*ic);       
00276       
00277       if (first_daughter== 0) last_daughter = 0;
00278       da1=first_daughter;
00279       da2=last_daughter;
00280     } else {
00281       da1=0;
00282       da2=0;
00283     }
00284   } 
00285 } 
00286 
00287 
00288 //-------------------------------------------------------------------------
00289 int HepMCFileReader::find_in_map( const std::map<HepMC::GenParticle*,int>& m, 
00290                                   HepMC::GenParticle *p) const
00291 {
00292   std::map<HepMC::GenParticle*,int>::const_iterator iter = m.find(p);
00293   return (iter == m.end()) ? 0 : iter->second;
00294 }
00295