HepMC Reference Documentation

HepMC

testHepMCIteration.cc.in

Use Matt's example_EventSelection along with example_UsingIterators to check HepMC iteration. Apply an event selection to the events in testHepMC.input Events containing a photon of pT > 25 GeV pass the selection. Use iterators on these events.

00001 
00002 // testHepMCIteration.cc.in
00003 //
00004 // garren@fnal.gov, May 2007
00005 // Use Matt's example_EventSelection along with example_UsingIterators 
00006 // to check HepMC iteration.
00007 // Apply an event selection to the events in testHepMC.input
00008 // Events containing a photon of pT > 25 GeV pass the selection.
00009 // Use iterators on these events.
00011 
00012 #include <list>
00013 
00014 #include "HepMC/IO_GenEvent.h"
00015 #include "HepMC/IO_AsciiParticles.h"
00016 #include "HepMC/GenEvent.h"
00017 
00018 // define methods and classes used by this test
00019 #include "IsGoodEvent.h"
00020 #include "testHepMCIteration.h"
00021 
00022 bool findW( HepMC::GenEvent* evt, std::ofstream& os);
00023 bool simpleIter( HepMC::GenEvent* evt, std::ofstream& os );
00024 
00025 int main() { 
00026     // declare an input strategy to read the data produced with the 
00027     // example_MyPythia
00028     HepMC::IO_GenEvent ascii_in("@srcdir@/testIOGenEvent.input",std::ios::in);
00029     // declare an instance of the event selection predicate
00030     IsGoodEvent is_good_event;
00031     // define an output stream
00032     std::ofstream os( "testHepMCIteration.out" );
00033     //........................................EVENT LOOP
00034     int icount=0;
00035     int num_good_events=0;
00036     HepMC::GenEvent* evt = ascii_in.read_next_event();
00037     HepMC::GenEvent* evcopy;
00038     while ( evt ) {
00039         icount++;
00040         if ( icount%50==1 ) std::cout << "Processing Event Number " << icount
00041                                       << " its # " << evt->event_number() 
00042                                       << std::endl;
00043         // icount of 100 should be the last event
00044         if ( icount==100 ) std::cout << "Processing Event Number " << icount
00045                                       << " its # " << evt->event_number() 
00046                                       << std::endl;
00047         evcopy = evt;
00048         if ( is_good_event(evcopy) ) {
00049             os << "Event " << evcopy->event_number() << " is good " << std::endl;
00050             ++num_good_events;
00051             // test iterators
00052             simpleIter( evcopy, os );
00053             findW( evcopy, os );
00054         }
00055         evcopy->clear();
00056         
00057         // clean up and get next event
00058         delete evt;
00059         evt = ascii_in.read_next_event();
00060     }
00061     //........................................PRINT RESULT
00062     std::cout << num_good_events << " out of " << icount 
00063               << " processed events passed the cuts. Finished." << std::endl;
00064 }
00065 
00066 bool simpleIter( HepMC::GenEvent* evt, std::ofstream& os )
00067 {
00068         // use GenEvent::vertex_iterator to fill a list of all 
00069         // vertices in the event
00070         std::list<HepMC::GenVertex*> allvertices;
00071         for ( HepMC::GenEvent::vertex_iterator v = evt->vertices_begin();
00072               v != evt->vertices_end(); ++v ) {
00073             allvertices.push_back(*v);
00074         }
00075 
00076         // we could do the same thing with the STL algorithm copy
00077         std::list<HepMC::GenVertex*> allvertices2;
00078         copy( evt->vertices_begin(), evt->vertices_end(), 
00079               back_inserter(allvertices2) );
00080 
00081         // fill a list of all final state particles in the event, by requiring
00082         // that each particle satisfyies the IsFinalState predicate
00083         IsFinalState isfinal;
00084         std::list<HepMC::GenParticle*> finalstateparticles;
00085         for ( HepMC::GenEvent::particle_iterator p = evt->particles_begin();
00086               p != evt->particles_end(); ++p ) {
00087             if ( isfinal(*p) ) finalstateparticles.push_back(*p);
00088         }
00089 
00090         // an STL-like algorithm called HepMC::copy_if is provided in the
00091         // GenEvent.h header to do this sort of operation more easily,
00092         // you could get the identical results as above by using:
00093         std::list<HepMC::GenParticle*> finalstateparticles2;
00094         HepMC::copy_if( evt->particles_begin(), evt->particles_end(), 
00095                         back_inserter(finalstateparticles2), IsFinalState() );
00096 
00097         // print all photons in the event that satisfy the IsPhoton criteria
00098         os << "photons in event " << evt->event_number() << ":" << std::endl;
00099         for ( HepMC::GenEvent::particle_iterator p = evt->particles_begin();
00100               p != evt->particles_end(); ++p ) {
00101             if ( IsPhoton(*p) ) (*p)->print( os );
00102         }
00103         return true;
00104 }
00105 
00106 bool findW( HepMC::GenEvent* evt, std::ofstream& os )
00107 {
00108     int num_W=0;
00109         // use GenEvent::particle_iterator to find all W's in the event,
00110         // then 
00111         // (1) for each W user the GenVertex::particle_iterator with a range of
00112         //     parents to return and print the immediate mothers of these W's.
00113         // (2) for each W user the GenVertex::particle_iterator with a range of
00114         //     descendants to return and print all descendants of these W's.
00115         for ( HepMC::GenEvent::particle_iterator p = evt->particles_begin();
00116               p != evt->particles_end(); ++p ) {
00117             if ( IsWBoson(*p) ) {
00118                 ++num_W;
00119                 os << "A W boson has been found in event: " << evt->event_number() << std::endl;
00120                 (*p)->print( os );
00121                 // return all parents
00122                 // we do this by pointing to the production vertex of the W 
00123                 // particle and asking for all particle parents of that vertex
00124                 os << "\t Its parents are: " << std::endl;
00125                 if ( (*p)->production_vertex() ) {
00126                     for ( HepMC::GenVertex::particle_iterator mother 
00127                               = (*p)->production_vertex()->
00128                               particles_begin(HepMC::parents);
00129                           mother != (*p)->production_vertex()->
00130                               particles_end(HepMC::parents); 
00131                           ++mother ) {
00132                         os << "\t";
00133                         (*mother)->print( os );
00134                     }
00135                 }
00136 
00137                 // return immediate children
00138                 os << "\t\t" << "Its children are: " << std::endl;
00139                 if ( (*p)->end_vertex() ) {
00140                     for ( HepMC::GenVertex::particle_iterator child = 
00141                           (*p)->end_vertex()->particles_begin(HepMC::children); 
00142                           child != (*p)->end_vertex()->particles_end(HepMC::children); 
00143                           ++child ) {
00144                         // make a copy
00145                         HepMC::GenVertex::particle_iterator cp = child;
00146                         // use the copy and the original
00147                         os << "\t\t\t (id,barcode,status) " 
00148                            << (*cp)->pdg_id() << " " 
00149                            << (*child)->barcode() << " "
00150                            << (*cp)->status() << std::endl;
00151                     }
00152                 }
00153 
00154                 // return all descendants
00155                 // we do this by pointing to the end vertex of the W 
00156                 // particle and asking for all particle descendants of that vertex
00157                 os << "\t\t Its descendants are: " << std::endl;
00158                 if ( (*p)->end_vertex() ) {
00159                     for ( HepMC::GenVertex::particle_iterator des 
00160                               =(*p)->end_vertex()->
00161                               particles_begin(HepMC::descendants);
00162                           des != (*p)->end_vertex()->
00163                               particles_end(HepMC::descendants);
00164                           ++des ) {
00165                         os << "\t\t";
00166                         (*des)->print( os );
00167                     }
00168                 }
00169             }   // if IsWBoson
00170         }       // end particle loop
00171         return true;
00172 }
00173 

Generated on Thu Jan 7 13:10:15 2010 for HepMC by  doxygen 1.4.7