File indexing completed on 2024-06-13 03:24:04
0001 #include "RecoHGCal/TICL/interface/TracksterLinkingAlgoBase.h"
0002 #include "fastjet/ClusterSequence.hh"
0003 #include "DataFormats/Math/interface/deltaR.h"
0004 #include "RecoHGCal/TICL/plugins/TracksterLinkingbyFastJet.h"
0005
0006 using namespace ticl;
0007
0008 void TracksterLinkingbyFastJet::linkTracksters(
0009 const Inputs& input,
0010 std::vector<Trackster>& resultTracksters,
0011 std::vector<std::vector<unsigned int>>& linkedResultTracksters,
0012 std::vector<std::vector<unsigned int>>& linkedTracksterIdToInputTracksterId) {
0013
0014 std::vector<fastjet::PseudoJet> fjInputs;
0015 for (size_t i = 0; i < input.tracksters.size(); ++i) {
0016
0017 fastjet::PseudoJet pj(input.tracksters[i].barycenter().x(),
0018 input.tracksters[i].barycenter().y(),
0019 input.tracksters[i].barycenter().z(),
0020 input.tracksters[i].raw_energy());
0021 pj.set_user_index(i);
0022 fjInputs.push_back(pj);
0023 }
0024
0025
0026 fastjet::ClusterSequence sequence(fjInputs, fastjet::JetDefinition(algorithm_, radius_));
0027 auto jets = fastjet::sorted_by_pt(sequence.inclusive_jets(0));
0028 linkedTracksterIdToInputTracksterId.resize(jets.size());
0029
0030 for (unsigned int i = 0; i < jets.size(); ++i) {
0031 const auto& jet = jets[i];
0032
0033 std::vector<unsigned int> linkedTracksters;
0034 Trackster outTrackster;
0035 if (!jet.constituents().empty()) {
0036
0037 for (const auto& constituent : jet.constituents()) {
0038 auto tracksterIndex = constituent.user_index();
0039 linkedTracksterIdToInputTracksterId[i].push_back(tracksterIndex);
0040 outTrackster.mergeTracksters(input.tracksters[tracksterIndex]);
0041 }
0042 linkedTracksters.push_back(resultTracksters.size());
0043 resultTracksters.push_back(outTrackster);
0044
0045 linkedResultTracksters.push_back(linkedTracksters);
0046 }
0047 }
0048 }