StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
FastJet3.h
1 // FastJet3.h is a part of the PYTHIA event generator.
2 // Copyright (C) 2020 Torbjorn Sjostrand.
3 // PYTHIA is licenced under the GNU GPL v2 or later, see COPYING for details.
4 // Please respect the MCnet Guidelines, see GUIDELINES for details.
5 
6 // This header file written by Gavin Salam.
7 
8 #ifndef Pythia8_FastJet3_H
9 #define Pythia8_FastJet3_H
10 
11 //----------------------------------------------------------------------
76 // ----------------------------------------------------------------------
77 // Copyright 2011 by Matteo Cacciari, Gavin Salam and Gregory
78 // Soyez. Permission is granted to redistribute this file and modify
79 // it, as long as this notice is retained and any changes are clearly
80 // marked. No warranties are provided!
81 // ----------------------------------------------------------------------
82 
83 #include "fastjet/config.h" // will allow a test for FJ3
84 #include "fastjet/ClusterSequence.hh" // also gives PseudoJet & JetDefinition
85 #include "fastjet/Selector.hh"
86 #include "Pythia8/Event.h" // this is what we need from Pythia8
87 
88 // FASTJET_VERSION is only defined from version 3 onwards so we can
89 // use it to test that we have a sufficiently recent version
90 #ifndef FASTJET_VERSION
91 #error "FastJet3 is required in order to obtain the features of this interface"
92 #endif
93 
94 FASTJET_BEGIN_NAMESPACE // place the code here inside the FJ namespace
95 
101 class Py8Particle: public Pythia8::Particle,
102  public PseudoJet::UserInfoBase {
103 public:
104  Py8Particle(const Pythia8::Particle & particle) : Particle(particle),
105  mIndex(particle.index()) {}
106  virtual int index() const override {return mIndex;}
107 private:
108  int mIndex;
109 };
110 
113 template<>
114 inline PseudoJet::PseudoJet(const Pythia8::Particle & particle) {
115  reset(particle.px(),particle.py(),particle.pz(), particle.e());
116  set_user_index( particle.index() );
117  set_user_info(new Py8Particle(particle));
118 }
119 
122 template<>
123 inline PseudoJet::PseudoJet(const Pythia8::Vec4 & particle) {
124  reset(particle.px(),particle.py(),particle.pz(), particle.e());
125 }
126 
127 
139 template<class T> class SelectorWorkerPy8 : public SelectorWorker {
140 public:
142  typedef T (Pythia8::Particle::*Py8ParticleFnPtr)() const;
143 
146  SelectorWorkerPy8(Py8ParticleFnPtr member_fn_ptr, T value) :
147  _member_fn_ptr(member_fn_ptr), _value(value) {};
148 
153  bool pass(const PseudoJet & p) const {
154  const Pythia8::Particle * py8_particle
155  = dynamic_cast<const Pythia8::Particle *>(p.user_info_ptr());
156  if (py8_particle == 0) {
157  return false; // no info, so false
158  } else {
159  return (py8_particle->*_member_fn_ptr)() == _value;
160  }
161  }
162 private:
163  Py8ParticleFnPtr _member_fn_ptr;
164  T _value;
165 };
166 
174 inline Selector SelectorIsFinal () {return
175  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isFinal , true));}
176 inline Selector SelectorIsCharged () {return
177  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isCharged , true));}
178 inline Selector SelectorIsNeutral () {return
179  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isNeutral , true));}
180 inline Selector SelectorIsResonance() {return
181  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isResonance,true));}
182 inline Selector SelectorIsVisible () {return
183  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isVisible , true));}
184 inline Selector SelectorIsLepton () {return
185  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isLepton , true));}
186 inline Selector SelectorIsQuark () {return
187  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isQuark , true));}
188 inline Selector SelectorIsGluon () {return
189  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isGluon , true));}
190 inline Selector SelectorIsDiquark () {return
191  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isDiquark , true));}
192 inline Selector SelectorIsParton () {return
193  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isParton , true));}
194 inline Selector SelectorIsHadron () {return
195  Selector(new SelectorWorkerPy8<bool>(&Pythia8::Particle::isHadron , true));}
197 
205 inline Selector SelectorId (int i) {return
206  Selector(new SelectorWorkerPy8<int>(&Pythia8::Particle::id , i));}
207 inline Selector SelectorIdAbs (int i) {return
208  Selector(new SelectorWorkerPy8<int>(&Pythia8::Particle::idAbs , i));}
209 inline Selector SelectorStatus (int i) {return
210  Selector(new SelectorWorkerPy8<int>(&Pythia8::Particle::status , i));}
211 inline Selector SelectorStatusAbs(int i) {return
212  Selector(new SelectorWorkerPy8<int>(&Pythia8::Particle::statusAbs, i));}
214 
215 
216 FASTJET_END_NAMESPACE
217 
218 #endif // Pythia8_FastJet3_H
SelectorWorkerPy8(Py8ParticleFnPtr member_fn_ptr, T value)
Definition: FastJet3.h:146
bool pass(const PseudoJet &p) const
Definition: FastJet3.h:153
T(Pythia8::Particle::* Py8ParticleFnPtr)() const
the typedef helps with the notation for member function pointers
Definition: FastJet3.h:137