CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
ticl::TracksterLinkingbyFastJet Class Reference

#include <TracksterLinkingbyFastJet.h>

Inheritance diagram for ticl::TracksterLinkingbyFastJet:
ticl::TracksterLinkingAlgoBase

Public Member Functions

void initialize (const HGCalDDDConstants *hgcons, const hgcal::RecHitTools rhtools, const edm::ESHandle< MagneticField > bfieldH, const edm::ESHandle< Propagator > propH) override
 
void linkTracksters (const Inputs &input, std::vector< Trackster > &resultTracksters, std::vector< std::vector< unsigned int >> &linkedResultTracksters, std::vector< std::vector< unsigned int >> &linkedTracksterIdToInputTracksterId) override
 
 TracksterLinkingbyFastJet (const edm::ParameterSet &conf, edm::ConsumesCollector iC)
 
 ~TracksterLinkingbyFastJet () override
 
- Public Member Functions inherited from ticl::TracksterLinkingAlgoBase
 TracksterLinkingAlgoBase (const edm::ParameterSet &conf, edm::ConsumesCollector)
 
virtual ~TracksterLinkingAlgoBase ()
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &iDesc)
 
- Static Public Member Functions inherited from ticl::TracksterLinkingAlgoBase
static void fillPSetDescription (edm::ParameterSetDescription &desc)
 

Private Attributes

fastjet::JetAlgorithm algorithm_
 
const float radius_
 

Additional Inherited Members

- Protected Attributes inherited from ticl::TracksterLinkingAlgoBase
int algo_verbosity_
 

Detailed Description

Definition at line 14 of file TracksterLinkingbyFastJet.h.

Constructor & Destructor Documentation

◆ TracksterLinkingbyFastJet()

ticl::TracksterLinkingbyFastJet::TracksterLinkingbyFastJet ( const edm::ParameterSet conf,
edm::ConsumesCollector  iC 
)
inline

Definition at line 16 of file TracksterLinkingbyFastJet.h.

References algorithm_, Exception, and edm::ParameterSet::getParameter().

17  : TracksterLinkingAlgoBase(conf, iC), radius_(conf.getParameter<double>("radius")) {
18  // Cluster tracksters into jets using FastJet with configurable algorithm
19  auto algo = conf.getParameter<int>("jet_algorithm");
20 
21  switch (algo) {
22  case 0:
23  algorithm_ = fastjet::kt_algorithm;
24  break;
25  case 1:
26  algorithm_ = fastjet::cambridge_algorithm;
27  break;
28  case 2:
29  algorithm_ = fastjet::antikt_algorithm;
30  break;
31  default:
32  throw cms::Exception("BadConfig") << "FastJet jet clustering algorithm not set correctly.";
33  }
34  }
TracksterLinkingAlgoBase(const edm::ParameterSet &conf, edm::ConsumesCollector)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307

◆ ~TracksterLinkingbyFastJet()

ticl::TracksterLinkingbyFastJet::~TracksterLinkingbyFastJet ( )
inlineoverride

Definition at line 36 of file TracksterLinkingbyFastJet.h.

36 {}

Member Function Documentation

◆ fillPSetDescription()

static void ticl::TracksterLinkingbyFastJet::fillPSetDescription ( edm::ParameterSetDescription iDesc)
inlinestatic

Definition at line 48 of file TracksterLinkingbyFastJet.h.

References edm::ParameterSetDescription::add().

48  {
49  iDesc.add<int>("algo_verbosity", 0);
50  iDesc.add<int>("jet_algorithm", 2)
51  ->setComment("FastJet jet clustering algorithm: 0 = kt, 1 = Cambridge/Aachen, 2 = anti-kt");
52  iDesc.add<double>("radius", 0.1);
53  }
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ initialize()

void ticl::TracksterLinkingbyFastJet::initialize ( const HGCalDDDConstants hgcons,
const hgcal::RecHitTools  rhtools,
const edm::ESHandle< MagneticField bfieldH,
const edm::ESHandle< Propagator propH 
)
inlineoverridevirtual

Implements ticl::TracksterLinkingAlgoBase.

Definition at line 43 of file TracksterLinkingbyFastJet.h.

46  {};

◆ linkTracksters()

void TracksterLinkingbyFastJet::linkTracksters ( const Inputs input,
std::vector< Trackster > &  resultTracksters,
std::vector< std::vector< unsigned int >> &  linkedResultTracksters,
std::vector< std::vector< unsigned int >> &  linkedTracksterIdToInputTracksterId 
)
overridevirtual

Implements ticl::TracksterLinkingAlgoBase.

Definition at line 8 of file TracksterLinkingbyFastJet.cc.

References algorithm_, mps_fire::i, input, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, ticl::Trackster::mergeTracksters(), and radius_.

12  {
13  // Create jets of tracksters using FastJet
14  std::vector<fastjet::PseudoJet> fjInputs;
15  for (size_t i = 0; i < input.tracksters.size(); ++i) {
16  // Convert Trackster information to PseudoJet
17  fastjet::PseudoJet pj(input.tracksters[i].barycenter().x(),
18  input.tracksters[i].barycenter().y(),
19  input.tracksters[i].barycenter().z(),
20  input.tracksters[i].raw_energy());
21  pj.set_user_index(i);
22  fjInputs.push_back(pj);
23  }
24 
25  // Cluster tracksters into jets using FastJet
26  fastjet::ClusterSequence sequence(fjInputs, fastjet::JetDefinition(algorithm_, radius_));
27  auto jets = fastjet::sorted_by_pt(sequence.inclusive_jets(0));
28  linkedTracksterIdToInputTracksterId.resize(jets.size());
29  // Link tracksters based on which ones are components of the same jet
30  for (unsigned int i = 0; i < jets.size(); ++i) {
31  const auto& jet = jets[i];
32 
33  std::vector<unsigned int> linkedTracksters;
34  Trackster outTrackster;
35  if (!jet.constituents().empty()) {
36  // Check if a trackster is a component of the current jet
37  for (const auto& constituent : jet.constituents()) {
38  auto tracksterIndex = constituent.user_index();
39  linkedTracksterIdToInputTracksterId[i].push_back(tracksterIndex);
40  outTrackster.mergeTracksters(input.tracksters[tracksterIndex]);
41  }
42  linkedTracksters.push_back(resultTracksters.size());
43  resultTracksters.push_back(outTrackster);
44  // Store the linked tracksters
45  linkedResultTracksters.push_back(linkedTracksters);
46  }
47  }
48 }
static std::string const input
Definition: EdmProvDump.cc:50
void mergeTracksters(const Trackster &other)
Definition: Trackster.h:81

Member Data Documentation

◆ algorithm_

fastjet::JetAlgorithm ticl::TracksterLinkingbyFastJet::algorithm_
private

Definition at line 56 of file TracksterLinkingbyFastJet.h.

Referenced by linkTracksters(), and TracksterLinkingbyFastJet().

◆ radius_

const float ticl::TracksterLinkingbyFastJet::radius_
private

Definition at line 57 of file TracksterLinkingbyFastJet.h.

Referenced by linkTracksters().