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/GenRanges.h"
00009 #include "HepMC/IO_GenEvent.h"
00010 #include "Rivet/Tools/RivetSTL.hh"
00011 #include "Rivet/Tools/Exceptions.hh"
00012 
00013 namespace Rivet {
00014 
00015 
00016   using HepMC::GenEvent;
00017   using HepMC::GenParticle;
00018   using HepMC::GenVertex;
00019 
00020   #if HEPMC_VERSION_CODE >= 2007000
00021   using HepMC::GenEventPtr;
00022   using HepMC::GenParticlePtr;
00023   using HepMC::GenVertexPtr;
00024   #else
00025   typedef GenEvent* GenEventPtr;
00026   typedef GenParticle* GenParticlePtr;
00027   typedef GenVertex* GenVertexPtr;
00028   #endif
00029 
00030 
00031   /// @todo Use mcutils?
00032 
00033 
00034   inline std::vector<GenParticle const *> particles(const GenEvent* ge) {
00035     assert(ge);
00036     std::vector<const GenParticle*> rtn;
00037     for (GenEvent::particle_const_iterator pi = ge->particles_begin(); pi != ge->particles_end(); ++pi)
00038       rtn.push_back(*pi);
00039     return rtn;
00040   }
00041 
00042   inline std::vector<GenParticlePtr> particles(GenEvent* ge) {
00043     assert(ge);
00044     std::vector<GenParticle*> rtn;
00045     for (GenEvent::particle_iterator pi = ge->particles_begin(); pi != ge->particles_end(); ++pi)
00046       rtn.push_back(*pi);
00047     return rtn;
00048   }
00049 
00050 
00051   inline std::vector<const GenVertex*> vertices(const GenEvent* ge) {
00052     std::vector<GenVertex const *> rtn;
00053     for (GenEvent::vertex_const_iterator vi = ge->vertices_begin(); vi != ge->vertices_end(); ++vi)
00054       rtn.push_back(*vi);
00055     return rtn;
00056   }
00057 
00058   inline std::vector<GenVertex*> vertices(GenEvent* ge) {
00059     std::vector<GenVertex*> rtn;
00060     for (GenEvent::vertex_iterator vi = ge->vertices_begin(); vi != ge->vertices_end(); ++vi)
00061       rtn.push_back(*vi);
00062     return rtn;
00063   }
00064 
00065 
00066   //////////////////////////
00067 
00068 
00069   inline std::vector<const GenParticle*> particles(const GenVertex* gv, HepMC::IteratorRange range=HepMC::relatives) {
00070     std::vector<GenParticle const *> rtn;
00071     /// @todo A particle_const_iterator on GenVertex would be nice...
00072     // Before HepMC 2.7.0 there were no GV::particles_const_iterators and constness consistency was all screwed up :-/
00073     #if HEPMC_VERSION_CODE >= 2007000
00074     // for (GenVertex::particle_iterator pi = gv->particles_const_begin(range); pi != gv->particles_const_end(range); ++pi)
00075     for (GenVertex::particle_iterator pi = gv->particles_begin(range); pi != gv->particles_end(range); ++pi)
00076       rtn.push_back(*pi);
00077     #else
00078     GenVertex* gv2 = const_cast<GenVertex*>(gv);
00079     for (GenVertex::particle_iterator pi = gv2->particles_begin(range); pi != gv2->particles_end(range); ++pi)
00080       rtn.push_back(const_cast<const GenParticle*>(*pi));
00081     #endif
00082     return rtn;
00083   }
00084 
00085   inline std::vector<GenParticle*> particles(GenVertex* gv, HepMC::IteratorRange range=HepMC::relatives) {
00086     std::vector<GenParticle*> rtn;
00087     for (GenVertex::particle_iterator pi = gv->particles_begin(range); pi != gv->particles_end(range); ++pi)
00088       rtn.push_back(*pi);
00089     return rtn;
00090   }
00091 
00092 
00093 
00094   // Get iterator ranges as wrapped begin/end pairs
00095   /// @note GenVertex _in and _out iterators are actually, secretly the same types *sigh*
00096   struct GenVertexIterRangeC {
00097     GenVertexIterRangeC(const GenVertex::particles_in_const_iterator& begin, const GenVertex::particles_in_const_iterator& end)
00098       : _begin(begin), _end(end) {  }
00099     const GenVertex::particles_in_const_iterator& begin() { return _begin; }
00100     const GenVertex::particles_in_const_iterator& end() { return _end; }
00101   private:
00102     const GenVertex::particles_in_const_iterator _begin, _end;
00103   };
00104 
00105   inline GenVertexIterRangeC particles_in(const GenVertex* gv) {
00106     return GenVertexIterRangeC(gv->particles_in_const_begin(), gv->particles_in_const_end());
00107   }
00108 
00109   inline GenVertexIterRangeC particles_out(const GenVertex* gv) {
00110     return GenVertexIterRangeC(gv->particles_out_const_begin(), gv->particles_out_const_end());
00111   }
00112 
00113 
00114 
00115   #if HEPMC_VERSION_CODE >= 2007000
00116 
00117   // Get iterator ranges as wrapped begin/end pairs
00118   /// @note GenVertex _in and _out iterators are actually, secretly the same types *sigh*
00119   struct GenVertexIterRange {
00120     GenVertexIterRange(const GenVertex::particles_in_iterator& begin, const GenVertex::particles_in_iterator& end)
00121       : _begin(begin), _end(end) {  }
00122     const GenVertex::particles_in_iterator& begin() { return _begin; }
00123     const GenVertex::particles_in_iterator& end() { return _end; }
00124   private:
00125     const GenVertex::particles_in_iterator _begin, _end;
00126   };
00127 
00128   inline GenVertexIterRange particles_in(const GenVertex* gv) {
00129     return GenVertexIterRange(gv->particles_in_begin(), gv->particles_in_end());
00130   }
00131 
00132   inline GenVertexIterRange particles_out(const GenVertex* gv) {
00133     return GenVertexIterRange(gv->particles_out_begin(), gv->particles_out_end());
00134   }
00135 
00136   #endif
00137 
00138 
00139   //////////////////////////
00140 
00141 
00142   /// Get the direct parents or all-ancestors of GenParticle @a gp
00143   inline std::vector<const GenParticle*> particles_in(const GenParticle* gp, HepMC::IteratorRange range=HepMC::ancestors) {
00144     if (range != HepMC::parents && range != HepMC::ancestors)
00145       throw UserError("Requested particles_in(GenParticle*) with a non-'in' iterator range");
00146     if (!gp->production_vertex()) return std::vector<const GenParticle*>();
00147     #if HEPMC_VERSION_CODE >= 2007000
00148     return particles(gp->production_vertex(), range);
00149     #else
00150     // Before HepMC 2.7.0 the constness consistency of methods and their return types was all screwed up :-/
00151     std::vector<const GenParticle*> rtn;
00152     foreach (GenParticle* gp2, particles(gp->production_vertex(), range))
00153       rtn.push_back( const_cast<const GenParticle*>(gp2) );
00154     return rtn;
00155     #endif
00156   }
00157 
00158   /// Get the direct parents or all-ancestors of GenParticle @a gp
00159   inline std::vector<GenParticle*> particles_in(GenParticle* gp, HepMC::IteratorRange range=HepMC::ancestors) {
00160     if (range != HepMC::parents && range != HepMC::ancestors)
00161       throw UserError("Requested particles_in(GenParticle*) with a non-'in' iterator range");
00162     return (gp->production_vertex()) ? particles(gp->production_vertex(), range) : std::vector<GenParticle*>();
00163   }
00164 
00165 
00166   /// Get the direct children or all-descendents of GenParticle @a gp
00167   inline std::vector<const GenParticle*> particles_out(const GenParticle* gp, HepMC::IteratorRange range=HepMC::descendants) {
00168     if (range != HepMC::children && range != HepMC::descendants)
00169       throw UserError("Requested particles_out(GenParticle*) with a non-'out' iterator range");
00170     if (!gp->end_vertex()) return std::vector<const GenParticle*>();
00171     #if HEPMC_VERSION_CODE >= 2007000
00172     return particles(gp->end_vertex(), range);
00173     #else
00174     // Before HepMC 2.7.0 the constness consistency of methods and their return types was all screwed up :-/
00175     std::vector<const GenParticle*> rtn;
00176     foreach (GenParticle* gp2, particles(gp->end_vertex(), range))
00177       rtn.push_back( const_cast<const GenParticle*>(gp2) );
00178     return rtn;
00179     #endif
00180   }
00181 
00182   /// Get the direct children or all-descendents of GenParticle @a gp
00183   inline std::vector<GenParticle*> particles_out(GenParticle* gp, HepMC::IteratorRange range=HepMC::descendants) {
00184     if (range != HepMC::children && range != HepMC::descendants)
00185       throw UserError("Requested particles_out(GenParticle*) with a non-'out' iterator range");
00186     return (gp->end_vertex()) ? particles(gp->end_vertex(), range) : std::vector<GenParticle*>();
00187   }
00188 
00189 
00190   /// Get any relatives of GenParticle @a gp
00191   inline std::vector<const GenParticle*> particles(const GenParticle* gp, HepMC::IteratorRange range=HepMC::ancestors) {
00192     if (range == HepMC::parents || range == HepMC::ancestors)
00193       return particles_in(gp, range);
00194     if (range == HepMC::children || range == HepMC::descendants)
00195       return particles_in(gp, range);
00196     throw UserError("Requested particles(GenParticle*) with an unsupported iterator range");
00197   }
00198 
00199 
00200 }
00201 
00202 #endif