HepMC Reference Documentation

HepMC

example_MyPythia.cc

Go to the documentation of this file.
00001 
00002 // Matt.Dobbs@Cern.CH, December 1999
00003 // November 2000, updated to use Pythia 6.1
00004 // 
00046 
00047 
00048 #include <iostream>
00049 #include "HepMC/PythiaWrapper.h"
00050 #include "HepMC/IO_HEPEVT.h"
00051 #include "HepMC/IO_GenEvent.h"
00052 #include "HepMC/IO_AsciiParticles.h"
00053 #include "HepMC/GenEvent.h"
00054 #include "PythiaHelper.h"
00055 
00057 
00061 class IsGoodEventMyPythia {
00062 public:
00064     bool operator()( const HepMC::GenEvent* evt ) { 
00065         for ( HepMC::GenEvent::particle_const_iterator p 
00066                   = evt->particles_begin(); p != evt->particles_end(); ++p ){
00067             if ( (*p)->pdg_id() == 22 && (*p)->momentum().perp() > 25. ) {
00068                 //std::cout << "Event " << evt->event_number()
00069                 //     << " is a good event." << std::endl;
00070                 //(*p)->print();
00071                 return 1;
00072             }
00073         }
00074         return 0;
00075     }
00076 };
00077     
00078 
00079 void pythia_out();
00080 void pythia_in();
00081 void pythia_in_out();
00082 void event_selection();
00083 void pythia_particle_out();
00084 
00085 int main() { 
00086     // example to generate events and write output
00087     pythia_out();
00088     // example to generate events and perform simple event selection
00089     event_selection();
00090     // example to read the file written by pythia_out
00091     pythia_in();
00092     // example to generate events, write them, and read them back
00093     pythia_in_out();
00094 
00095     return 0;
00096 }
00097 
00098 
00099 void pythia_out()
00100 {
00101     std::cout << std::endl;
00102     std::cout << "Begin pythia_out()" << std::endl;
00103    //........................................HEPEVT
00104     // Pythia 6.1 uses HEPEVT with 4000 entries and 8-byte floating point
00105     //  numbers. We need to explicitly pass this information to the 
00106     //  HEPEVT_Wrapper.
00107     //
00108     HepMC::HEPEVT_Wrapper::set_max_number_entries(4000);
00109     HepMC::HEPEVT_Wrapper::set_sizeof_real(8);
00110     //
00111     //........................................PYTHIA INITIALIZATIONS
00112     initPythia();
00113 
00114     //........................................HepMC INITIALIZATIONS
00115     //
00116     // Instantiate an IO strategy for reading from HEPEVT.
00117     HepMC::IO_HEPEVT hepevtio;
00118     //
00119     { // begin scope of ascii_io
00120         // Instantiate an IO strategy to write the data to file 
00121         HepMC::IO_GenEvent ascii_io("example_MyPythia.dat",std::ios::out);
00122         //
00123         //........................................EVENT LOOP
00124         for ( int i = 1; i <= 100; i++ ) {
00125             if ( i%50==1 ) std::cout << "Processing Event Number " 
00126                                      << i << std::endl;
00127             call_pyevnt();      // generate one event with Pythia
00128             // pythia pyhepc routine converts common PYJETS in common HEPEVT
00129             call_pyhepc( 1 );
00130             HepMC::GenEvent* evt = hepevtio.read_next_event();
00131             // define the units (Pythia uses GeV and mm)
00132             evt->use_units(HepMC::Units::GEV, HepMC::Units::MM);
00133             // add some information to the event
00134             evt->set_event_number(i);
00135             evt->set_signal_process_id(20);
00136             // set number of multi parton interactions
00137             evt->set_mpi( pypars.msti[31-1] );
00138             // set cross section information
00139             evt->set_cross_section( getPythiaCrossSection() );
00140             // write the event out to the ascii files
00141             ascii_io << evt;
00142             // we also need to delete the created event from memory
00143             delete evt;
00144         }
00145         //........................................TERMINATION
00146         // write out some information from Pythia to the screen
00147         call_pystat( 1 );    
00148     } // end scope of ascii_io
00149 }
00150 
00151  
00152 void event_selection()
00153 {
00154     std::cout << std::endl;
00155     std::cout << "Begin event_selection()" << std::endl;
00156     //........................................HEPEVT
00157     // Pythia 6.1 uses HEPEVT with 4000 entries and 8-byte floating point
00158     //  numbers. We need to explicitly pass this information to the 
00159     //  HEPEVT_Wrapper.
00160     //
00161     HepMC::HEPEVT_Wrapper::set_max_number_entries(4000);
00162     HepMC::HEPEVT_Wrapper::set_sizeof_real(8);
00163     //  
00164     //........................................PYTHIA INITIALIZATIONS
00165     initPythia();
00166     //
00167     //........................................HepMC INITIALIZATIONS
00168     // Instantiate an IO strategy for reading from HEPEVT.
00169     HepMC::IO_HEPEVT hepevtio;
00170     // declare an instance of the event selection predicate
00171     IsGoodEventMyPythia is_good_event;
00172     //........................................EVENT LOOP
00173     int icount=0;
00174     int num_good_events=0;
00175     for ( int i = 1; i <= 100; i++ ) {
00176         icount++;
00177         if ( i%50==1 ) std::cout << "Processing Event Number " 
00178                                  << i << std::endl;
00179         call_pyevnt(); // generate one event with Pythia
00180         // pythia pyhepc routine convert common PYJETS in common HEPEVT
00181         call_pyhepc( 1 );
00182         HepMC::GenEvent* evt = hepevtio.read_next_event();
00183         // define the units (Pythia uses GeV and mm)
00184         evt->use_units(HepMC::Units::GEV, HepMC::Units::MM);
00185         // set number of multi parton interactions
00186         evt->set_mpi( pypars.msti[31-1] );
00187         // set cross section information
00188         evt->set_cross_section( getPythiaCrossSection() );
00189         // do event selection
00190         if ( is_good_event(evt) ) {
00191             std::cout << "Good Event Number " << i << std::endl;
00192             ++num_good_events;
00193         }
00194         // we also need to delete the created event from memory
00195         delete evt;
00196     }
00197     //........................................TERMINATION
00198     // write out some information from Pythia to the screen
00199     call_pystat( 1 );    
00200     //........................................PRINT RESULTS
00201     std::cout << num_good_events << " out of " << icount 
00202               << " processed events passed the cuts. Finished." << std::endl;
00203 }
00204 
00205 void pythia_in()
00206 {
00207     std::cout << std::endl;
00208     std::cout << "Begin pythia_in()" << std::endl;
00209     std::cout << "reading example_MyPythia.dat" << std::endl;
00210     //........................................define an input scope
00211     {
00212         // open input stream
00213         std::ifstream istr( "example_MyPythia.dat" );
00214         if( !istr ) {
00215           std::cerr << "example_ReadMyPythia: cannot open example_MyPythia.dat" << std::endl;
00216           exit(-1);
00217         }
00218         HepMC::IO_GenEvent ascii_in(istr);
00219         // open output stream (alternate method)
00220         HepMC::IO_GenEvent ascii_out("example_MyPythia2.dat",std::ios::out);
00221         // now read the file
00222         int icount=0;
00223         HepMC::GenEvent* evt = ascii_in.read_next_event();
00224         while ( evt ) {
00225             icount++;
00226             if ( icount%50==1 ) std::cout << "Processing Event Number " << icount
00227                                           << " its # " << evt->event_number() 
00228                                           << std::endl;
00229             // write the event out to the ascii file
00230             ascii_out << evt;
00231             delete evt;
00232             ascii_in >> evt;
00233         }
00234         //........................................PRINT RESULT
00235         std::cout << icount << " events found. Finished." << std::endl;
00236     } // ascii_out and istr destructors are called here
00237 }
00238 
00239 void pythia_in_out()
00240 {
00241     std::cout << std::endl;
00242     std::cout << "Begin pythia_in_out()" << std::endl;
00243     //........................................HEPEVT
00244     // Pythia 6.3 uses HEPEVT with 4000 entries and 8-byte floating point
00245     //  numbers. We need to explicitly pass this information to the 
00246     //  HEPEVT_Wrapper.
00247     //
00248     HepMC::HEPEVT_Wrapper::set_max_number_entries(4000);
00249     HepMC::HEPEVT_Wrapper::set_sizeof_real(8);
00250     //
00251     //........................................PYTHIA INITIALIZATIONS
00252     initPythia();
00253 
00254     //........................................HepMC INITIALIZATIONS
00255     //
00256     // Instantiate an IO strategy for reading from HEPEVT.
00257     HepMC::IO_HEPEVT hepevtio;
00258     //
00259     //........................................define the output scope
00260     {
00261         // Instantial an IO strategy to write the data to file - it uses the 
00262         //  same ParticleDataTable
00263         HepMC::IO_GenEvent ascii_io("example_MyPythiaRead.dat",std::ios::out);
00264         //
00265         //........................................EVENT LOOP
00266         for ( int i = 1; i <= 100; i++ ) {
00267             if ( i%50==1 ) std::cout << "Processing Event Number " 
00268                                      << i << std::endl;
00269             call_pyevnt();      // generate one event with Pythia
00270             // pythia pyhepc routine converts common PYJETS in common HEPEVT
00271             call_pyhepc( 1 );
00272             HepMC::GenEvent* evt = hepevtio.read_next_event();
00273             // define the units (Pythia uses GeV and mm)
00274             evt->use_units(HepMC::Units::GEV, HepMC::Units::MM);
00275             // set cross section information
00276             evt->set_cross_section( getPythiaCrossSection() );
00277             // add some information to the event
00278             evt->set_event_number(i);
00279             evt->set_signal_process_id(20);
00280             // write the event out to the ascii file
00281             ascii_io << evt;
00282             // we also need to delete the created event from memory
00283             delete evt;
00284         }
00285         //........................................TERMINATION
00286         // write out some information from Pythia to the screen
00287         call_pystat( 1 );    
00288     }  // ascii_io destructor is called here
00289     //
00290     //........................................define an input scope
00291     {
00292         // now read the file we wrote
00293         HepMC::IO_GenEvent ascii_in("example_MyPythiaRead.dat",std::ios::in);
00294         HepMC::IO_GenEvent ascii_io2("example_MyPythiaRead2.dat",std::ios::out);
00295         int icount=0;
00296         HepMC::GenEvent* evt = ascii_in.read_next_event();
00297         while ( evt ) {
00298             icount++;
00299             if ( icount%50==1 ) std::cout << "Processing Event Number " << icount
00300                                           << " its # " << evt->event_number() 
00301                                           << std::endl;
00302             // write the event out to the ascii file
00303             ascii_io2 << evt;
00304             delete evt;
00305             ascii_in >> evt;
00306         }
00307         //........................................PRINT RESULT
00308         std::cout << icount << " events found. Finished." << std::endl;
00309     } // ascii_io2 and ascii_in destructors are called here
00310 }
00311 
00312 void pythia_particle_out()
00313 {
00314     std::cout << std::endl;
00315     std::cout << "Begin pythia_particle_out()" << std::endl;
00316     //........................................HEPEVT
00317     // Pythia 6.1 uses HEPEVT with 4000 entries and 8-byte floating point
00318     //  numbers. We need to explicitly pass this information to the 
00319     //  HEPEVT_Wrapper.
00320     //
00321     HepMC::HEPEVT_Wrapper::set_max_number_entries(4000);
00322     HepMC::HEPEVT_Wrapper::set_sizeof_real(8);
00323     //
00324     //........................................PYTHIA INITIALIZATIONS
00325     initPythia();
00326 
00327     //........................................HepMC INITIALIZATIONS
00328     //
00329     // Instantiate an IO strategy for reading from HEPEVT.
00330     HepMC::IO_HEPEVT hepevtio;
00331     //
00332     { // begin scope of ascii_io
00333         // Instantiate an IO strategy to write the data to file 
00334         HepMC::IO_AsciiParticles ascii_io("example_PythiaParticle.dat",std::ios::out);
00335         //
00336         //........................................EVENT LOOP
00337         for ( int i = 1; i <= 100; i++ ) {
00338             if ( i%50==1 ) std::cout << "Processing Event Number " 
00339                                      << i << std::endl;
00340             call_pyevnt();      // generate one event with Pythia
00341             // pythia pyhepc routine converts common PYJETS in common HEPEVT
00342             call_pyhepc( 1 );
00343             HepMC::GenEvent* evt = hepevtio.read_next_event();
00344             // define the units (Pythia uses GeV and mm)
00345             evt->use_units(HepMC::Units::GEV, HepMC::Units::MM);
00346             // set cross section information
00347             evt->set_cross_section( getPythiaCrossSection() );
00348             // add some information to the event
00349             evt->set_event_number(i);
00350             evt->set_signal_process_id(20);
00351             // write the event out to the ascii file
00352             ascii_io << evt;
00353             // we also need to delete the created event from memory
00354             delete evt;
00355         }
00356         //........................................TERMINATION
00357         // write out some information from Pythia to the screen
00358         call_pystat( 1 );    
00359     } // end scope of ascii_io
00360 }
00361 

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