#include <HepMCConverter.h>
Classes | |
struct | Vertex |
Public Types | |
typedef map< tcColinePtr, long > | FlowMap |
typedef Traits::EventT | GenEvent |
typedef Traits::ParticleT | GenParticle |
typedef Traits::VertexT | GenVertex |
typedef map< const Vertex *, GenVertex * > | GenVertexMap |
typedef map< tcPPtr, GenParticle * > | ParticleMap |
typedef Traits::PdfInfoT | PdfInfo |
typedef map< tcPPtr, Vertex * > | VertexMap |
Static Public Member Functions | |
static GenEvent * | convert (const Event &ev, bool nocopies=false, Energy eunit=Traits::defaultEnergyUnit(), Length lunit=Traits::defaultLengthUnit()) |
static void | convert (const Event &ev, GenEvent &gev, bool nocopies=false) |
static void | convert (const Event &ev, GenEvent &gev, bool nocopies, Energy eunit, Length lunit) |
Private Member Functions | |
GenParticle * | createParticle (tcPPtr p) const |
GenVertex * | createVertex (Vertex *v) |
HepMCConverter (const Event &ev, GenEvent &gev, bool nocopies, Energy eunit, Length lunit) | |
HepMCConverter () | |
HepMCConverter (const HepMCConverter &) | |
HepMCConverter (const Event &ev, bool nocopies, Energy eunit, Length lunit) | |
void | init (const Event &ev, bool nocopies) |
void | join (tcPPtr parent, tcPPtr child) |
HepMCConverter & | operator= (const HepMCConverter &) |
void | setPdfInfo (const Event &e) |
void | sortTopologically (tcPVector &allv) |
Private Attributes | |
VertexMap | decv |
tcEHPtr | eh |
Energy | energyUnit |
FlowMap | flowmap |
GenEvent * | geneve |
Length | lengthUnit |
ParticleMap | pmap |
VertexMap | prov |
vector< Vertex > | vertices |
GenVertexMap | vmap |
The HepMCConverter defines only one public static function which converts a ThePEG::Event object to a HepMC::GenEvent
. All mother-daughter relationships and colour information is preserved.
Definition at line 32 of file HepMCConverter.h.
typedef map<tcColinePtr,long> ThePEG::HepMCConverter< HepMCEventT, Traits >::FlowMap |
Map ThePEG colour lines to HepMC colour indices.
Definition at line 63 of file HepMCConverter.h.
typedef Traits::EventT ThePEG::HepMCConverter< HepMCEventT, Traits >::GenEvent |
Forward typedefs from Traits class.
Definition at line 55 of file HepMCConverter.h.
typedef Traits::ParticleT ThePEG::HepMCConverter< HepMCEventT, Traits >::GenParticle |
Forward typedefs from Traits class.
Definition at line 53 of file HepMCConverter.h.
typedef Traits::VertexT ThePEG::HepMCConverter< HepMCEventT, Traits >::GenVertex |
Forward typedefs from Traits class.
Definition at line 57 of file HepMCConverter.h.
typedef map<const Vertex *, GenVertex*> ThePEG::HepMCConverter< HepMCEventT, Traits >::GenVertexMap |
Map vertices to GenVertex
Definition at line 67 of file HepMCConverter.h.
typedef map<tcPPtr,GenParticle*> ThePEG::HepMCConverter< HepMCEventT, Traits >::ParticleMap |
Map ThePEG particles to HepMC particles.
Definition at line 61 of file HepMCConverter.h.
typedef Traits::PdfInfoT ThePEG::HepMCConverter< HepMCEventT, Traits >::PdfInfo |
Forward typedefs from Traits class.
Definition at line 59 of file HepMCConverter.h.
typedef map<tcPPtr,Vertex*> ThePEG::HepMCConverter< HepMCEventT, Traits >::VertexMap |
Map ThePEG particles to vertices.
Definition at line 65 of file HepMCConverter.h.
ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter | ( | const Event & | ev, |
bool | nocopies, | ||
Energy | eunit, | ||
Length | lunit | ||
) | [private] |
The proper constructors are private. The class is only instantiated within the convert method.
Definition at line 65 of file HepMCConverter.cc.
References ThePEG::HepMCConverter< HepMCEventT, Traits >::geneve, and ThePEG::HepMCConverter< HepMCEventT, Traits >::init().
: energyUnit(eunit), lengthUnit(lunit) { geneve = Traits::newEvent(ev.number(), ev.weight()); init(ev, nocopies); }
ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter | ( | const Event & | ev, |
GenEvent & | gev, | ||
bool | nocopies, | ||
Energy | eunit, | ||
Length | lunit | ||
) | [private] |
The proper constructors are private. The class is only instantiated within the convert method.
Definition at line 76 of file HepMCConverter.cc.
References ThePEG::HepMCConverter< HepMCEventT, Traits >::geneve, and ThePEG::HepMCConverter< HepMCEventT, Traits >::init().
: energyUnit(eunit), lengthUnit(lunit) { geneve = &gev; init(ev, nocopies); }
ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter | ( | ) | [private] |
Default constructor is unimplemented and private and should never be used.
ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter | ( | const HepMCConverter< HepMCEventT, Traits > & | ) | [private] |
Copy constructor is unimplemented and private and should never be used.
HepMCConverter< HepMCEventT, Traits >::GenEvent * ThePEG::HepMCConverter< HepMCEventT, Traits >::convert | ( | const Event & | ev, |
bool | nocopies = false , |
||
Energy | eunit = Traits::defaultEnergyUnit() , |
||
Length | lunit = Traits::defaultLengthUnit() |
||
) | [static] |
Convert a ThePEG::Event to a HepMC::GenEvent. The caller is responsible for deleting the constructed GenEvent object. If nocopies is true, only final copies of particles connected with Particle::previous() and Particle::next() will be entered in the HepMC::GenEvent. In the GenEvent object, the energy/momentum variables will be in units of eunit and lengths variables in units of lunit.
Definition at line 43 of file HepMCConverter.cc.
References ThePEG::HepMCConverter< HepMCEventT, Traits >::geneve.
{
HepMCConverter<HepMCEventT,Traits> converter(ev, nocopies, eunit, lunit);
return converter.geneve;
}
void ThePEG::HepMCConverter< HepMCEventT, Traits >::convert | ( | const Event & | ev, |
GenEvent & | gev, | ||
bool | nocopies, | ||
Energy | eunit, | ||
Length | lunit | ||
) | [static] |
Convert a ThePEG::Event to a HepMC::GenEvent. The caller supplies a GenEvent object, gev, which will be filled. If nocopies is true, only final copies of particles connected with Particle::previous() and Particle::next() will be entered in the HepMC::GenEvent. In the GenEvent object, the energy/momentum variables will be in units of eunit and lengths variables in units of lunit.
Definition at line 58 of file HepMCConverter.cc.
{ HepMCConverter<HepMCEventT,Traits> converter(ev, gev, nocopies, eunit, lunit); }
void ThePEG::HepMCConverter< HepMCEventT, Traits >::convert | ( | const Event & | ev, |
GenEvent & | gev, | ||
bool | nocopies = false |
||
) | [static] |
Convert a ThePEG::Event to a HepMC::GenEvent. The caller supplies a GenEvent object, gev, which will be filled. If nocopies is true, only final copies of particles connected with Particle::previous() and Particle::next() will be entered in the HepMC::GenEvent. In the GenEvent object, the energy/momentum variables will be in units of eunit and lengths variables in units of lunit.
Definition at line 50 of file HepMCConverter.cc.
{ HepMCConverter<HepMCEventT,Traits> converter(ev, gev, nocopies, Traits::momentumUnit(gev), Traits::lengthUnit(gev)); }
HepMCConverter< HepMCEventT, Traits >::GenParticle * ThePEG::HepMCConverter< HepMCEventT, Traits >::createParticle | ( | tcPPtr | p | ) | const [private] |
Create a GenParticle from a ThePEG Particle.
Definition at line 225 of file HepMCConverter.cc.
References ntuplemaker::status, and launcher::step.
{ int status = 1; if ( !p->children().empty() || p->next() ) { tStepPtr step = p->birthStep(); if ((!step || (step && (!step->handler() || step->handler() == eh))) && p->id() != 82) status = 3; else status = 2; } GenParticle * gp = Traits::newParticle(p->momentum(), p->id(), status, energyUnit); if ( p->spinInfo() && p->spinInfo()->hasPolarization() ) { DPair pol = p->spinInfo()->polarization(); Traits::setPolarization(*gp, pol.first, pol.second); } return gp; }
HepMCConverter< HepMCEventT, Traits >::GenVertex * ThePEG::HepMCConverter< HepMCEventT, Traits >::createVertex | ( | Vertex * | v | ) | [private] |
Create a GenVertex from a temporary Vertex.
Definition at line 268 of file HepMCConverter.cc.
References ThePEG::HepMCConverter< HepMCEventT, Traits >::Vertex::in, ThePEG::HepMCConverter< HepMCEventT, Traits >::Vertex::out, and AlCaHLTBitMon_ParallelJobs::p.
{ if ( !v ) throw HepMCConverterException() << "Found internal null Vertex." << Exception::abortnow; GenVertex * gv = new GenVertex(); // We assume that the vertex position is the average of the decay // vertices of all incoming and the creation vertices of all // outgoing particles in the lab. Note that this will probably not // be useful information for very small distances. LorentzPoint p; for ( tcParticleSet::iterator it = v->in.begin(); it != v->in.end(); ++it ) { p += (**it).labDecayVertex(); Traits::addIncoming(*gv, pmap[*it]); } for ( tcParticleSet::iterator it = v->out.begin(); it != v->out.end(); ++it ) { p += (**it).labVertex(); Traits::addOutgoing(*gv, pmap[*it]); } p /= double(v->in.size() + v->out.size()); Traits::setPosition(*gv, p, lengthUnit); return gv; }
void ThePEG::HepMCConverter< HepMCEventT, Traits >::init | ( | const Event & | ev, |
bool | nocopies | ||
) | [private] |
Common init function used by the constructors.
Definition at line 93 of file HepMCConverter.cc.
References cond::ecalcond::all, i, join(), prof2calltree::l, N, AlCaHLTBitMon_ParallelJobs::p, and createTree::pp.
Referenced by ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter().
{ if ( Traits::hasUnits() ) { if ( lengthUnit != millimeter && lengthUnit != centimeter ) throw HepMCConverterException() << "Length unit used for HepMC::GenEvent was not MM nor CM." << Exception::runerror; if ( energyUnit != GeV && energyUnit != MeV ) throw HepMCConverterException() << "Momentum unit used for HepMC::GenEvent was not GEV nor MEV." << Exception::runerror; } else { if ( lengthUnit != millimeter ) throw HepMCConverterException() << "Length unit used for HepMC::GenEvent was not MM and the " << "HepMCTraits class claims that the used version of HepMC " << "cannot handle user-defined units." << Exception::runerror; if ( energyUnit != GeV ) throw HepMCConverterException() << "Momentum unit used for HepMC::GenEvent was not GEV and the " << "HepMCTraits class claims that the used version of HepMC " << "cannot handle user-defined units." << Exception::runerror; } Traits::setUnits(*geneve, energyUnit, lengthUnit); if ( ev.primaryCollision() && ( eh = dynamic_ptr_cast<tcEHPtr>(ev.primaryCollision()->handler()) ) ) { // Get general event info if present. } // Extract all particles and order them. tcPVector all; ev.select(back_inserter(all), SelectAll()); vertices.reserve(all.size()*2); sortTopologically(all); // Create GenParticle's and map them to the ThePEG particles. for ( int i = 0, N = all.size(); i < N; ++i ) { tcPPtr p = all[i]; if ( nocopies && p->next() ) continue; if ( pmap.find(p) != pmap.end() ) continue; GenParticle * gp = pmap[p] = createParticle(p); if ( p->hasColourInfo() ) { // Check if the particle is connected to colour lines, in which // case the lines are mapped to an integer and set in the // GenParticle's Flow info. tcColinePtr l; if ( l = p->colourLine() ) { if ( !member(flowmap, l) ) flowmap[l] = flowmap.size() + 500; Traits::setColourLine(*gp, 1, flowmap[l]); } if ( l = p->antiColourLine() ) { if ( !member(flowmap, l) ) flowmap[l] = flowmap.size() + 500; Traits::setColourLine(*gp, 2, flowmap[l]); } } if ( !p->children().empty() || p->next() ) { // If the particle has children it should have a decay vertex: vertices.push_back(Vertex()); decv[p] = &vertices.back(); vertices.back().in.insert(p); } if ( !p->parents().empty() || p->previous() || (p->children().empty() && !p->next()) ) { // If the particle has parents it should have a production // vertex. If neither parents or children it should still have a // dummy production vertex. vertices.push_back(Vertex()); prov[p] = &vertices.back(); vertices.back().out.insert(p); } } // Now go through the the particles again, and join the vertices. for ( int i = 0, N = all.size(); i < N; ++i ) { tcPPtr p = all[i]; if ( nocopies ) { if ( p->next() ) continue; for ( int i = 0, N = p->children().size(); i < N; ++i ) join(p, p->children()[i]->final()); tcPPtr pp = p; while ( pp->parents().empty() && pp->previous() ) pp = pp->previous(); for ( int i = 0, N = pp->parents().size(); i < N; ++i ) join(pp->parents()[i]->final(), p); } else { for ( int i = 0, N = p->children().size(); i < N; ++i ) join(p, p->children()[i]); if ( p->next() ) join(p, p->next()); for ( int i = 0, N = p->parents().size(); i < N; ++i ) join(p->parents()[i], p); if ( p->previous() ) join(p->previous(), p); } } // Time to create the GenVertex's for ( typename VertexMap::iterator it = prov.begin(); it != prov.end(); ++it ) if ( !member(vmap, it->second) ) vmap[it->second] = createVertex(it->second); for ( typename VertexMap::iterator it = decv.begin(); it != decv.end(); ++it ) if ( !member(vmap, it->second) ) vmap[it->second] = createVertex(it->second); // Add the vertices. for ( typename GenVertexMap::iterator it = vmap.begin(); it != vmap.end(); ++it ) Traits::addVertex(*geneve, it->second); // Now find the primary signal process vertex defined to be the // decay vertex of the first parton coming into the primary hard // sub-collision. tSubProPtr sub = ev.primarySubProcess(); if ( sub && sub->incoming().first ) { const Vertex * prim = decv[sub->incoming().first]; Traits::setSignalProcessVertex(*geneve, vmap[prim]); vmap.erase(prim); } // and the incoming beam particles geneve->set_beam_particles(pmap[ev.incoming().first], pmap[ev.incoming().second]); // and the PDF info }
void ThePEG::HepMCConverter< HepMCEventT, Traits >::join | ( | tcPPtr | parent, |
tcPPtr | child | ||
) | [private] |
Join the decay vertex of the parent with the decay vertex of the child.
Definition at line 247 of file HepMCConverter.cc.
References ThePEG::HepMCConverter< HepMCEventT, Traits >::Vertex::in, ThePEG::HepMCConverter< HepMCEventT, Traits >::Vertex::out, and dbtoconf::parent.
{ Vertex * dec = decv[parent]; Vertex * pro = prov[child]; if ( !pro || !dec ) throw HepMCConverterException() << "Found a reference to a ThePEG::Particle which was not in the Event." << Exception::eventerror; if ( pro == dec ) return; while ( !pro->in.empty() ) { dec->in.insert(*(pro->in.begin())); decv[*(pro->in.begin())] = dec; pro->in.erase(pro->in.begin()); } while ( !pro->out.empty() ) { dec->out.insert(*(pro->out.begin())); prov[*(pro->out.begin())] = dec; pro->out.erase(pro->out.begin()); } }
HepMCConverter& ThePEG::HepMCConverter< HepMCEventT, Traits >::operator= | ( | const HepMCConverter< HepMCEventT, Traits > & | ) | [private] |
Assignment is unimplemented and private and should never be used.
void ThePEG::HepMCConverter< HepMCEventT, Traits >::setPdfInfo | ( | const Event & | e | ) | [private] |
Create and set a PdfInfo object for the event
Definition at line 297 of file HepMCConverter.cc.
References pileupReCalc_HLTpaths::scale, and mathSSE::sqrt().
{ // ids of the partons going into the primary sub process tSubProPtr sub = e.primarySubProcess(); int id1 = sub->incoming().first ->id(); int id2 = sub->incoming().second->id(); // get the event handler tcEHPtr eh = dynamic_ptr_cast<tcEHPtr>(e.handler()); // get the values of x double x1 = eh->lastX1(); double x2 = eh->lastX2(); // get the pdfs pair<PDF,PDF> pdfs; pdfs.first = eh->pdf<PDF>(sub->incoming().first ); pdfs.second = eh->pdf<PDF>(sub->incoming().second); // get the scale Energy2 scale = eh->lastScale(); // get the values of the pdfs double xf1 = pdfs.first.xfx(sub->incoming().first->dataPtr(), scale, x1); double xf2 = pdfs.second.xfx(sub->incoming().second->dataPtr(), scale, x2); Traits::setPdfInfo(*geneve, id1, id2, x1, x2, sqrt(scale/GeV2), xf1, xf2); }
void ThePEG::HepMCConverter< HepMCEventT, Traits >::sortTopologically | ( | tcPVector & | allv | ) | [private] |
Sort vertices topologically
Definition at line 322 of file HepMCConverter.cc.
References parents, and swap().
{ // order the particles topologically std::set<tcPPtr> visited; for(unsigned int head = 0; head < all.size(); head++) { bool vetoed = true; for(unsigned int cur = head; cur < all.size(); cur++) { vetoed = false; for(tParticleVector::const_iterator iter = all[cur]->parents().begin(); iter != all[cur]->parents().end(); ++iter) { if (visited.find(*iter) == visited.end()) { vetoed = true; break; } } if (!vetoed) { if (cur != head) std::swap(all[head], all[cur]); break; } } visited.insert(all[head]); } visited.clear(); }
VertexMap ThePEG::HepMCConverter< HepMCEventT, Traits >::decv [private] |
The mapping of particles to their decy vertices.
Definition at line 203 of file HepMCConverter.h.
tcEHPtr ThePEG::HepMCConverter< HepMCEventT, Traits >::eh [private] |
The primary event handler
Definition at line 223 of file HepMCConverter.h.
Energy ThePEG::HepMCConverter< HepMCEventT, Traits >::energyUnit [private] |
The energy unit to be used in the GenEvent.
Definition at line 213 of file HepMCConverter.h.
FlowMap ThePEG::HepMCConverter< HepMCEventT, Traits >::flowmap [private] |
The translation table between ThePEG ColourLine objects and HepMC Flow indices.
Definition at line 188 of file HepMCConverter.h.
GenEvent* ThePEG::HepMCConverter< HepMCEventT, Traits >::geneve [private] |
The constructed GenEvent.
Definition at line 176 of file HepMCConverter.h.
Referenced by ThePEG::HepMCConverter< HepMCEventT, Traits >::convert(), and ThePEG::HepMCConverter< HepMCEventT, Traits >::HepMCConverter().
Length ThePEG::HepMCConverter< HepMCEventT, Traits >::lengthUnit [private] |
The length unit to be used in the GenEvent.
Definition at line 218 of file HepMCConverter.h.
ParticleMap ThePEG::HepMCConverter< HepMCEventT, Traits >::pmap [private] |
The translation table between the ThePEG particles and the GenParticles.
Definition at line 182 of file HepMCConverter.h.
VertexMap ThePEG::HepMCConverter< HepMCEventT, Traits >::prov [private] |
The mapping of particles to their production vertices.
Definition at line 198 of file HepMCConverter.h.
vector<Vertex> ThePEG::HepMCConverter< HepMCEventT, Traits >::vertices [private] |
All temporary vertices created.
Definition at line 193 of file HepMCConverter.h.
GenVertexMap ThePEG::HepMCConverter< HepMCEventT, Traits >::vmap [private] |
The mapping between temporary vertices and the created GenVertex Objects.
Definition at line 208 of file HepMCConverter.h.