rivet is hosted by Hepforge, IPPP Durham
RivetHepMC.hh
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #ifndef RIVET_RivetHepMC_HH
00003 #define RIVET_RivetHepMC_HH
00004 
00005 #include "HepMC/GenEvent.h"
00006 #include "HepMC/GenParticle.h"
00007 #include "HepMC/GenVertex.h"
00008 #include "HepMC/IO_GenEvent.h"
00009 #include "Rivet/Tools/RivetSTL.hh"
00010 #include "Rivet/Tools/RivetBoost.hh"
00011 #include "Rivet/Exceptions.hh"
00012 
00013 namespace Rivet {
00014 
00015 
00016   using HepMC::GenEvent;
00017   using HepMC::GenParticle;
00018   using HepMC::GenVertex;
00019 
00020 
00021   /// @todo Add functions from mcutils
00022 
00023 
00024   inline std::vector<GenParticle*> particles(const GenEvent& ge) {
00025     std::vector<GenParticle*> rtn;
00026     for (GenEvent::particle_const_iterator pi = ge.particles_begin(); pi != ge.particles_end(); ++pi) {
00027       rtn.push_back(*pi);
00028     }
00029     return rtn;
00030   }
00031   inline std::vector<GenParticle*> particles(const GenEvent* ge) {
00032     assert(ge);
00033     return particles(*ge);
00034   }
00035   //  inline std::pair<GenEvent::particle_const_iterator, GenEvent::particle_const_iterator> particles(const GenEvent& ge) {
00036   //    return make_pair(ge.particles_begin(), ge.particles_end());
00037   //  }
00038   //  inline std::pair<GenEvent::particle_const_iterator, GenEvent::particle_const_iterator> particles(const GenEvent* ge) {
00039   //    assert(ge);
00040   //    return particles(*ge);
00041   //  }
00042 
00043 
00044 
00045   inline std::vector<GenVertex*> vertices(const GenEvent& ge) {
00046     std::vector<GenVertex*> rtn;
00047     for (GenEvent::vertex_const_iterator vi = ge.vertices_begin(); vi != ge.vertices_end(); ++vi) {
00048       rtn.push_back(*vi);
00049     }
00050     return rtn;
00051   }
00052   inline std::vector<GenVertex*> vertices(const GenEvent* ge) {
00053     assert(ge);
00054     return vertices(*ge);
00055   }
00056   //  inline std::pair<GenEvent::vertex_const_iterator, GenEvent::vertex_const_iterator> vertices(const GenEvent& ge) {
00057   //    return make_pair(ge.vertices_begin(), ge.vertices_end());
00058   //  }
00059   //  inline std::pair<GenEvent::vertex_const_iterator, GenEvent::vertex_const_iterator> vertices(const GenEvent* ge) {
00060   //    assert(ge);
00061   //    return vertices(*ge);
00062   //  }
00063 
00064 
00065   ////////////////////////////////////
00066 
00067 
00068   inline std::pair<GenVertex::particles_in_const_iterator, GenVertex::particles_in_const_iterator>
00069   particles_in(const GenVertex* gv) {
00070     return make_pair(gv->particles_in_const_begin(), gv->particles_in_const_end());
00071   }
00072 
00073 
00074   inline std::pair<GenVertex::particles_out_const_iterator, GenVertex::particles_out_const_iterator>
00075   particles_out(const GenVertex* gv) {
00076     return make_pair(gv->particles_out_const_begin(), gv->particles_out_const_end());
00077   }
00078 
00079 
00080   inline std::vector<GenParticle*> particles(GenVertex* gv, HepMC::IteratorRange range=HepMC::relatives) {
00081     std::vector<GenParticle*> rtn;
00082     for (GenVertex::particle_iterator pi = gv->particles_begin(range); pi != gv->particles_end(range); ++pi) {
00083       rtn.push_back(*pi);
00084     }
00085     return rtn;
00086   }
00087   //  inline std::pair<GenVertex::particle_iterator, GenVertex::particle_iterator> particles(GenVertex* gv, HepMC::IteratorRange range=HepMC::relatives) {
00088   //    return make_pair(gv->particles_begin(range), gv->particles_end(range));
00089   //  }
00090 
00091 
00092   //////////////////////////
00093 
00094 
00095   /// Get the direct parents or all-ancestors of GenParticle @a gp
00096   inline std::vector<GenParticle*> particles_in(const GenParticle* gp, HepMC::IteratorRange range=HepMC::ancestors) {
00097     if (range != HepMC::parents && range != HepMC::ancestors)
00098       throw UserError("Requested particles_in(GenParticle*) with a non-'in' iterator range");
00099     return (gp->production_vertex()) ? particles(gp->production_vertex(), range) : std::vector<GenParticle*>();
00100   }
00101 
00102 
00103   /// Get the direct children or all-descendents of GenParticle @a gp
00104   inline std::vector<GenParticle*> particles_out(const GenParticle* gp, HepMC::IteratorRange range=HepMC::descendants) {
00105     if (range != HepMC::children && range != HepMC::descendants)
00106       throw UserError("Requested particles_out(GenParticle*) with a non-'out' iterator range");
00107     return (gp->end_vertex()) ? particles(gp->end_vertex(), range) : std::vector<GenParticle*>();
00108   }
00109 
00110 
00111 }
00112 
00113 #endif