JetAlg.hh
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #ifndef RIVET_JetAlg_HH 00003 #define RIVET_JetAlg_HH 00004 00005 #include "Rivet/Projection.hh" 00006 #include "Rivet/Projections/FinalState.hh" 00007 #include "Rivet/Projections/VisibleFinalState.hh" 00008 #include "Rivet/Particle.hh" 00009 #include "Rivet/Jet.hh" 00010 00011 namespace Rivet { 00012 00013 00014 /// Abstract base class for projections which can return a set of {@link Jet}s. 00015 class JetAlg : public Projection { 00016 public: 00017 00018 /// Constructor 00019 JetAlg(const FinalState& fs); 00020 00021 JetAlg() {}; 00022 00023 /// Clone on the heap. 00024 virtual const Projection* clone() const = 0; 00025 00026 /// Destructor 00027 virtual ~JetAlg() { } 00028 00029 /// @brief Include invisible particles in jet construction. 00030 /// The default behaviour is that jets are only constructed from visible 00031 /// (i.e. charged under an SM gauge group) particles. Some jet studies, 00032 /// including those from ATLAS, use a definition in which neutrinos from hadron 00033 /// decays are included (via MC correction) in the experimental jet definition. 00034 /// Setting this flag to true avoids the automatic restriction to a VisibleFinalState. 00035 void useInvisibles(bool useinvis=true) { 00036 _useInvisibles = useinvis; 00037 } 00038 00039 /// Get jets in no guaranteed order, with optional cuts on \f$ p_\perp \f$ and rapidity. 00040 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00041 virtual Jets jets(double ptmin=0.0, double ptmax=MAXDOUBLE, 00042 double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, 00043 RapScheme rapscheme=PSEUDORAPIDITY) const { 00044 const Jets rawjets = _jets(ptmin); 00045 Jets rtn; 00046 MSG_DEBUG("Raw jet size (with pTmin cut = " << ptmin/GeV << " GeV) = " << rawjets.size()); 00047 foreach (const Jet& j, rawjets) { 00048 const FourMomentum pj = j.momentum(); 00049 if (!inRange(pj.pT(), ptmin, ptmax)) continue; 00050 if (rapscheme == PSEUDORAPIDITY && !inRange(pj.eta(), rapmin, rapmax)) continue; 00051 if (rapscheme == RAPIDITY && !inRange(pj.rapidity(), rapmin, rapmax)) continue; 00052 rtn += j; 00053 } 00054 return rtn; 00055 } 00056 00057 /// Get the jets, ordered by supplied sorting function object, with optional cuts on \f$ p_\perp \f$ and rapidity. 00058 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00059 template <typename F> 00060 Jets jets(F sorter, double ptmin, double ptmax, 00061 double rapmin, double rapmax, 00062 RapScheme rapscheme) const { 00063 Jets js = jets(ptmin, ptmax, rapmin, rapmax, rapscheme); 00064 if (sorter != 0) { 00065 std::sort(js.begin(), js.end(), sorter); 00066 } 00067 return js; 00068 } 00069 00070 /// Get the jets, ordered by \f$ p_T \f$, with optional cuts on \f$ p_\perp \f$ and rapidity. 00071 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00072 Jets jetsByPt(double ptmin=0.0, double ptmax=MAXDOUBLE, 00073 double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, 00074 RapScheme rapscheme=PSEUDORAPIDITY) const { 00075 return jets(cmpMomByPt, ptmin, ptmax, rapmin, rapmax, rapscheme); 00076 } 00077 00078 /// Get the jets, ordered by \f$ |p| \f$, with optional cuts on \f$ p_\perp \f$ and rapidity. 00079 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00080 Jets jetsByP(double ptmin=0.0, double ptmax=MAXDOUBLE, 00081 double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, 00082 RapScheme rapscheme=PSEUDORAPIDITY) const { 00083 return jets(cmpMomByP, ptmin, ptmax, rapmin, rapmax, rapscheme); 00084 } 00085 00086 /// Get the jets, ordered by \f$ E \f$, with optional cuts on \f$ p_\perp \f$ and rapidity. 00087 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00088 Jets jetsByE(double ptmin=0.0, double ptmax=MAXDOUBLE, 00089 double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, 00090 RapScheme rapscheme=PSEUDORAPIDITY) const { 00091 return jets(cmpMomByE, ptmin, ptmax, rapmin, rapmax, rapscheme); 00092 } 00093 00094 /// Get the jets, ordered by \f$ E_T \f$, with optional cuts on \f$ p_\perp \f$ and rapidity. 00095 /// @todo Introduce MomentumFilter objects for pT, ET, eta, y, etc. filtering, to avoid double-arg ambiguities 00096 Jets jetsByEt(double ptmin=0.0, double ptmax=MAXDOUBLE, 00097 double rapmin=-MAXDOUBLE, double rapmax=MAXDOUBLE, 00098 RapScheme rapscheme=PSEUDORAPIDITY) const { 00099 return jets(cmpMomByEt, ptmin, ptmax, rapmin, rapmax, rapscheme); 00100 } 00101 00102 00103 protected: 00104 00105 /// @brief Internal pure virtual method for getting jets in no guaranteed order. 00106 /// An optional cut on min \f$ p_\perp \f$ is applied in this function, since that is 00107 /// directly supported by FastJet and it seems a shame to not make use of that. But 00108 /// all other jet cuts are applied at the @c ::jets() function level. 00109 virtual Jets _jets(double ptmin) const = 0; 00110 00111 00112 public: 00113 00114 /// Number of jets. 00115 virtual size_t size() const = 0; 00116 /// Determine if the jet collection is empty. 00117 bool empty() const { return size() != 0; } 00118 00119 /// Clear the projection. 00120 virtual void reset() = 0; 00121 00122 typedef Jet entity_type; 00123 typedef Jets collection_type; 00124 00125 /// Template-usable interface common to FinalState. 00126 collection_type entities() const { return jets(); } 00127 00128 /// Do the calculation locally (no caching). 00129 virtual void calc(const Particles& ps) = 0; 00130 00131 00132 protected: 00133 00134 /// Perform the projection on the Event. 00135 virtual void project(const Event& e) = 0; 00136 00137 /// Compare projections. 00138 virtual int compare(const Projection& p) const = 0; 00139 00140 00141 protected: 00142 00143 /// Flag to determine whether or not the VFS wrapper is to be used. 00144 bool _useInvisibles; 00145 00146 }; 00147 00148 00149 } 00150 00151 #endif Generated on Tue May 13 2014 11:38:28 for The Rivet MC analysis system by ![]() |