ATLAS Offline Software
Loading...
Searching...
No Matches
Qjets.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef jetsubstructureutils_qjets_header
6#define jetsubstructureutils_qjets_header
7
8#include <queue>
9#include <vector>
10#include <list>
11#include <algorithm>
12#include "fastjet/JetDefinition.hh"
13#include "fastjet/PseudoJet.hh"
14#include "fastjet/ClusterSequence.hh"
15#include "TRandom3.h"
16
17
18namespace JetSubStructureUtils {
20 double dij;
21 int j1;
22 int j2;
23 };
24
26 public:
28 bool operator() (const jet_distance& lhs, const jet_distance &rhs) const {return lhs.dij > rhs.dij;};
29 };
30
31 class Qjets {
32 private:
33 TRandom3 m_rand;
35 std::map<int,bool> m_merged_jets;
36 std::priority_queue <jet_distance, std::vector<jet_distance>, JetDistanceCompare> m_distances;
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(jet_distance& jd,fastjet::ClusterSequence & cs);
43 bool JetsUnmerged(const jet_distance& jd) const;
44 bool JetUnmerged(int num) const;
45 void ComputeNewDistanceMeasures(fastjet::ClusterSequence & cs, int new_jet);
46 void ComputeAllDistances(const std::vector<fastjet::PseudoJet>& inp);
48 double ComputeNormalization(double dmin);
50 bool Same(const jet_distance& lhs, const jet_distance& rhs);
51
52 public:
53 Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr);
54 void Cluster(fastjet::ClusterSequence & cs);
55 void SetRandSeed(unsigned int seed); /* In case you want reproducible behavior */
56 };
57}
58#endif
bool operator()(const jet_distance &lhs, const jet_distance &rhs) const
Definition Qjets.h:28
void ComputeNewDistanceMeasures(fastjet::ClusterSequence &cs, int new_jet)
Definition Qjets.cxx:159
double d_ij(const fastjet::PseudoJet &v1, const fastjet::PseudoJet &v2) const
Definition Qjets.cxx:171
void SetRandSeed(unsigned int seed)
Definition Qjets.cxx:22
void ComputeAllDistances(const std::vector< fastjet::PseudoJet > &inp)
Definition Qjets.cxx:144
bool Same(const jet_distance &lhs, const jet_distance &rhs)
bool JetsUnmerged(const jet_distance &jd) const
Definition Qjets.cxx:30
std::priority_queue< jet_distance, std::vector< jet_distance >, JetDistanceCompare > m_distances
Definition Qjets.h:36
bool JetUnmerged(int num) const
Definition Qjets.cxx:26
double ComputeNormalization(double dmin)
jet_distance GetNextDistance()
Definition Qjets.cxx:34
std::map< int, bool > m_merged_jets
Definition Qjets.h:35
void ComputeDCut(fastjet::ClusterSequence &cs)
Definition Qjets.cxx:126
void Cluster(fastjet::ClusterSequence &cs)
Definition Qjets.cxx:80
bool Prune(jet_distance &jd, fastjet::ClusterSequence &cs)
Definition Qjets.cxx:134
Qjets(double zcut, double dcut_fctr, double exp_min, double exp_max, double rigidity, double truncation_fctr)
Definition Qjets.cxx:10