CMS 3D CMS Logo

Qjets.h
Go to the documentation of this file.
1 #ifndef RecoJets_JetAlgorithms_QJets_h
2 #define RecoJets_JetAlgorithms_QJets_h
3 #include <queue>
4 #include <vector>
5 #include <list>
6 #include <algorithm>
7 #include "fastjet/JetDefinition.hh"
8 #include "fastjet/PseudoJet.hh"
9 #include "fastjet/ClusterSequence.hh"
14 #include "CLHEP/Random/RandomEngine.h"
15 
16 struct JetDistance {
17  double dij;
18  int j1;
19  int j2;
20 };
21 
23 public:
25  bool operator()(const JetDistance& lhs, const JetDistance& rhs) const { return lhs.dij > rhs.dij; };
26 };
27 
28 class Qjets {
29 private:
30  double omega;
32  unsigned int _seed;
34  std::map<int, bool> _merged_jets;
35  std::priority_queue<JetDistance, std::vector<JetDistance>, JetDistanceCompare> _distances;
36  CLHEP::HepRandomEngine* _rnEngine;
37 
38  double d_ij(const fastjet::PseudoJet& v1, const fastjet::PseudoJet& v2) const;
39  void computeDCut(fastjet::ClusterSequence& cs);
40 
41  double Rand();
42  bool Prune(JetDistance& jd, fastjet::ClusterSequence& cs);
43  bool JetsUnmerged(const JetDistance& jd) const;
44  bool JetUnmerged(int num) const;
45  void ComputeNewDistanceMeasures(fastjet::ClusterSequence& cs, unsigned int new_jet);
46  void ComputeAllDistances(const std::vector<fastjet::PseudoJet>& inp);
47  double ComputeMinimumDistance();
48  double ComputeNormalization(double dmin);
50  bool Same(const JetDistance& lhs, const JetDistance& rhs);
51 
52 public:
53  Qjets(double zcut,
54  double dcut_fctr,
55  double exp_min,
56  double exp_max,
57  double rigidity,
58  double truncation_fctr,
59  CLHEP::HepRandomEngine* rnEngine)
61  _zcut(zcut),
62  _dcut(-1.),
63  _dcut_fctr(dcut_fctr),
64  _exp_min(exp_min),
65  _exp_max(exp_max),
67  _truncation_fctr(truncation_fctr),
68  _rnEngine(rnEngine){};
69 
70  void Cluster(fastjet::ClusterSequence& cs);
71  void SetRandSeed(unsigned int seed); /* In case you want reproducible behavior */
72 };
73 
75 public:
76  QjetsBaseExtras() : _wij(-1.) {}
77  ~QjetsBaseExtras() override {}
78  virtual double weight() const { return _wij; }
79  friend class Qjets;
80 
81 protected:
82  double _wij;
83 };
84 
85 #endif
QjetsBaseExtras::_wij
double _wij
Definition: Qjets.h:82
Qjets::ComputeNormalization
double ComputeNormalization(double dmin)
JetDistanceCompare::operator()
bool operator()(const JetDistance &lhs, const JetDistance &rhs) const
Definition: Qjets.h:25
Qjets::SetRandSeed
void SetRandSeed(unsigned int seed)
Definition: Qjets.cc:5
MessageLogger.h
Qjets::Prune
bool Prune(JetDistance &jd, fastjet::ClusterSequence &cs)
Definition: Qjets.cc:120
funct::false
false
Definition: Factorize.h:29
fwrapper::cs
unique_ptr< ClusterSequence > cs
Definition: fastjetfortran_madfks.cc:47
Qjets::JetsUnmerged
bool JetsUnmerged(const JetDistance &jd) const
Definition: Qjets.cc:12
Qjets::_dcut
double _dcut
Definition: Qjets.h:33
RandomNumberGenerator.h
JetDistanceCompare::JetDistanceCompare
JetDistanceCompare()
Definition: Qjets.h:24
align::Extras
Definition: StructureType.h:94
Qjets::_exp_max
double _exp_max
Definition: Qjets.h:33
Qjets::_zcut
double _zcut
Definition: Qjets.h:33
Qjets::_dcut_fctr
double _dcut_fctr
Definition: Qjets.h:33
Qjets::_seed
unsigned int _seed
Definition: Qjets.h:32
Qjets::computeDCut
void computeDCut(fastjet::ClusterSequence &cs)
Definition: Qjets.cc:112
Qjets::_merged_jets
std::map< int, bool > _merged_jets
Definition: Qjets.h:34
fileCollector.seed
seed
Definition: fileCollector.py:127
Qjets::omega
double omega
Definition: Qjets.h:30
qjetsadder_cfi.rigidity
rigidity
Definition: qjetsadder_cfi.py:9
Service.h
JetDistance::j2
int j2
Definition: Qjets.h:19
JetDistanceCompare
Definition: Qjets.h:22
JetDistance::dij
double dij
Definition: Qjets.h:17
Qjets::GetNextDistance
JetDistance GetNextDistance()
Definition: Qjets.cc:14
Qjets::_rand_seed_set
bool _rand_seed_set
Definition: Qjets.h:31
Qjets::_rigidity
double _rigidity
Definition: Qjets.h:33
Qjets::d_ij
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
Definition: Qjets.cc:157
Qjets::Qjets
Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr, CLHEP::HepRandomEngine *rnEngine)
Definition: Qjets.h:53
Qjets::_exp_min
double _exp_min
Definition: Qjets.h:33
Qjets::_truncation_fctr
double _truncation_fctr
Definition: Qjets.h:33
QjetsBaseExtras
Definition: Qjets.h:74
QjetsBaseExtras::weight
virtual double weight() const
Definition: Qjets.h:78
JetDistance::j1
int j1
Definition: Qjets.h:18
QjetsBaseExtras::~QjetsBaseExtras
~QjetsBaseExtras() override
Definition: Qjets.h:77
JetDistance
Definition: Qjets.h:16
Qjets::ComputeMinimumDistance
double ComputeMinimumDistance()
Qjets::_rnEngine
CLHEP::HepRandomEngine * _rnEngine
Definition: Qjets.h:36
HLT_FULL_cff.zcut
zcut
Definition: HLT_FULL_cff.py:8647
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:33
QjetsBaseExtras::QjetsBaseExtras
QjetsBaseExtras()
Definition: Qjets.h:76
isFinite.h
Qjets
Definition: Qjets.h:28
Qjets::JetUnmerged
bool JetUnmerged(int num) const
Definition: Qjets.cc:10
Qjets::Rand
double Rand()
Definition: Qjets.cc:166
Qjets::ComputeNewDistanceMeasures
void ComputeNewDistanceMeasures(fastjet::ClusterSequence &cs, unsigned int new_jet)
Definition: Qjets.cc:145
Qjets::ComputeAllDistances
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
Definition: Qjets.cc:130
Qjets::Cluster
void Cluster(fastjet::ClusterSequence &cs)
Definition: Qjets.cc:60
Qjets::_distances
std::priority_queue< JetDistance, std::vector< JetDistance >, JetDistanceCompare > _distances
Definition: Qjets.h:35
Qjets::Same
bool Same(const JetDistance &lhs, const JetDistance &rhs)