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